MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
structarraymodule Module Reference

This module contains the StructArrayModule. More...

Data Types

type  structarraytype
 type for structured array More...
 

Functions/Subroutines

type(structarraytype) function, pointer, public constructstructarray (mf6_input, ncol, nrow, blocknum, mempath, component_mempath)
 constructor for a struct_array More...
 
subroutine, public destructstructarray (struct_array)
 destructor for a struct_array More...
 
subroutine mem_create_vector (this, icol, idt)
 create new vector in StructArrayType More...
 
integer(i4b) function count (this)
 
subroutine set_pointer (sv, sv_target)
 
type(structvectortype) function, pointer get (this, idx)
 
subroutine allocate_int_type (this, sv)
 allocate integer input type More...
 
subroutine allocate_dbl_type (this, sv)
 allocate double input type More...
 
subroutine allocate_charstr_type (this, sv)
 allocate charstr input type More...
 
subroutine allocate_int1d_type (this, sv)
 allocate int1d input type More...
 
subroutine allocate_dbl1d_type (this, sv)
 allocate dbl1d input type More...
 
subroutine load_deferred_vector (this, icol)
 
subroutine memload_vectors (this)
 load deferred vectors into managed memory More...
 
subroutine log_structarray_vars (this, iout)
 log information about the StructArrayType More...
 
subroutine check_reallocate (this)
 reallocate local memory for deferred vectors if necessary More...
 
subroutine read_param (this, parser, sv_col, irow, timeseries, iout, auxcol)
 
integer(i4b) function read_from_parser (this, parser, timeseries, iout, input_name)
 read from the block parser to fill the StructArrayType More...
 
integer(i4b) function read_from_parser_keystring (this, parser, timeseries, nleading, iout, input_name)
 read keystring period block into the StructArrayType More...
 
integer(i4b) function read_from_binary (this, inunit, iout)
 read from binary input to fill the StructArrayType More...
 
subroutine ts_update (this, tsmanager, subcomp_name, iprpak, input_name, auxname_cst, clear_strlocs)
 link time-series strings in this struct array to a tsmanager More...
 

Detailed Description

This module contains the routines for reading a structured list, which consists of a separate vector for each column in the list.

Function/Subroutine Documentation

◆ allocate_charstr_type()

subroutine structarraymodule::allocate_charstr_type ( class(structarraytype this,
type(structvectortype), intent(inout)  sv 
)
private
Parameters
thisStructArrayType

Definition at line 260 of file StructArray.f90.

261  class(StructArrayType) :: this !< StructArrayType
262  type(StructVectorType), intent(inout) :: sv
263  type(CharacterStringType), dimension(:), pointer, contiguous :: charstr1d
264  integer(I4B) :: j
265 
266  if (this%deferred_shape) then
267  allocate (charstr1d(this%deferred_size_init))
268  else
269  call mem_allocate(charstr1d, linelength, this%nrow, &
270  sv%idt%mf6varname, this%mempath)
271  end if
272 
273  do j = 1, this%nrow
274  charstr1d(j) = ''
275  end do
276 
277  sv%memtype = mtype_str
278  sv%charstr1d => charstr1d

◆ allocate_dbl1d_type()

subroutine structarraymodule::allocate_dbl1d_type ( class(structarraytype this,
type(structvectortype), intent(inout)  sv 
)
Parameters
thisStructArrayType

Definition at line 357 of file StructArray.f90.

358  use memorymanagermodule, only: get_isize
359  class(StructArrayType) :: this !< StructArrayType
360  type(StructVectorType), intent(inout) :: sv
361  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
362  integer(I4B), pointer :: naux, nseg, nseg_1
363  integer(I4B) :: nseg1_isize, n, m
364 
365  if (sv%idt%shape == 'NAUX') then
366  call mem_setptr(naux, sv%idt%shape, this%mempath)
367  call mem_allocate(dbl2d, naux, this%nrow, sv%idt%mf6varname, this%mempath)
368 
369  ! initialize
370  do m = 1, this%nrow
371  do n = 1, naux
372  dbl2d(n, m) = dzero
373  end do
374  end do
375 
376  sv%memtype = mtype_dbl2d
377  sv%dbl2d => dbl2d
378  sv%intshape => naux
379  else if (sv%idt%shape == 'NSEG-1') then
380  call mem_setptr(nseg, 'NSEG', this%mempath)
381  call get_isize('NSEG_1', this%mempath, nseg1_isize)
382 
383  if (nseg1_isize < 0) then
384  call mem_allocate(nseg_1, 'NSEG_1', this%mempath)
385  nseg_1 = nseg - 1
386  else
387  call mem_setptr(nseg_1, 'NSEG_1', this%mempath)
388  end if
389 
390  ! allocate
391  call mem_allocate(dbl2d, nseg_1, this%nrow, sv%idt%mf6varname, this%mempath)
392 
393  ! initialize
394  do m = 1, this%nrow
395  do n = 1, nseg_1
396  dbl2d(n, m) = dzero
397  end do
398  end do
399 
400  sv%memtype = mtype_dbl2d
401  sv%dbl2d => dbl2d
402  sv%intshape => nseg_1
403  else
404  errmsg = 'IDM unimplemented. StructArray::allocate_dbl1d_type &
405  & unsupported shape "'//trim(sv%idt%shape)//'".'
406  call store_error(errmsg, terminate=.true.)
407  end if
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
Here is the call graph for this function:

◆ allocate_dbl_type()

subroutine structarraymodule::allocate_dbl_type ( class(structarraytype this,
type(structvectortype), intent(inout)  sv 
)
private
Parameters
thisStructArrayType

Definition at line 233 of file StructArray.f90.

234  class(StructArrayType) :: this !< StructArrayType
235  type(StructVectorType), intent(inout) :: sv
236  real(DP), dimension(:), pointer, contiguous :: dbl1d
237  integer(I4B) :: j, nrow
238 
239  if (this%deferred_shape) then
240  ! shape not known, allocate locally
241  nrow = this%deferred_size_init
242  allocate (dbl1d(this%deferred_size_init))
243  else
244  ! shape known, allocate in managed memory
245  nrow = this%nrow
246  call mem_allocate(dbl1d, this%nrow, sv%idt%mf6varname, this%mempath)
247  end if
248 
249  ! initialize
250  do j = 1, nrow
251  dbl1d(j) = dzero
252  end do
253 
254  sv%memtype = mtype_dbl
255  sv%dbl1d => dbl1d

◆ allocate_int1d_type()

subroutine structarraymodule::allocate_int1d_type ( class(structarraytype this,
type(structvectortype), intent(inout)  sv 
)
private
Parameters
thisStructArrayType

Definition at line 283 of file StructArray.f90.

284  use constantsmodule, only: lenmodelname
287  class(StructArrayType) :: this !< StructArrayType
288  type(StructVectorType), intent(inout) :: sv
289  integer(I4B), dimension(:, :), pointer, contiguous :: int2d
290  type(STLVecInt), pointer :: intvector
291  integer(I4B), pointer :: ncelldim, exgid
292  character(len=LENMEMPATH) :: input_mempath
293  character(len=LENMODELNAME) :: mname
294  type(CharacterStringType), dimension(:), contiguous, &
295  pointer :: charstr1d
296  integer(I4B) :: nrow, n, m
297 
298  if (sv%idt%shape == 'NCELLDIM') then
299  ! if EXCHANGE set to NCELLDIM of appropriate model
300  if (this%mf6_input%component_type == 'EXG') then
301  ! set pointer to EXGID
302  call mem_setptr(exgid, 'EXGID', this%mf6_input%mempath)
303  ! set pointer to appropriate exchange model array
304  input_mempath = create_mem_path('SIM', 'NAM', idm_context)
305  if (sv%idt%tagname == 'CELLIDM1') then
306  call mem_setptr(charstr1d, 'EXGMNAMEA', input_mempath)
307  else if (sv%idt%tagname == 'CELLIDM2') then
308  call mem_setptr(charstr1d, 'EXGMNAMEB', input_mempath)
309  end if
310 
311  ! set the model name
312  mname = charstr1d(exgid)
313 
314  ! set ncelldim pointer
315  input_mempath = create_mem_path(component=mname, context=idm_context)
316  call mem_setptr(ncelldim, sv%idt%shape, input_mempath)
317  else
318  call mem_setptr(ncelldim, sv%idt%shape, this%component_mempath)
319  end if
320 
321  if (this%deferred_shape) then
322  ! shape not known, allocate locally
323  nrow = this%deferred_size_init
324  allocate (int2d(ncelldim, this%deferred_size_init))
325  else
326  ! shape known, allocate in managed memory
327  nrow = this%nrow
328  call mem_allocate(int2d, ncelldim, this%nrow, &
329  sv%idt%mf6varname, this%mempath)
330  end if
331 
332  ! initialize
333  do m = 1, nrow
334  do n = 1, ncelldim
335  int2d(n, m) = izero
336  end do
337  end do
338 
339  sv%memtype = mtype_int2d
340  sv%int2d => int2d
341  sv%intshape => ncelldim
342  else
343  ! allocate intvector object
344  allocate (intvector)
345  ! initialize STLVecInt
346  call intvector%init()
347  sv%memtype = mtype_intvec
348  sv%intvector => intvector
349  sv%size = -1
350  ! set pointer to dynamic shape
351  call mem_setptr(sv%intvector_shape, sv%idt%shape, this%mempath)
352  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenmodelname
maximum length of the model name
Definition: Constants.f90:22
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
Here is the call graph for this function:

◆ allocate_int_type()

subroutine structarraymodule::allocate_int_type ( class(structarraytype this,
type(structvectortype), intent(inout)  sv 
)
private
Parameters
thisStructArrayType

Definition at line 206 of file StructArray.f90.

207  class(StructArrayType) :: this !< StructArrayType
208  type(StructVectorType), intent(inout) :: sv
209  integer(I4B), dimension(:), pointer, contiguous :: int1d
210  integer(I4B) :: j, nrow
211 
212  if (this%deferred_shape) then
213  ! shape not known, allocate locally
214  nrow = this%deferred_size_init
215  allocate (int1d(this%deferred_size_init))
216  else
217  ! shape known, allocate in managed memory
218  nrow = this%nrow
219  call mem_allocate(int1d, this%nrow, sv%idt%mf6varname, this%mempath)
220  end if
221 
222  ! initialize vector values
223  do j = 1, nrow
224  int1d(j) = izero
225  end do
226 
227  sv%memtype = mtype_int
228  sv%int1d => int1d

◆ check_reallocate()

subroutine structarraymodule::check_reallocate ( class(structarraytype this)
private
Parameters
thisStructArrayType

Definition at line 681 of file StructArray.f90.

682  class(StructArrayType) :: this !< StructArrayType
683  integer(I4B) :: i, j, k, newsize
684  integer(I4B), dimension(:), pointer, contiguous :: p_int1d
685  integer(I4B), dimension(:, :), pointer, contiguous :: p_int2d
686  real(DP), dimension(:), pointer, contiguous :: p_dbl1d
687  type(CharacterStringType), dimension(:), pointer, contiguous :: p_charstr1d
688  integer(I4B) :: reallocate_mult
689 
690  ! set growth rate
691  reallocate_mult = 2
692 
693  do j = 1, this%ncol
694  ! reallocate based on memtype
695  select case (this%struct_vectors(j)%memtype)
696  case (mtype_int)
697  ! check if more space needed
698  if (this%nrow > this%struct_vectors(j)%size) then
699  ! calculate new size
700  newsize = this%struct_vectors(j)%size * reallocate_mult
701  ! allocate new vector
702  allocate (p_int1d(newsize))
703 
704  ! copy from old to new
705  do i = 1, this%struct_vectors(j)%size
706  p_int1d(i) = this%struct_vectors(j)%int1d(i)
707  end do
708 
709  ! deallocate old vector
710  deallocate (this%struct_vectors(j)%int1d)
711 
712  ! update struct array object
713  this%struct_vectors(j)%int1d => p_int1d
714  this%struct_vectors(j)%size = newsize
715  end if
716  case (mtype_dbl)
717  if (this%nrow > this%struct_vectors(j)%size) then
718  newsize = this%struct_vectors(j)%size * reallocate_mult
719  allocate (p_dbl1d(newsize))
720 
721  do i = 1, this%struct_vectors(j)%size
722  p_dbl1d(i) = this%struct_vectors(j)%dbl1d(i)
723  end do
724 
725  deallocate (this%struct_vectors(j)%dbl1d)
726 
727  this%struct_vectors(j)%dbl1d => p_dbl1d
728  this%struct_vectors(j)%size = newsize
729  end if
730  !
731  case (mtype_str)
732  if (this%nrow > this%struct_vectors(j)%size) then
733  newsize = this%struct_vectors(j)%size * reallocate_mult
734  allocate (p_charstr1d(newsize))
735 
736  do i = 1, this%struct_vectors(j)%size
737  p_charstr1d(i) = this%struct_vectors(j)%charstr1d(i)
738  call this%struct_vectors(j)%charstr1d(i)%destroy()
739  end do
740 
741  deallocate (this%struct_vectors(j)%charstr1d)
742 
743  this%struct_vectors(j)%charstr1d => p_charstr1d
744  this%struct_vectors(j)%size = newsize
745  end if
746  case (mtype_int2d)
747  if (this%nrow > this%struct_vectors(j)%size) then
748  newsize = this%struct_vectors(j)%size * reallocate_mult
749  allocate (p_int2d(this%struct_vectors(j)%intshape, newsize))
750 
751  do i = 1, this%struct_vectors(j)%size
752  do k = 1, this%struct_vectors(j)%intshape
753  p_int2d(k, i) = this%struct_vectors(j)%int2d(k, i)
754  end do
755  end do
756 
757  deallocate (this%struct_vectors(j)%int2d)
758 
759  this%struct_vectors(j)%int2d => p_int2d
760  this%struct_vectors(j)%size = newsize
761  end if
762  ! TODO: case (6)
763  case default
764  errmsg = 'IDM unimplemented. StructArray::check_reallocate &
765  &unsupported memtype.'
766  call store_error(errmsg, terminate=.true.)
767  end select
768  end do
Here is the call graph for this function:

◆ constructstructarray()

type(structarraytype) function, pointer, public structarraymodule::constructstructarray ( type(modflowinputtype), intent(in)  mf6_input,
integer(i4b), intent(in)  ncol,
integer(i4b), intent(in)  nrow,
integer(i4b), intent(in)  blocknum,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  component_mempath 
)
Parameters
[in]ncolnumber of columns in the StructArrayType
[in]nrownumber of rows in the StructArrayType
[in]blocknumvalid block number or 0
[in]mempathmemory path for storing the vector
Returns
new StructArrayType

Definition at line 78 of file StructArray.f90.

80  type(ModflowInputType), intent(in) :: mf6_input
81  integer(I4B), intent(in) :: ncol !< number of columns in the StructArrayType
82  integer(I4B), intent(in) :: nrow !< number of rows in the StructArrayType
83  integer(I4B), intent(in) :: blocknum !< valid block number or 0
84  character(len=*), intent(in) :: mempath !< memory path for storing the vector
85  character(len=*), intent(in) :: component_mempath
86  type(StructArrayType), pointer :: struct_array !< new StructArrayType
87 
88  ! allocate StructArrayType
89  allocate (struct_array)
90 
91  ! set description of input
92  struct_array%mf6_input = mf6_input
93 
94  ! set number of arrays
95  struct_array%ncol = ncol
96 
97  ! set rows if known or set deferred
98  struct_array%nrow = nrow
99  if (struct_array%nrow == -1) then
100  struct_array%nrow = 0
101  struct_array%deferred_shape = .true.
102  end if
103 
104  ! set blocknum
105  if (blocknum > 0) then
106  struct_array%blocknum = blocknum
107  else
108  struct_array%blocknum = 0
109  end if
110 
111  ! set mempath
112  struct_array%mempath = mempath
113  struct_array%component_mempath = component_mempath
114 
115  ! allocate StructVectorType objects
116  allocate (struct_array%struct_vectors(ncol))
117  allocate (struct_array%startidx(ncol))
118  allocate (struct_array%numcols(ncol))
Here is the caller graph for this function:

◆ count()

integer(i4b) function structarraymodule::count ( class(structarraytype this)
private
Parameters
thisStructArrayType

Definition at line 185 of file StructArray.f90.

186  class(StructArrayType) :: this !< StructArrayType
187  integer(I4B) :: count
188  count = size(this%struct_vectors)

◆ destructstructarray()

subroutine, public structarraymodule::destructstructarray ( type(structarraytype), intent(inout), pointer  struct_array)
Parameters
[in,out]struct_arrayStructArrayType to destroy

Definition at line 123 of file StructArray.f90.

124  type(StructArrayType), pointer, intent(inout) :: struct_array !< StructArrayType to destroy
125  deallocate (struct_array%struct_vectors)
126  deallocate (struct_array%startidx)
127  deallocate (struct_array%numcols)
128  deallocate (struct_array)
129  nullify (struct_array)
Here is the caller graph for this function:

◆ get()

type(structvectortype) function, pointer structarraymodule::get ( class(structarraytype this,
integer(i4b), intent(in)  idx 
)
private
Parameters
thisStructArrayType

Definition at line 197 of file StructArray.f90.

198  class(StructArrayType) :: this !< StructArrayType
199  integer(I4B), intent(in) :: idx
200  type(StructVectorType), pointer :: sv
201  call set_pointer(sv, this%struct_vectors(idx))
Here is the call graph for this function:

◆ load_deferred_vector()

subroutine structarraymodule::load_deferred_vector ( class(structarraytype this,
integer(i4b), intent(in)  icol 
)
Parameters
thisStructArrayType

Definition at line 410 of file StructArray.f90.

411  use memorymanagermodule, only: get_isize
412  class(StructArrayType) :: this !< StructArrayType
413  integer(I4B), intent(in) :: icol
414  integer(I4B) :: i, j, isize
415  integer(I4B), dimension(:), pointer, contiguous :: p_int1d
416  integer(I4B), dimension(:, :), pointer, contiguous :: p_int2d
417  real(DP), dimension(:), pointer, contiguous :: p_dbl1d
418  type(CharacterStringType), dimension(:), pointer, contiguous :: p_charstr1d
419  character(len=LENVARNAME) :: varname
420  logical(LGP) :: overwrite
421 
422  overwrite = .true.
423  if (this%struct_vectors(icol)%idt%blockname == 'SOLUTIONGROUP') &
424  overwrite = .false.
425 
426  ! set varname
427  varname = this%struct_vectors(icol)%idt%mf6varname
428  ! check if already mem managed variable
429  call get_isize(varname, this%mempath, isize)
430 
431  ! allocate and load based on memtype
432  select case (this%struct_vectors(icol)%memtype)
433  case (mtype_int)
434  if (isize > -1) then
435  ! variable exists, reallocate and append
436  call mem_setptr(p_int1d, varname, this%mempath)
437 
438  if (overwrite) then
439  ! overwrite existing array
440  if (this%nrow > isize) then
441  ! reallocate
442  call mem_reallocate(p_int1d, this%nrow, varname, this%mempath)
443  end if
444 
445  ! write new data
446  do i = 1, this%nrow
447  p_int1d(i) = this%struct_vectors(icol)%int1d(i)
448  end do
449 
450  if (isize > this%nrow) then
451  ! initialize excess space
452  do i = this%nrow + 1, isize
453  p_int1d(i) = izero
454  end do
455  end if
456  else
457  ! reallocate to new size
458  call mem_reallocate(p_int1d, this%nrow + isize, varname, this%mempath)
459 
460  ! write new data after existing
461  do i = 1, this%nrow
462  p_int1d(isize + i) = this%struct_vectors(icol)%int1d(i)
463  end do
464  end if
465  else
466  ! allocate memory manager vector
467  call mem_allocate(p_int1d, this%nrow, varname, this%mempath)
468 
469  ! load local vector to managed memory
470  do i = 1, this%nrow
471  p_int1d(i) = this%struct_vectors(icol)%int1d(i)
472  end do
473  end if
474 
475  ! deallocate local memory
476  deallocate (this%struct_vectors(icol)%int1d)
477 
478  ! update structvector
479  this%struct_vectors(icol)%int1d => p_int1d
480  this%struct_vectors(icol)%size = this%nrow
481  case (mtype_dbl)
482  if (isize > -1) then
483  call mem_setptr(p_dbl1d, varname, this%mempath)
484 
485  if (overwrite) then
486  if (this%nrow > isize) then
487  call mem_reallocate(p_dbl1d, this%nrow, varname, this%mempath)
488  end if
489 
490  do i = 1, this%nrow
491  p_dbl1d(i) = this%struct_vectors(icol)%dbl1d(i)
492  end do
493 
494  if (isize > this%nrow) then
495  do i = this%nrow + 1, isize
496  p_dbl1d(i) = dzero
497  end do
498  end if
499  else
500  call mem_reallocate(p_dbl1d, this%nrow + isize, varname, &
501  this%mempath)
502  do i = 1, this%nrow
503  p_dbl1d(isize + i) = this%struct_vectors(icol)%dbl1d(i)
504  end do
505  end if
506  else
507  call mem_allocate(p_dbl1d, this%nrow, varname, this%mempath)
508 
509  do i = 1, this%nrow
510  p_dbl1d(i) = this%struct_vectors(icol)%dbl1d(i)
511  end do
512  end if
513 
514  deallocate (this%struct_vectors(icol)%dbl1d)
515 
516  this%struct_vectors(icol)%dbl1d => p_dbl1d
517  this%struct_vectors(icol)%size = this%nrow
518  !
519  case (mtype_str)
520  if (isize > -1) then
521  call mem_setptr(p_charstr1d, varname, this%mempath)
522 
523  if (overwrite) then
524  if (this%nrow > isize) then
525  call mem_reallocate(p_charstr1d, linelength, this%nrow, varname, &
526  this%mempath)
527  end if
528 
529  do i = 1, this%nrow
530  p_charstr1d(i) = this%struct_vectors(icol)%charstr1d(i)
531  end do
532 
533  if (isize > this%nrow) then
534  do i = this%nrow + 1, isize
535  p_charstr1d(i) = ''
536  end do
537  end if
538  else
539  call mem_reallocate(p_charstr1d, linelength, this%nrow + isize, &
540  varname, this%mempath)
541  do i = 1, this%nrow
542  p_charstr1d(isize + i) = this%struct_vectors(icol)%charstr1d(i)
543  end do
544  end if
545  else
546  call mem_allocate(p_charstr1d, linelength, this%nrow, varname, &
547  this%mempath)
548  do i = 1, this%nrow
549  p_charstr1d(i) = this%struct_vectors(icol)%charstr1d(i)
550  call this%struct_vectors(icol)%charstr1d(i)%destroy()
551  end do
552  end if
553 
554  deallocate (this%struct_vectors(icol)%charstr1d)
555 
556  this%struct_vectors(icol)%charstr1d => p_charstr1d
557  this%struct_vectors(icol)%size = this%nrow
558  case (mtype_intvec) ! intvector reallocate unimplemented
559  errmsg = 'StructArray::load_deferred_vector &
560  &intvector reallocate unimplemented.'
561  call store_error(errmsg, terminate=.true.)
562  case (mtype_int2d)
563  if (isize > -1) then
564  errmsg = 'StructArray::load_deferred_vector &
565  &int2d reallocate unimplemented.'
566  call store_error(errmsg, terminate=.true.)
567  else
568  call mem_allocate(p_int2d, this%struct_vectors(icol)%intshape, &
569  this%nrow, varname, this%mempath)
570  do i = 1, this%nrow
571  do j = 1, this%struct_vectors(icol)%intshape
572  p_int2d(j, i) = this%struct_vectors(icol)%int2d(j, i)
573  end do
574  end do
575  end if
576 
577  deallocate (this%struct_vectors(icol)%int2d)
578 
579  this%struct_vectors(icol)%int2d => p_int2d
580  this%struct_vectors(icol)%size = this%nrow
581  case (mtype_dbl2d)
582  errmsg = 'StructArray::load_deferred_vector &
583  &dbl2d reallocate unimplemented.'
584  call store_error(errmsg, terminate=.true.)
585  case default
586  end select
Here is the call graph for this function:

◆ log_structarray_vars()

subroutine structarraymodule::log_structarray_vars ( class(structarraytype this,
integer(i4b), intent(in)  iout 
)
private
Parameters
thisStructArrayType
[in]ioutunit number for output

Definition at line 629 of file StructArray.f90.

630  class(StructArrayType) :: this !< StructArrayType
631  integer(I4B), intent(in) :: iout !< unit number for output
632  integer(I4B) :: j, nts
633  integer(I4B), dimension(:), pointer, contiguous :: int1d
634  character(len=LINELENGTH) :: ts_count_str
635 
636  ! idm variable logging
637  do j = 1, this%ncol
638  ! log based on memtype
639  select case (this%struct_vectors(j)%memtype)
640  case (mtype_int)
641  call idm_log_var(this%struct_vectors(j)%int1d, &
642  this%struct_vectors(j)%idt%tagname, &
643  this%mempath, iout)
644  case (mtype_dbl)
645  nts = this%struct_vectors(j)%ts_strlocs%count()
646  if (nts > 0) then
647  write (ts_count_str, '(i0, " time-series bound entries")') nts
648  call idm_log_var(this%struct_vectors(j)%idt%tagname, &
649  this%mempath, iout, .false., trim(ts_count_str))
650  else
651  call idm_log_var(this%struct_vectors(j)%dbl1d, &
652  this%struct_vectors(j)%idt%tagname, &
653  this%mempath, iout)
654  end if
655  case (mtype_intvec)
656  call mem_setptr(int1d, this%struct_vectors(j)%idt%mf6varname, &
657  this%mempath)
658  call idm_log_var(int1d, this%struct_vectors(j)%idt%tagname, &
659  this%mempath, iout)
660  case (mtype_int2d)
661  call idm_log_var(this%struct_vectors(j)%int2d, &
662  this%struct_vectors(j)%idt%tagname, &
663  this%mempath, iout)
664  case (mtype_dbl2d)
665  nts = this%struct_vectors(j)%ts_strlocs%count()
666  if (nts > 0) then
667  write (ts_count_str, '(i0, " time-series bound entries")') nts
668  call idm_log_var(this%struct_vectors(j)%idt%tagname, &
669  this%mempath, iout, .false., trim(ts_count_str))
670  else
671  call idm_log_var(this%struct_vectors(j)%dbl2d, &
672  this%struct_vectors(j)%idt%tagname, &
673  this%mempath, iout)
674  end if
675  end select
676  end do

◆ mem_create_vector()

subroutine structarraymodule::mem_create_vector ( class(structarraytype this,
integer(i4b), intent(in)  icol,
type(inputparamdefinitiontype), pointer  idt 
)
private
Parameters
thisStructArrayType
[in]icolcolumn to create

Definition at line 134 of file StructArray.f90.

135  class(StructArrayType) :: this !< StructArrayType
136  integer(I4B), intent(in) :: icol !< column to create
137  type(InputParamDefinitionType), pointer :: idt
138  type(StructVectorType) :: sv
139  integer(I4B) :: numcol
140 
141  ! initialize
142  numcol = 1
143  sv%idt => idt
144  sv%icol = icol
145 
146  ! set size
147  if (this%deferred_shape) then
148  sv%size = this%deferred_size_init
149  else
150  sv%size = this%nrow
151  end if
152 
153  ! allocate array memory for StructVectorType
154  select case (idt%datatype)
155  case ('INTEGER')
156  call this%allocate_int_type(sv)
157  case ('DOUBLE')
158  call this%allocate_dbl_type(sv)
159  case ('STRING', 'KEYWORD')
160  call this%allocate_charstr_type(sv)
161  case ('INTEGER1D')
162  call this%allocate_int1d_type(sv)
163  if (sv%memtype == mtype_int2d) then
164  numcol = sv%intshape
165  end if
166  case ('DOUBLE1D')
167  call this%allocate_dbl1d_type(sv)
168  numcol = sv%intshape
169  case default
170  errmsg = 'IDM unimplemented. StructArray::mem_create_vector &
171  &type='//trim(idt%datatype)
172  call store_error(errmsg, .true.)
173  end select
174 
175  ! set the object in the Struct Array
176  this%struct_vectors(icol) = sv
177  this%numcols(icol) = numcol
178  if (icol == 1) then
179  this%startidx(icol) = 1
180  else
181  this%startidx(icol) = this%startidx(icol - 1) + this%numcols(icol - 1)
182  end if
Here is the call graph for this function:

◆ memload_vectors()

subroutine structarraymodule::memload_vectors ( class(structarraytype this)
Parameters
thisStructArrayType

Definition at line 591 of file StructArray.f90.

592  class(StructArrayType) :: this !< StructArrayType
593  integer(I4B) :: icol, j
594  integer(I4B), dimension(:), pointer, contiguous :: p_intvector
595  character(len=LENVARNAME) :: varname
596 
597  do icol = 1, this%ncol
598  ! set varname
599  varname = this%struct_vectors(icol)%idt%mf6varname
600 
601  if (this%struct_vectors(icol)%memtype == mtype_intvec) then
602  ! intvectors always need to be loaded
603  ! size intvector to number of values read
604  call this%struct_vectors(icol)%intvector%shrink_to_fit()
605 
606  ! allocate memory manager vector
607  call mem_allocate(p_intvector, &
608  this%struct_vectors(icol)%intvector%size, &
609  varname, this%mempath)
610 
611  ! load local vector to managed memory
612  do j = 1, this%struct_vectors(icol)%intvector%size
613  p_intvector(j) = this%struct_vectors(icol)%intvector%at(j)
614  end do
615 
616  ! cleanup local memory
617  call this%struct_vectors(icol)%intvector%destroy()
618  deallocate (this%struct_vectors(icol)%intvector)
619  nullify (this%struct_vectors(icol)%intvector_shape)
620  else if (this%deferred_shape) then
621  ! load as shape wasn't known
622  call this%load_deferred_vector(icol)
623  end if
624  end do

◆ read_from_binary()

integer(i4b) function structarraymodule::read_from_binary ( class(structarraytype this,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)
Parameters
thisStructArrayType
[in]inunitunit number for binary input
[in]ioutunit number for output

Definition at line 1046 of file StructArray.f90.

1047  class(StructArrayType) :: this !< StructArrayType
1048  integer(I4B), intent(in) :: inunit !< unit number for binary input
1049  integer(I4B), intent(in) :: iout !< unit number for output
1050  integer(I4B) :: irow, ierr
1051  integer(I4B) :: j, k
1052  integer(I4B) :: intval, numval
1053  character(len=LINELENGTH) :: fname
1054  character(len=*), parameter :: fmtlsterronly = &
1055  "('Error reading LIST from file: ',&
1056  &1x,a,1x,' on UNIT: ',I0)"
1057 
1058  ! set error and exit if deferred shape
1059  if (this%deferred_shape) then
1060  errmsg = 'IDM unimplemented. StructArray::read_from_binary deferred shape &
1061  &not supported for binary inputs.'
1062  call store_error(errmsg, terminate=.true.)
1063  end if
1064  ! initialize
1065  irow = 0
1066  ierr = 0
1067  readloop: do
1068  ! update irow index
1069  irow = irow + 1
1070  ! handle line reads by column memtype
1071  do j = 1, this%ncol
1072  select case (this%struct_vectors(j)%memtype)
1073  case (mtype_int)
1074  read (inunit, iostat=ierr) this%struct_vectors(j)%int1d(irow)
1075  case (mtype_dbl)
1076  read (inunit, iostat=ierr) this%struct_vectors(j)%dbl1d(irow)
1077  case (mtype_str)
1078  errmsg = 'List style binary inputs not supported &
1079  &for text columns, tag='// &
1080  trim(this%struct_vectors(j)%idt%tagname)//'.'
1081  call store_error(errmsg, terminate=.true.)
1082  case (mtype_intvec)
1083  ! get shape for this row
1084  numval = this%struct_vectors(j)%intvector_shape(irow)
1085  ! read and store row values
1086  do k = 1, numval
1087  if (ierr == 0) then
1088  read (inunit, iostat=ierr) intval
1089  call this%struct_vectors(j)%intvector%push_back(intval)
1090  end if
1091  end do
1092  case (mtype_int2d)
1093  ! read and store row values
1094  do k = 1, this%struct_vectors(j)%intshape
1095  if (ierr == 0) then
1096  read (inunit, iostat=ierr) this%struct_vectors(j)%int2d(k, irow)
1097  end if
1098  end do
1099  case (mtype_dbl2d)
1100  do k = 1, this%struct_vectors(j)%intshape
1101  if (ierr == 0) then
1102  read (inunit, iostat=ierr) this%struct_vectors(j)%dbl2d(k, irow)
1103  end if
1104  end do
1105  end select
1106 
1107  ! handle error cases
1108  select case (ierr)
1109  case (0)
1110  ! no error
1111  case (:-1)
1112  ! End of block was encountered
1113  irow = irow - 1
1114  exit readloop
1115  case (1:)
1116  ! Error
1117  inquire (unit=inunit, name=fname)
1118  write (errmsg, fmtlsterronly) trim(adjustl(fname)), inunit
1119  call store_error(errmsg, terminate=.true.)
1120  case default
1121  end select
1122  end do
1123  if (irow == this%nrow) exit readloop
1124  end do readloop
1125 
1126  ! Stop if errors were detected
1127  !if (count_errors() > 0) then
1128  ! call store_error_unit(inunit)
1129  !end if
1130 
1131  ! if deferred shape vectors were read, load to input path
1132  call this%memload_vectors()
1133 
1134  ! log loaded variables
1135  if (iout > 0) then
1136  call this%log_structarray_vars(iout)
1137  end if
Here is the call graph for this function:

◆ read_from_parser()

integer(i4b) function structarraymodule::read_from_parser ( class(structarraytype this,
type(blockparsertype parser,
logical(lgp), intent(in)  timeseries,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  input_name 
)
private
Parameters
thisStructArrayType
parserblock parser to read from
[in]ioutunit number for output
[in]input_nameinput filename for error messages

Definition at line 852 of file StructArray.f90.

854  class(StructArrayType) :: this !< StructArrayType
855  type(BlockParserType) :: parser !< block parser to read from
856  logical(LGP), intent(in) :: timeseries
857  integer(I4B), intent(in) :: iout !< unit number for output
858  character(len=*), intent(in) :: input_name !< input filename for error messages
859  integer(I4B) :: irow, j
860  logical(LGP) :: endOfBlock
861 
862  ! initialize index irow
863  irow = 0
864 
865  ! reset nrow if deferred shape
866  if (this%deferred_shape) then
867  this%nrow = 0
868  end if
869 
870  ! read entire block
871  do
872  ! read next line
873  call parser%GetNextLine(endofblock)
874  if (endofblock) then
875  ! no more lines
876  exit
877  else if (this%deferred_shape) then
878  ! shape unknown, track lines read
879  this%nrow = this%nrow + 1
880  ! check and update memory allocation
881  call this%check_reallocate()
882  end if
883  ! update irow index
884  irow = irow + 1
885  if (this%deferred_shape) then
886  else
887  ! check allocated array size against user bound
888  if (irow > this%nrow) then
889  write (errmsg, '(a,i0,a)') &
890  'Input error: line count exceeds input dimension. Expected rows=', &
891  this%nrow, '.'
892  call store_error(errmsg)
893  call store_error_filename(input_name)
894  end if
895  end if
896  ! handle line reads by column memtype
897  do j = 1, this%ncol
898  call this%read_param(parser, j, irow, timeseries, iout)
899  end do
900  end do
901  ! if deferred shape vectors were read, load to input path
902  call this%memload_vectors()
903  ! log loaded variables
904  if (iout > 0) then
905  call this%log_structarray_vars(iout)
906  end if
Here is the call graph for this function:

◆ read_from_parser_keystring()

integer(i4b) function structarraymodule::read_from_parser_keystring ( class(structarraytype this,
type(blockparsertype parser,
logical(lgp), intent(in)  timeseries,
integer(i4b), intent(in)  nleading,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  input_name 
)
private

Each input line contains nleading fixed columns followed by a dispatch keyword. Two dispatch modes are supported:

Simple dispatch: keyword matches a DOUBLE/STRING/INTEGER column. One value token is read from the parser into that column. All other member columns receive their sentinel for that row.

Compound dispatch: keyword matches a KEYWORD-type column (e.g. FLOWING_WELL). No parser token is read for the KEYWORD column — the dispatch keyword itself is stored. Subsequent non-KEYWORD columns (the compound sub-members, e.g. FWELEV/FWCOND/FWRLEN) are read from the parser in order until the next KEYWORD column or the end of the member columns.

No-value KEYWORD dispatch: a KEYWORD column with no sub-members immediately following stores the dispatch keyword and reads nothing further.

Parameters
thisStructArrayType
parserblock parser to read from
[in]timeseries.true. when TS files loaded
[in]nleadingnumber of leading (fixed) columns
[in]ioutunit number for output
[in]input_nameinput filename for error messages

Definition at line 930 of file StructArray.f90.

932  use inputoutputmodule, only: upcase
934  class(StructArrayType) :: this !< StructArrayType
935  type(BlockParserType) :: parser !< block parser to read from
936  logical(LGP), intent(in) :: timeseries !< .true. when TS files loaded
937  integer(I4B), intent(in) :: nleading !< number of leading (fixed) columns
938  integer(I4B), intent(in) :: iout !< unit number for output
939  character(len=*), intent(in) :: input_name !< input filename for error messages
940  integer(I4B) :: irow
941  logical(LGP) :: endOfBlock, is_keyword_dispatch
942  character(len=LINELENGTH) :: keyword
943  integer(I4B) :: icol, found_col, last_set_col
944 
945  irow = 0
946 
947  ! reset nrow if deferred shape
948  if (this%deferred_shape) then
949  this%nrow = 0
950  end if
951 
952  do
953  call parser%GetNextLine(endofblock)
954  if (endofblock) exit
955 
956  if (this%deferred_shape) then
957  this%nrow = this%nrow + 1
958  call this%check_reallocate()
959  end if
960  irow = irow + 1
961 
962  ! bounds check for pre-allocated (non-deferred) arrays
963  if (.not. this%deferred_shape) then
964  if (irow > this%nrow) then
965  write (errmsg, '(a,i0,a)') &
966  'Input error: keystring row count exceeds pre-allocated maxbound=', &
967  this%nrow, '.'
968  call store_error(errmsg)
969  call store_error_filename(input_name)
970  exit
971  end if
972  end if
973 
974  ! read leading fixed columns (e.g. CELLID, IFNO)
975  do icol = 1, nleading
976  call this%read_param(parser, icol, irow, .false., iout)
977  end do
978 
979  ! read dispatch keyword
980  call parser%GetString(keyword)
981  call upcase(keyword)
982 
983  ! find the matching keystring-member column
984  found_col = 0
985  do icol = nleading + 1, this%ncol
986  if (trim(this%struct_vectors(icol)%idt%tagname) == trim(keyword)) then
987  found_col = icol
988  exit
989  end if
990  end do
991 
992  if (found_col < 1) then
993  write (errmsg, '(a,a,a)') &
994  'Unrecognized keystring keyword "', trim(keyword), &
995  '" in PERIOD block.'
996  call store_error(errmsg)
997  call store_error_filename(input_name)
998  cycle
999  end if
1000 
1001  ! determine dispatch mode and set/read matched column(s)
1002  is_keyword_dispatch = &
1003  (this%struct_vectors(found_col)%idt%datatype == 'KEYWORD')
1004 
1005  if (is_keyword_dispatch) then
1006  ! Compound or no-value KEYWORD dispatch:
1007  ! the dispatch keyword token was already consumed; store it in the
1008  ! KEYWORD column (charstr) without reading another token.
1009  this%struct_vectors(found_col)%charstr1d(irow) = trim(keyword)
1010  last_set_col = found_col
1011  ! read exactly nsubmembers compound sub-member columns that
1012  ! immediately follow in the struct array
1013  do icol = found_col + 1, &
1014  found_col + this%struct_vectors(found_col)%nsubmembers
1015  if (icol > this%ncol) exit
1016  call this%read_param(parser, icol, irow, timeseries, iout)
1017  last_set_col = icol
1018  end do
1019  else
1020  ! Simple single-value dispatch: read one value for the matched column
1021  call this%read_param(parser, found_col, irow, timeseries, iout)
1022  last_set_col = found_col
1023  end if
1024 
1025  ! fill sentinels for all non-matched member columns
1026  do icol = nleading + 1, this%ncol
1027  if (icol >= found_col .and. icol <= last_set_col) cycle
1028  select case (this%struct_vectors(icol)%memtype)
1029  case (mtype_dbl) ! DOUBLE: use DNODATA sentinel
1030  this%struct_vectors(icol)%dbl1d(irow) = dnodata
1031  case (mtype_str) ! STRING or KEYWORD: use empty string sentinel
1032  this%struct_vectors(icol)%charstr1d(irow) = ''
1033  end select
1034  end do
1035  end do
1036 
1037  call this%memload_vectors()
1038 
1039  if (iout > 0) then
1040  call this%log_structarray_vars(iout)
1041  end if
subroutine, public upcase(word)
Convert to upper case.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
Here is the call graph for this function:

◆ read_param()

subroutine structarraymodule::read_param ( class(structarraytype this,
type(blockparsertype), intent(inout)  parser,
integer(i4b), intent(in)  sv_col,
integer(i4b), intent(in)  irow,
logical(lgp), intent(in)  timeseries,
integer(i4b), intent(in)  iout,
integer(i4b), intent(in), optional  auxcol 
)
private
Parameters
thisStructArrayType
[in,out]parserblock parser to read from
[in]ioutunit number for output

Definition at line 771 of file StructArray.f90.

772  class(StructArrayType) :: this !< StructArrayType
773  type(BlockParserType), intent(inout) :: parser !< block parser to read from
774  integer(I4B), intent(in) :: sv_col
775  integer(I4B), intent(in) :: irow
776  logical(LGP), intent(in) :: timeseries
777  integer(I4B), intent(in) :: iout !< unit number for output
778  integer(I4B), optional, intent(in) :: auxcol
779  integer(I4B) :: n, intval, numval, icol
780  character(len=LINELENGTH) :: str
781  character(len=:), allocatable :: line
782  logical(LGP) :: preserve_case
783 
784  select case (this%struct_vectors(sv_col)%memtype)
785  case (mtype_int)
786  ! if reloadable block and first col, store blocknum
787  if (sv_col == 1 .and. this%blocknum > 0) then
788  ! store blocknum
789  this%struct_vectors(sv_col)%int1d(irow) = this%blocknum
790  else
791  ! read and store int
792  this%struct_vectors(sv_col)%int1d(irow) = parser%GetInteger()
793  end if
794  case (mtype_dbl)
795  if (this%struct_vectors(sv_col)%idt%timeseries .and. timeseries) then
796  call parser%GetString(str)
797  if (present(auxcol)) then
798  icol = auxcol
799  else
800  icol = 1
801  end if
802  this%struct_vectors(sv_col)%dbl1d(irow) = &
803  this%struct_vectors(sv_col)%read_token(str, this%startidx(sv_col), &
804  icol, irow)
805  else
806  this%struct_vectors(sv_col)%dbl1d(irow) = parser%GetDouble()
807  end if
808  case (mtype_str)
809  if (this%struct_vectors(sv_col)%idt%shape /= '') then
810  ! if last column with any shape, store rest of line
811  if (sv_col == this%ncol) then
812  call parser%GetRemainingLine(line)
813  this%struct_vectors(sv_col)%charstr1d(irow) = line
814  deallocate (line)
815  end if
816  else
817  ! read string token
818  preserve_case = (.not. this%struct_vectors(sv_col)%idt%preserve_case)
819  call parser%GetString(str, preserve_case)
820  this%struct_vectors(sv_col)%charstr1d(irow) = str
821  end if
822  case (mtype_intvec)
823  ! get shape for this row
824  numval = this%struct_vectors(sv_col)%intvector_shape(irow)
825  ! read and store row values
826  do n = 1, numval
827  intval = parser%GetInteger()
828  call this%struct_vectors(sv_col)%intvector%push_back(intval)
829  end do
830  case (mtype_int2d)
831  ! read and store row values
832  do n = 1, this%struct_vectors(sv_col)%intshape
833  this%struct_vectors(sv_col)%int2d(n, irow) = parser%GetInteger()
834  end do
835  case (mtype_dbl2d)
836  ! read and store row values
837  do n = 1, this%struct_vectors(sv_col)%intshape
838  if (this%struct_vectors(sv_col)%idt%timeseries .and. timeseries) then
839  call parser%GetString(str)
840  icol = this%startidx(sv_col) + n - 1
841  this%struct_vectors(sv_col)%dbl2d(n, irow) = &
842  this%struct_vectors(sv_col)%read_token(str, icol, n, irow)
843  else
844  this%struct_vectors(sv_col)%dbl2d(n, irow) = parser%GetDouble()
845  end if
846  end do
847  end select

◆ set_pointer()

subroutine structarraymodule::set_pointer ( type(structvectortype), pointer  sv,
type(structvectortype), target  sv_target 
)
private

Definition at line 191 of file StructArray.f90.

192  type(StructVectorType), pointer :: sv
193  type(StructVectorType), target :: sv_target
194  sv => sv_target
Here is the caller graph for this function:

◆ ts_update()

subroutine structarraymodule::ts_update ( class(structarraytype), intent(inout)  this,
type(timeseriesmanagertype), intent(inout), pointer  tsmanager,
character(len=*), intent(in)  subcomp_name,
integer(i4b), intent(in)  iprpak,
character(len=*), intent(in)  input_name,
type(characterstringtype), dimension(:), intent(in), optional, pointer  auxname_cst,
logical(lgp), intent(in), optional  clear_strlocs 
)
private

Iterates over struct vectors that carry deferred TS tokens and registers each with the supplied tsmanager. Handles both BND (MTYPE_DBL / dbl1d) and AUX (MTYPE_DBL2D / dbl2d) columns. Pass auxname_cst only when AUX columns may carry time series.

Parameters
[in]clear_strlocsif .false. strlocs are preserved for re-registration (default .true.)

Definition at line 1148 of file StructArray.f90.

1150  class(StructArrayType), intent(inout) :: this
1151  type(TimeSeriesManagerType), pointer, intent(inout) :: tsmanager
1152  character(len=*), intent(in) :: subcomp_name
1153  integer(I4B), intent(in) :: iprpak
1154  character(len=*), intent(in) :: input_name
1155  type(CharacterStringType), dimension(:), pointer, intent(in), &
1156  optional :: auxname_cst
1157  logical(LGP), optional, intent(in) :: clear_strlocs !< if .false. strlocs are preserved for re-registration (default .true.)
1158  type(TSStringLocType), pointer :: ts_strloc
1159  type(TimeSeriesLinkType), pointer :: tsLink
1160  real(DP), pointer :: bndElem
1161  character(len=LENBOUNDNAME) :: boundname
1162  integer(I4B) :: m, n, iboundname
1163  logical(LGP) :: do_clear
1164 
1165  do_clear = .true.
1166  if (present(clear_strlocs)) do_clear = clear_strlocs
1167 
1168  ! find BOUNDNAME column (0 = none)
1169  iboundname = 0
1170  do m = 1, this%ncol
1171  if (this%struct_vectors(m)%idt%mf6varname == 'BOUNDNAME') then
1172  iboundname = m
1173  exit
1174  end if
1175  end do
1176 
1177  do m = 1, this%ncol
1178  if (.not. this%struct_vectors(m)%idt%timeseries) cycle
1179  do n = 1, this%struct_vectors(m)%ts_strlocs%count()
1180  ts_strloc => this%struct_vectors(m)%get_ts_strloc(n)
1181  nullify (tslink)
1182  select case (this%struct_vectors(m)%memtype)
1183  case (mtype_dbl) ! dbl1d (BND)
1184  bndelem => this%struct_vectors(m)%dbl1d(ts_strloc%row)
1185  call read_value_or_time_series(ts_strloc%token, ts_strloc%row, &
1186  ts_strloc%structarray_col, bndelem, &
1187  subcomp_name, 'BND', tsmanager, &
1188  iprpak, tslink)
1189  if (associated(tslink)) then
1190  tslink%Text = this%struct_vectors(m)%idt%mf6varname
1191  if (iboundname > 0) then
1192  boundname = &
1193  this%struct_vectors(iboundname)%charstr1d(ts_strloc%row)
1194  tslink%BndName = boundname
1195  end if
1196  end if
1197  case (mtype_dbl2d) ! dbl2d (AUX)
1198  if (.not. present(auxname_cst)) cycle
1199  if (.not. associated(auxname_cst)) cycle
1200  bndelem => this%struct_vectors(m)%dbl2d(ts_strloc%col, ts_strloc%row)
1201  call read_value_or_time_series(ts_strloc%token, ts_strloc%row, &
1202  ts_strloc%structarray_col, bndelem, &
1203  subcomp_name, 'AUX', tsmanager, &
1204  iprpak, tslink)
1205  if (associated(tslink)) then
1206  tslink%Text = auxname_cst(ts_strloc%col)
1207  if (iboundname > 0) then
1208  boundname = &
1209  this%struct_vectors(iboundname)%charstr1d(ts_strloc%row)
1210  tslink%BndName = boundname
1211  end if
1212  end if
1213  end select
1214  end do
1215  if (do_clear) call this%struct_vectors(m)%clear()
1216  end do
1217 
1218  if (count_errors() > 0) then
1219  call store_error_filename(input_name)
1220  end if
Here is the call graph for this function: