MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwf-rch.f90
Go to the documentation of this file.
1 module rchmodule
2  !
3  use kindmodule, only: dp, i4b, lgp
7  use bndmodule, only: bndtype
8  use bndextmodule, only: bndexttype
10  use simvariablesmodule, only: errmsg
16  use geomutilmodule, only: get_node
17  !
18  implicit none
19  !
20  private
21  public :: rch_create
22  !
23  character(len=LENFTYPE) :: ftype = 'RCH'
24  character(len=LENPACKAGENAME) :: text = ' RCH'
25  character(len=LENPACKAGENAME) :: texta = ' RCHA'
26  !
27  type, extends(bndexttype) :: rchtype
28  real(dp), dimension(:), pointer, contiguous :: recharge => null() !< boundary recharge array
29  integer(I4B), dimension(:), pointer, contiguous :: nodesontop => null() ! User provided cell numbers; nodelist is cells where recharge is applied)
30  logical, pointer, private :: fixed_cell
31  logical, pointer, private :: read_as_arrays
32 
33  contains
34 
35  procedure :: rch_allocate_scalars
36  procedure :: allocate_arrays => rch_allocate_arrays
37  procedure :: source_options => rch_source_options
38  procedure :: source_dimensions => rch_source_dimensions
39  procedure :: log_rch_options
40  procedure :: read_initial_attr => rch_read_initial_attr
41  procedure :: bnd_rp => rch_rp
42  procedure :: bnd_cf => rch_cf
43  procedure :: bnd_fc => rch_fc
44  procedure :: bnd_da => rch_da
45  procedure :: set_nodesontop
46  procedure :: define_listlabel => rch_define_listlabel
47  procedure :: bound_value => rch_bound_value
48  procedure, private :: default_nodelist
49  ! -- for observations
50  procedure, public :: bnd_obs_supported => rch_obs_supported
51  procedure, public :: bnd_df_obs => rch_df_obs
52 
53  end type rchtype
54 
55 contains
56 
57  !> @brief Create a New Recharge Package
58  !!
59  !! Create new RCH package and point packobj to the new package
60  !<
61  subroutine rch_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
62  mempath)
63  ! -- dummy
64  class(bndtype), pointer :: packobj
65  integer(I4B), intent(in) :: id
66  integer(I4B), intent(in) :: ibcnum
67  integer(I4B), intent(in) :: inunit
68  integer(I4B), intent(in) :: iout
69  character(len=*), intent(in) :: namemodel
70  character(len=*), intent(in) :: pakname
71  character(len=*), intent(in) :: mempath
72  ! -- local
73  type(rchtype), pointer :: rchobj
74  !
75  ! -- allocate recharge object and scalar variables
76  allocate (rchobj)
77  packobj => rchobj
78  !
79  ! -- create name and memory path
80  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
81  packobj%text = text
82  !
83  ! -- allocate scalars
84  call rchobj%rch_allocate_scalars()
85  !
86  ! -- initialize package
87  call packobj%pack_initialize()
88  !
89  packobj%inunit = inunit
90  packobj%iout = iout
91  packobj%id = id
92  packobj%ibcnum = ibcnum
93  packobj%ictMemPath = create_mem_path(namemodel, 'NPF')
94  !
95  ! -- Return
96  return
97  end subroutine rch_create
98 
99  !> @brief Allocate scalar members
100  !<
101  subroutine rch_allocate_scalars(this)
102  ! -- dummy
103  class(rchtype), intent(inout) :: this
104  !
105  ! -- allocate base scalars
106  call this%BndExtType%allocate_scalars()
107  !
108  ! -- allocate internal members
109  allocate (this%fixed_cell)
110  allocate (this%read_as_arrays)
111  !
112  ! -- Set values
113  this%fixed_cell = .false.
114  this%read_as_arrays = .false.
115  !
116  ! -- Return
117  return
118  end subroutine rch_allocate_scalars
119 
120  !> @brief Allocate package arrays
121  !<
122  subroutine rch_allocate_arrays(this, nodelist, auxvar)
123  ! -- modules
125  ! -- dummy
126  class(rchtype) :: this
127  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
128  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
129  !
130  ! -- allocate base arrays
131  call this%BndExtType%allocate_arrays(nodelist, auxvar)
132  !
133  ! -- set recharge input context pointer
134  call mem_setptr(this%recharge, 'RECHARGE', this%input_mempath)
135  !
136  ! -- checkin recharge input context pointer
137  call mem_checkin(this%recharge, 'RECHARGE', this%memoryPath, &
138  'RECHARGE', this%input_mempath)
139  !
140  ! -- Return
141  return
142  end subroutine rch_allocate_arrays
143 
144  !> @brief Source options specific to RchType
145  !<
146  subroutine rch_source_options(this)
147  ! -- modules
149  implicit none
150  ! -- dummy
151  class(rchtype), intent(inout) :: this
152  ! -- local
153  logical(LGP) :: found_fixed_cell = .false.
154  logical(LGP) :: found_readasarrays = .false.
155  !
156  ! -- source common bound options
157  call this%BndExtType%source_options()
158  !
159  ! -- update defaults with idm sourced values
160  call mem_set_value(this%fixed_cell, 'FIXED_CELL', this%input_mempath, &
161  found_fixed_cell)
162  call mem_set_value(this%read_as_arrays, 'READASARRAYS', this%input_mempath, &
163  found_readasarrays)
164  !
165  if (found_readasarrays) then
166  if (this%dis%supports_layers()) then
167  this%text = texta
168  else
169  errmsg = 'READASARRAYS option is not compatible with selected'// &
170  ' discretization type.'
171  call store_error(errmsg)
172  call store_error_filename(this%input_fname)
173  end if
174  end if
175  !
176  ! -- log rch params
177  call this%log_rch_options(found_fixed_cell, found_readasarrays)
178  !
179  ! -- Return
180  return
181  end subroutine rch_source_options
182 
183  !> @brief Log options specific to RchType
184  !<
185  subroutine log_rch_options(this, found_fixed_cell, found_readasarrays)
186  implicit none
187  ! -- dummy
188  class(rchtype), intent(inout) :: this
189  logical(LGP), intent(in) :: found_fixed_cell
190  logical(LGP), intent(in) :: found_readasarrays
191  ! -- formats
192  character(len=*), parameter :: fmtfixedcell = &
193  &"(4x, 'RECHARGE WILL BE APPLIED TO SPECIFIED CELL.')"
194  character(len=*), parameter :: fmtreadasarrays = &
195  &"(4x, 'RECHARGE INPUT WILL BE READ AS ARRAY(S).')"
196  !
197  ! -- log found options
198  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
199  //' OPTIONS'
200  !
201  if (found_fixed_cell) then
202  write (this%iout, fmtfixedcell)
203  end if
204  !
205  if (found_readasarrays) then
206  write (this%iout, fmtreadasarrays)
207  end if
208  !
209  ! -- close logging block
210  write (this%iout, '(1x,a)') &
211  'END OF '//trim(adjustl(this%text))//' OPTIONS'
212  !
213  ! -- Return
214  return
215  end subroutine log_rch_options
216 
217  !> @brief Source the dimensions for this package
218  !!
219  !! Dimensions block is not required if:
220  !! (1) discretization is DIS or DISV, and
221  !! (2) READASARRAYS option has been specified.
222  !<
223  subroutine rch_source_dimensions(this)
224  ! -- dummy
225  class(rchtype), intent(inout) :: this
226  !
227  if (this%read_as_arrays) then
228  this%maxbound = this%dis%get_ncpl()
229  !
230  ! -- verify dimensions were set
231  if (this%maxbound <= 0) then
232  write (errmsg, '(a)') &
233  'MAXBOUND must be an integer greater than zero.'
234  call store_error(errmsg)
235  call store_error_filename(this%input_fname)
236  end if
237  !
238  else
239  !
240  ! -- source maxbound
241  call this%BndExtType%source_dimensions()
242  end if
243  !
244  ! -- Call define_listlabel to construct the list label that is written
245  ! when PRINT_INPUT option is used.
246  call this%define_listlabel()
247  !
248  ! -- Return
249  return
250  end subroutine rch_source_dimensions
251 
252  !> @brief Part of allocate and read
253  !<
254  subroutine rch_read_initial_attr(this)
255  ! -- dummy
256  class(rchtype), intent(inout) :: this
257  !
258  if (this%read_as_arrays) then
259  call this%default_nodelist()
260  end if
261  !
262  ! -- Return
263  return
264  end subroutine rch_read_initial_attr
265 
266  !> @brief Read and Prepare
267  !!
268  !! Read itmp and read new boundaries if itmp > 0
269  !<
270  subroutine rch_rp(this)
271  ! -- modules
272  use tdismodule, only: kper
273  implicit none
274  ! -- dummy
275  class(rchtype), intent(inout) :: this
276  !
277  if (this%iper /= kper) return
278  !
279  if (this%read_as_arrays) then
280  !
281  ! -- update nodelist based on IRCH input
282  call nodelist_update(this%nodelist, this%nbound, this%maxbound, &
283  this%dis, this%input_mempath)
284  !
285  else
286  !
287  call this%BndExtType%bnd_rp()
288  !
289  end if
290  !
291  ! -- copy nodelist to nodesontop if not fixed cell
292  if (.not. this%fixed_cell) call this%set_nodesontop()
293  !
294  ! -- Write the list to iout if requested
295  if (this%iprpak /= 0) then
296  call this%write_list()
297  end if
298  !
299  ! -- Return
300  return
301  end subroutine rch_rp
302 
303  !> @brief Store nodelist in nodesontop
304  !<
305  subroutine set_nodesontop(this)
306  implicit none
307  ! -- dummy
308  class(rchtype), intent(inout) :: this
309  ! -- local
310  integer(I4B) :: n
311  !
312  ! -- allocate if necessary
313  if (.not. associated(this%nodesontop)) then
314  allocate (this%nodesontop(this%maxbound))
315  end if
316  !
317  ! -- copy nodelist into nodesontop
318  do n = 1, this%nbound
319  this%nodesontop(n) = this%nodelist(n)
320  end do
321  !
322  ! -- Return
323  return
324  end subroutine set_nodesontop
325 
326  !> @brief Formulate the HCOF and RHS terms
327  !!
328  !! Skip if no recharge. Otherwise, calculate hcof and rhs
329  !<
330  subroutine rch_cf(this)
331  ! -- dummy
332  class(rchtype) :: this
333  ! -- local
334  integer(I4B) :: i, node
335  !
336  ! -- Return if no recharge
337  if (this%nbound == 0) return
338  !
339  ! -- Calculate hcof and rhs for each recharge entry
340  do i = 1, this%nbound
341  !
342  ! -- Find the node number
343  if (this%fixed_cell) then
344  node = this%nodelist(i)
345  else
346  node = this%nodesontop(i)
347  end if
348  !
349  ! -- cycle if nonexistent bound
350  if (node <= 0) then
351  this%hcof(i) = dzero
352  this%rhs(i) = dzero
353  cycle
354  end if
355  !
356  ! -- reset nodelist to highest active
357  if (.not. this%fixed_cell) then
358  if (this%ibound(node) == 0) &
359  call this%dis%highest_active(node, this%ibound)
360  this%nodelist(i) = node
361  end if
362  !
363  ! -- Set rhs and hcof
364  this%hcof(i) = dzero
365  if (this%iauxmultcol > 0) then
366  this%rhs(i) = -this%recharge(i) * this%dis%get_area(node) * &
367  this%auxvar(this%iauxmultcol, i)
368  else
369  this%rhs(i) = -this%recharge(i) * this%dis%get_area(node)
370  end if
371  if (this%ibound(node) <= 0) then
372  this%rhs(i) = dzero
373  cycle
374  end if
375  if (this%ibound(node) == iwetlake) then
376  this%rhs(i) = dzero
377  cycle
378  end if
379  end do
380  !
381  ! -- Return
382  return
383  end subroutine rch_cf
384 
385  !> @brief Copy rhs and hcof into solution rhs and amat
386  !<
387  subroutine rch_fc(this, rhs, ia, idxglo, matrix_sln)
388  ! -- dummy
389  class(rchtype) :: this
390  real(DP), dimension(:), intent(inout) :: rhs
391  integer(I4B), dimension(:), intent(in) :: ia
392  integer(I4B), dimension(:), intent(in) :: idxglo
393  class(matrixbasetype), pointer :: matrix_sln
394  ! -- local
395  integer(I4B) :: i, n, ipos
396  !
397  ! -- Copy package rhs and hcof into solution rhs and amat
398  do i = 1, this%nbound
399  n = this%nodelist(i)
400  if (n <= 0) cycle
401  ! -- reset hcof and rhs for excluded cells
402  if (this%ibound(n) == iwetlake) then
403  this%hcof(i) = dzero
404  this%rhs(i) = dzero
405  cycle
406  end if
407  rhs(n) = rhs(n) + this%rhs(i)
408  ipos = ia(n)
409  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
410  end do
411  !
412  ! -- Return
413  return
414  end subroutine rch_fc
415 
416  !> @brief Deallocate memory
417  !<
418  subroutine rch_da(this)
419  ! -- modules
421  ! -- dummy
422  class(rchtype) :: this
423  !
424  ! -- Deallocate parent package
425  call this%BndExtType%bnd_da()
426  !
427  ! -- scalars
428  deallocate (this%fixed_cell)
429  deallocate (this%read_as_arrays)
430  !
431  ! -- arrays
432  if (associated(this%nodesontop)) deallocate (this%nodesontop)
433  call mem_deallocate(this%recharge, 'RECHARGE', this%memoryPath)
434  !
435  ! -- Return
436  return
437  end subroutine rch_da
438 
439  !> @brief Define the list heading that is written to iout when PRINT_INPUT
440  !! option is used.
441  !<
442  subroutine rch_define_listlabel(this)
443  ! -- dummy
444  class(rchtype), intent(inout) :: this
445  !
446  ! -- create the header list label
447  this%listlabel = trim(this%filtyp)//' NO.'
448  if (this%dis%ndim == 3) then
449  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
450  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
451  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
452  elseif (this%dis%ndim == 2) then
453  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
454  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
455  else
456  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
457  end if
458  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'RECHARGE'
459 ! if(this%multindex > 0) &
460 ! write(this%listlabel, '(a, a16)') trim(this%listlabel), 'MULTIPLIER'
461  if (this%inamedbound == 1) then
462  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
463  end if
464  !
465  ! -- Return
466  return
467  end subroutine rch_define_listlabel
468 
469  !> @brief Assign default nodelist when READASARRAYS is specified.
470  !!
471  !! Equivalent to reading IRCH as CONSTANT 1
472  !<
473  subroutine default_nodelist(this)
474  ! -- dummy
475  class(rchtype) :: this
476  ! -- local
477  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nodeu, noder, ipos
478  !
479  ! -- set variables
480  if (this%dis%ndim == 3) then
481  nlay = this%dis%mshape(1)
482  nrow = this%dis%mshape(2)
483  ncol = this%dis%mshape(3)
484  elseif (this%dis%ndim == 2) then
485  nlay = this%dis%mshape(1)
486  nrow = 1
487  ncol = this%dis%mshape(2)
488  end if
489  !
490  ! -- Populate nodelist
491  ipos = 1
492  il = 1
493  do ir = 1, nrow
494  do ic = 1, ncol
495  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
496  noder = this%dis%get_nodenumber(nodeu, 0)
497  this%nodelist(ipos) = noder
498  ipos = ipos + 1
499  end do
500  end do
501  !
502  ! -- Assign nbound
503  this%nbound = ipos - 1
504  !
505  ! -- if fixed_cell option not set, then need to store nodelist
506  ! in the nodesontop array
507  if (.not. this%fixed_cell) call this%set_nodesontop()
508  !
509  ! -- Return
510  return
511  end subroutine default_nodelist
512 
513  ! -- Procedures related to observations
514 
515  !> @brief
516  !!
517  !! Overrides BndType%bnd_obs_supported()
518  !<
519  logical function rch_obs_supported(this)
520  implicit none
521  ! -- dummy
522  class(rchtype) :: this
523  !
524  rch_obs_supported = .true.
525  !
526  ! -- Return
527  return
528  end function rch_obs_supported
529 
530  !> @brief Implements bnd_df_obs
531  !!
532  !! Store observation type supported by RCH package. Overrides
533  !! BndType%bnd_df_obs
534  !<
535  subroutine rch_df_obs(this)
536  implicit none
537  ! -- dummy
538  class(rchtype) :: this
539  ! -- local
540  integer(I4B) :: indx
541  !
542  call this%obs%StoreObsType('rch', .true., indx)
543  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
544  !
545  ! -- Return
546  return
547  end subroutine rch_df_obs
548 
549  !> @brief Return requested boundary value
550  !<
551  function rch_bound_value(this, col, row) result(bndval)
552  ! -- modules
553  use constantsmodule, only: dzero
554  ! -- dummy
555  class(rchtype), intent(inout) :: this !< BndExtType object
556  integer(I4B), intent(in) :: col
557  integer(I4B), intent(in) :: row
558  ! -- result
559  real(dp) :: bndval
560  !
561  select case (col)
562  case (1)
563  if (this%iauxmultcol > 0) then
564  bndval = this%recharge(row) * this%auxvar(this%iauxmultcol, row)
565  else
566  bndval = this%recharge(row)
567  end if
568  case default
569  errmsg = 'Programming error. RCH bound value requested column '&
570  &'outside range of ncolbnd (1).'
571  call store_error(errmsg)
572  call store_error_filename(this%input_fname)
573  end select
574  !
575  ! -- Return
576  return
577  end function rch_bound_value
578 
579  !> @brief Update the nodelist based on IRCH input
580  !!
581  !! This is a module scoped routine to check for IRCH
582  !! input. If array input was provided, INIRCH and IRCH
583  !! will be allocated in the input context. If the read
584  !! state variable INIRCH is set to 1 during this period
585  !! update, IRCH input was read and is used here to update
586  !! the nodelist.
587  !!
588  !<
589  subroutine nodelist_update(nodelist, nbound, maxbound, &
590  dis, input_mempath)
591  ! -- modules
593  use basedismodule, only: disbasetype
594  ! -- dummy
595  integer(I4B), dimension(:), contiguous, &
596  pointer, intent(inout) :: nodelist
597  class(disbasetype), pointer, intent(in) :: dis
598  character(len=*), intent(in) :: input_mempath
599  integer(I4B), intent(inout) :: nbound
600  integer(I4B), intent(in) :: maxbound
601  character(len=24) :: aname = ' LAYER OR NODE INDEX'
602  ! -- local
603  integer(I4B), dimension(:), contiguous, &
604  pointer :: irch => null()
605  integer(I4B), pointer :: inirch => null()
606  !
607  ! -- set pointer to input context INIRCH
608  call mem_setptr(inirch, 'INIRCH', input_mempath)
609  !
610  ! -- check INIRCH read state
611  if (inirch == 1) then
612  ! -- irch was read this period
613  !
614  ! -- set pointer to input context IRCH
615  call mem_setptr(irch, 'IRCH', input_mempath)
616  !
617  ! -- update nodelist
618  call dis%nlarray_to_nodelist(irch, nodelist, &
619  maxbound, nbound, aname)
620  end if
621  !
622  ! -- Return
623  return
624  end subroutine nodelist_update
625 
626 end module rchmodule
627 
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:22
integer(i4b), parameter iwetlake
integer constant for a dry lake
Definition: Constants.f90:51
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:38
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:46
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:81
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:249
subroutine rch_allocate_arrays(this, nodelist, auxvar)
Allocate package arrays.
Definition: gwf-rch.f90:123
logical function rch_obs_supported(this)
Overrides BndTypebnd_obs_supported()
Definition: gwf-rch.f90:520
subroutine nodelist_update(nodelist, nbound, maxbound, dis, input_mempath)
Update the nodelist based on IRCH input.
Definition: gwf-rch.f90:591
subroutine rch_df_obs(this)
Implements bnd_df_obs.
Definition: gwf-rch.f90:536
subroutine default_nodelist(this)
Assign default nodelist when READASARRAYS is specified.
Definition: gwf-rch.f90:474
subroutine, public rch_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Recharge Package.
Definition: gwf-rch.f90:63
subroutine rch_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
Definition: gwf-rch.f90:388
real(dp) function rch_bound_value(this, col, row)
Return requested boundary value.
Definition: gwf-rch.f90:552
subroutine log_rch_options(this, found_fixed_cell, found_readasarrays)
Log options specific to RchType.
Definition: gwf-rch.f90:186
subroutine rch_allocate_scalars(this)
Allocate scalar members.
Definition: gwf-rch.f90:102
character(len=lenpackagename) texta
Definition: gwf-rch.f90:25
subroutine rch_read_initial_attr(this)
Part of allocate and read.
Definition: gwf-rch.f90:255
subroutine rch_rp(this)
Read and Prepare.
Definition: gwf-rch.f90:271
subroutine rch_define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
Definition: gwf-rch.f90:443
subroutine rch_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwf-rch.f90:331
subroutine rch_source_options(this)
Source options specific to RchType.
Definition: gwf-rch.f90:147
subroutine rch_da(this)
Deallocate memory.
Definition: gwf-rch.f90:419
character(len=lenpackagename) text
Definition: gwf-rch.f90:24
subroutine set_nodesontop(this)
Store nodelist in nodesontop.
Definition: gwf-rch.f90:306
character(len=lenftype) ftype
Definition: gwf-rch.f90:23
subroutine rch_source_dimensions(this)
Source the dimensions for this package.
Definition: gwf-rch.f90:224
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
@ brief BndType
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23