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

This module contains the LayerArrayLoadModule. More...

Data Types

type  layerarrayloadtype
 Ascii array layer dynamic loader type. More...
 

Functions/Subroutines

subroutine ainit (this, mf6_input, component_name, component_input_name, input_name, iperblock, parser, iout)
 
subroutine df (this)
 
subroutine ad (this)
 
subroutine rp (this, parser)
 
subroutine destroy (this)
 
subroutine reset (this)
 
subroutine init_charstr1d (this, varname, input_name)
 
subroutine params_alloc (this)
 
subroutine param_load (this, parser, idt, mempath, netcdf, iaux)
 
subroutine tas_arrays_alloc (this)
 
subroutine tas_links_create (this, inunit)
 

Detailed Description

This module contains the routines for reading period block array based input that is associated with a layer and an layer index array, such as with the EVTA and RCHA packages.

Function/Subroutine Documentation

◆ ad()

subroutine layerarrayloadmodule::ad ( class(layerarrayloadtype), intent(inout)  this)
private

Definition at line 121 of file Mf6FileLayerArray.f90.

122  class(LayerArrayLoadType), intent(inout) :: this
123  call this%tasmanager%ad()

◆ ainit()

subroutine layerarrayloadmodule::ainit ( class(layerarrayloadtype), intent(inout)  this,
type(modflowinputtype), intent(in)  mf6_input,
character(len=*), intent(in)  component_name,
character(len=*), intent(in)  component_input_name,
character(len=*), intent(in)  input_name,
integer(i4b), intent(in)  iperblock,
type(blockparsertype), intent(inout), pointer  parser,
integer(i4b), intent(in)  iout 
)
private

Definition at line 56 of file Mf6FileLayerArray.f90.

62  class(LayerArrayLoadType), intent(inout) :: this
63  type(ModflowInputType), intent(in) :: mf6_input
64  character(len=*), intent(in) :: component_name
65  character(len=*), intent(in) :: component_input_name
66  character(len=*), intent(in) :: input_name
67  integer(I4B), intent(in) :: iperblock
68  type(BlockParserType), pointer, intent(inout) :: parser
69  integer(I4B), intent(in) :: iout
70  type(LoadMf6FileType) :: loader
71  type(CharacterStringType), dimension(:), pointer, &
72  contiguous :: tas_fnames
73  character(len=LINELENGTH) :: fname
74  integer(I4B) :: tas6_size, n
75 
76  ! initialize base type
77  call this%DynamicPkgLoadType%init(mf6_input, component_name, &
78  component_input_name, &
79  input_name, iperblock, iout)
80  ! initialize
81  nullify (this%aux_tasnames)
82  nullify (this%param_tasnames)
83  this%tas_active = 0
84  this%iout = iout
85 
86  ! load static input
87  call loader%load(parser, mf6_input, this%nc_vars, this%input_name, iout)
88 
89  ! create tasmanager
90  allocate (this%tasmanager)
91  call tasmanager_cr(this%tasmanager, modelname=this%mf6_input%component_name, &
92  iout=this%iout)
93 
94  ! determine if TAS6 files were provided in OPTIONS block
95  call get_isize('TAS6_FILENAME', this%mf6_input%mempath, tas6_size)
96  if (tas6_size > 0) then
97  this%tas_active = 1
98  call mem_setptr(tas_fnames, 'TAS6_FILENAME', this%mf6_input%mempath)
99  ! add files to tasmanager
100  do n = 1, size(tas_fnames)
101  fname = tas_fnames(n)
102  call this%tasmanager%add_tasfile(fname)
103  end do
104  end if
105 
106  ! initialize input context memory
107  call this%ctx%init(mf6_input)
108 
109  ! allocate dfn params
110  call this%params_alloc()
111 
112  ! allocate memory for storing TAS strings
113  call this%tas_arrays_alloc()
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains the LoadMf6FileModule.
Definition: LoadMf6File.f90:8
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
Static parser based input loader.
Definition: LoadMf6File.f90:54
Here is the call graph for this function:

◆ destroy()

subroutine layerarrayloadmodule::destroy ( class(layerarrayloadtype), intent(inout)  this)

Definition at line 226 of file Mf6FileLayerArray.f90.

227  class(LayerArrayLoadType), intent(inout) :: this
228  !
229  ! nullify ctx pointers (including mshape) before deallocate
230  call this%ctx%destroy()
231  !
232  ! deallocate tasmanager
233  call this%tasmanager%da()
234  deallocate (this%tasmanager)
235  nullify (this%tasmanager)

◆ df()

subroutine layerarrayloadmodule::df ( class(layerarrayloadtype), intent(inout)  this)

Definition at line 116 of file Mf6FileLayerArray.f90.

117  class(LayerArrayLoadType), intent(inout) :: this
118  call this%tasmanager%tasmanager_df()

◆ init_charstr1d()

subroutine layerarrayloadmodule::init_charstr1d ( class(layerarrayloadtype this,
character(len=*), intent(in)  varname,
character(len=*), intent(in)  input_name 
)
private

Definition at line 263 of file Mf6FileLayerArray.f90.

265  class(LayerArrayLoadType) :: this
266  character(len=*), intent(in) :: varname
267  character(len=*), intent(in) :: input_name
268  type(CharacterStringType), dimension(:), pointer, &
269  contiguous :: charstr1d
270  integer(I4B) :: n
271  call mem_setptr(charstr1d, varname, this%mf6_input%mempath)
272  do n = 1, size(charstr1d)
273  charstr1d(n) = ''
274  end do

◆ param_load()

subroutine layerarrayloadmodule::param_load ( class(layerarrayloadtype), intent(inout)  this,
type(blockparsertype), intent(in)  parser,
type(inputparamdefinitiontype), intent(in)  idt,
character(len=*), intent(in)  mempath,
logical(lgp), intent(in)  netcdf,
integer(i4b), intent(in)  iaux 
)
private

Definition at line 301 of file Mf6FileLayerArray.f90.

302  use tdismodule, only: kper
304  use arrayhandlersmodule, only: ifind
311  use idmloggermodule, only: idm_log_var
312  class(LayerArrayLoadType), intent(inout) :: this
313  type(BlockParserType), intent(in) :: parser
314  type(InputParamDefinitionType), intent(in) :: idt
315  character(len=*), intent(in) :: mempath
316  logical(LGP), intent(in) :: netcdf
317  integer(I4B), intent(in) :: iaux
318  integer(I4B), dimension(:), pointer, contiguous :: int1d
319  real(DP), dimension(:), pointer, contiguous :: dbl1d
320  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
321  integer(I4B) :: iparam, n
322 
323  select case (idt%datatype)
324  case ('INTEGER1D')
325  call mem_setptr(int1d, idt%mf6varname, mempath)
326  if (netcdf) then
327  call netcdf_read_array(int1d, this%ctx%mshape, idt, &
328  this%mf6_input, this%nc_vars, this%input_name, &
329  this%iout, kper)
330  else
331  call read_int1d(parser, int1d, idt%mf6varname)
332  end if
333  call idm_log_var(int1d, idt%tagname, mempath, this%iout)
334  case ('DOUBLE1D')
335  call mem_setptr(dbl1d, idt%mf6varname, mempath)
336  if (netcdf) then
337  call netcdf_read_array(dbl1d, this%ctx%mshape, idt, &
338  this%mf6_input, this%nc_vars, this%input_name, &
339  this%iout, kper)
340  else
341  call read_dbl1d(parser, dbl1d, idt%mf6varname)
342  end if
343  call idm_log_var(dbl1d, idt%tagname, mempath, this%iout)
344  case ('DOUBLE2D')
345  call mem_setptr(dbl2d, idt%mf6varname, mempath)
346  allocate (dbl1d(this%ctx%ncpl))
347  if (netcdf) then
348  call netcdf_read_array(dbl1d, this%ctx%mshape, idt, &
349  this%mf6_input, this%nc_vars, this%input_name, &
350  this%iout, kper, iaux)
351  else
352  call read_dbl1d(parser, dbl1d, idt%mf6varname)
353  end if
354  do n = 1, this%ctx%ncpl
355  dbl2d(iaux, n) = dbl1d(n)
356  end do
357  call idm_log_var(dbl1d, idt%tagname, mempath, this%iout)
358  deallocate (dbl1d)
359  case default
360  errmsg = 'IDM unimplemented. LayerArrayLoad::param_load &
361  &datatype='//trim(idt%datatype)
362  call store_error(errmsg)
363  call store_error_filename(this%input_name)
364  end select
365 
366  ! if param is tracked set read state
367  iparam = ifind(this%param_names, idt%tagname)
368  if (iparam > 0) then
369  this%param_reads(iparam)%invar = 1
370  end if
This module contains the DefinitionSelectModule.
type(inputparamdefinitiontype) function, pointer, public get_param_definition_type(input_definition_types, component_type, subcomponent_type, blockname, tagname, filename, found)
Return parameter definition.
subroutine, public read_dbl1d(parser, dbl1d, aname)
subroutine, public read_dbl2d(parser, dbl2d, aname)
This module contains the Input Data Model Logger Module.
Definition: IdmLogger.f90:7
Input definition module.
subroutine, public read_int1d(parser, int1d, aname)
This module contains the LoadNCInputModule.
Definition: LoadNCInput.F90:7
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Input parameter definition. Describes an input parameter.
Here is the call graph for this function:

◆ params_alloc()

subroutine layerarrayloadmodule::params_alloc ( class(layerarrayloadtype), intent(inout)  this)

Definition at line 277 of file Mf6FileLayerArray.f90.

278  class(LayerArrayLoadType), intent(inout) :: this
279  character(len=LENVARNAME) :: rs_varname
280  integer(I4B), pointer :: intvar
281  integer(I4B) :: iparam
282 
283  ! set in scope param names
284  call this%ctx%tags(this%param_names, this%nparam, this%input_name, &
285  create=.true.)
286  call this%ctx%allocate_arrays()
287 
288  ! allocate and set param_reads pointer array
289  allocate (this%param_reads(this%nparam))
290 
291  ! store read state variable pointers
292  do iparam = 1, this%nparam
293  ! allocate and store name of read state variable
294  rs_varname = this%ctx%rsv_alloc(this%param_names(iparam))
295  call mem_setptr(intvar, rs_varname, this%mf6_input%mempath)
296  this%param_reads(iparam)%invar => intvar
297  this%param_reads(iparam)%invar = 0
298  end do

◆ reset()

subroutine layerarrayloadmodule::reset ( class(layerarrayloadtype), intent(inout)  this)
private

Definition at line 238 of file Mf6FileLayerArray.f90.

239  class(LayerArrayLoadType), intent(inout) :: this
240  integer(I4B) :: n, m
241 
242  if (this%tas_active /= 0) then
243  ! reset tasmanager
244  call this%tasmanager%reset(this%mf6_input%subcomponent_name)
245  ! reinitialize tas name arrays
246  call this%init_charstr1d('AUXTASNAME', this%input_name)
247  call this%init_charstr1d('PARAMTASNAME', this%input_name)
248  end if
249 
250  do n = 1, this%nparam
251  ! reset read state
252  this%param_reads(n)%invar = 0
253  end do
254 
255  ! explicitly reset auxvar array each period
256  do m = 1, this%ctx%ncpl
257  do n = 1, this%ctx%naux
258  this%ctx%auxvar(n, m) = dzero
259  end do
260  end do

◆ rp()

subroutine layerarrayloadmodule::rp ( class(layerarrayloadtype), intent(inout)  this,
type(blockparsertype), intent(inout), pointer  parser 
)
private

Definition at line 126 of file Mf6FileLayerArray.f90.

131  use arrayhandlersmodule, only: ifind
134  class(LayerArrayLoadType), intent(inout) :: this
135  type(BlockParserType), pointer, intent(inout) :: parser
136  logical(LGP) :: endOfBlock, netcdf
137  character(len=LINELENGTH) :: keyword, param_tag
138  type(InputParamDefinitionType), pointer :: idt
139  integer(I4B) :: iaux, iparam
140  character(len=LENTIMESERIESNAME) :: tas_name
141  integer(I4B), dimension(:), pointer, contiguous :: int1d
142 
143  ! reset for this period
144  call this%reset()
145 
146  ! log lst file header
147  call idm_log_header(this%mf6_input%component_name, &
148  this%mf6_input%subcomponent_name, this%iout)
149 
150  ! read array block
151  do
152  ! initialize
153  iaux = 0
154  netcdf = .false.
155 
156  ! read next line
157  call parser%GetNextLine(endofblock)
158  if (endofblock) exit
159  ! read param_tag
160  call parser%GetStringCaps(param_tag)
161 
162  ! is param tag an auxvar?
163  iaux = ifind_charstr(this%ctx%auxname_cst, param_tag)
164  ! any auvxar corresponds to the definition tag 'AUX'
165  if (iaux > 0) param_tag = 'AUX'
166 
167  ! set input definition
168  idt => get_param_definition_type(this%mf6_input%param_dfns, &
169  this%mf6_input%component_type, &
170  this%mf6_input%subcomponent_type, &
171  'PERIOD', param_tag, this%input_name)
172  ! look for TAS and NetCDF keywords
173  call parser%GetStringCaps(keyword)
174  if (keyword == 'TIMEARRAYSERIES') then
175  if (this%tas_active /= 0) then
176  call parser%GetStringCaps(tas_name)
177  if (param_tag == 'AUX') then
178  this%aux_tasnames(iaux) = tas_name
179  else
180  iparam = ifind(this%param_names, param_tag)
181  this%param_tasnames(iparam) = tas_name
182  this%param_reads(iparam)%invar = 2
183  end if
184  ! log variable
185  call idm_log_var(param_tag, this%mf6_input%mempath, this%iout, .true., &
186  trim(tas_name))
187  ! cycle to next input param
188  cycle
189  else
190  ! TODO: throw error
191  end if
192  else if (keyword == 'NETCDF') then
193  netcdf = .true.
194  end if
195 
196  ! read and load the parameter
197  call this%param_load(parser, idt, this%mf6_input%mempath, netcdf, iaux)
198  end do
199 
200  ! check if layer index variable was read; default to 1 if not provided.
201  ! the layer index variable follows the I<TYPE3> naming convention (e.g.
202  ! IEVT, IRCH), consistent with the in_scope check in LoadContext.
203  if (this%param_reads(1)%invar == 0) then
204  if (this%param_names(1) == &
205  'I'//trim(this%mf6_input%subcomponent_type(1:3))) then
206  idt => get_param_definition_type(this%mf6_input%param_dfns, &
207  this%mf6_input%component_type, &
208  this%mf6_input%subcomponent_type, &
209  'PERIOD', this%param_names(1), &
210  this%input_name)
211  ! set to default of 1 without updating invar
212  call mem_setptr(int1d, idt%mf6varname, this%mf6_input%mempath)
213  int1d = 1
214  end if
215  end if
216 
217  if (this%tas_active /= 0) then
218  call this%tas_links_create(parser%iuactive)
219  end if
220 
221  ! log lst file header
222  call idm_log_close(this%mf6_input%component_name, &
223  this%mf6_input%subcomponent_name, this%iout)
subroutine, public idm_log_close(component, subcomponent, iout)
@ brief log the closing message
Definition: IdmLogger.f90:56
subroutine, public idm_log_header(component, subcomponent, iout)
@ brief log a header message
Definition: IdmLogger.f90:44
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
integer(i4b) function, public ifind_charstr(array, str)
Here is the call graph for this function:

◆ tas_arrays_alloc()

subroutine layerarrayloadmodule::tas_arrays_alloc ( class(layerarrayloadtype), intent(inout)  this)

Definition at line 373 of file Mf6FileLayerArray.f90.

375  class(LayerArrayLoadType), intent(inout) :: this
376 
377  ! count params other than AUX
378  if (this%tas_active /= 0) then
379  call mem_allocate(this%aux_tasnames, lentimeseriesname, &
380  this%ctx%naux, 'AUXTASNAME', &
381  this%mf6_input%mempath)
382  call mem_allocate(this%param_tasnames, lentimeseriesname, this%nparam, &
383  'PARAMTASNAME', this%mf6_input%mempath)
384  call this%init_charstr1d('AUXTASNAME', this%input_name)
385  call this%init_charstr1d('PARAMTASNAME', this%input_name)
386  else
387  call mem_allocate(this%aux_tasnames, lentimeseriesname, 0, &
388  'AUXTASNAME', this%mf6_input%mempath)
389  call mem_allocate(this%param_tasnames, lentimeseriesname, 0, &
390  'PARAMTASNAME', this%mf6_input%mempath)
391  end if

◆ tas_links_create()

subroutine layerarrayloadmodule::tas_links_create ( class(layerarrayloadtype), intent(inout)  this,
integer(i4b), intent(in)  inunit 
)

Definition at line 395 of file Mf6FileLayerArray.f90.

398  class(LayerArrayLoadType), intent(inout) :: this
399  integer(I4B), intent(in) :: inunit
400  type(InputParamDefinitionType), pointer :: idt
401  ! non-contiguous because a slice of bound is passed
402  real(DP), dimension(:), pointer :: auxArrayPtr, bndArrayPtr
403  real(DP), dimension(:), pointer, contiguous :: bound
404  integer(I4B), dimension(:), pointer, contiguous :: nodelist
405  character(len=LENTIMESERIESNAME) :: tas_name
406  character(len=LENAUXNAME) :: aux_name
407  logical :: convertFlux
408  integer(I4B) :: n
409 
410  ! initialize
411  nullify (auxarrayptr)
412  nullify (bndarrayptr)
413  nullify (nodelist)
414  convertflux = .false.
415 
416  ! Create AUX Time Array Series links
417  do n = 1, this%ctx%naux
418  tas_name = this%aux_tasnames(n)
419  if (tas_name /= '') then
420  ! set auxvar pointer
421  auxarrayptr => this%ctx%auxvar(n, :)
422  aux_name = this%ctx%auxname_cst(n)
423  call this%tasmanager%MakeTasLink(this%mf6_input%subcomponent_name, &
424  auxarrayptr, this%ctx%iprpak, &
425  tas_name, aux_name, convertflux, &
426  nodelist, inunit)
427  end if
428  end do
429 
430  ! Create BND Time Array Series links
431  do n = 1, this%nparam
432  ! assign param definition pointer
433  idt => get_param_definition_type(this%mf6_input%param_dfns, &
434  this%mf6_input%component_type, &
435  this%mf6_input%subcomponent_type, &
436  'PERIOD', this%param_names(n), &
437  this%input_name)
438  if (idt%timeseries) then
439  if (this%param_reads(n)%invar == 2) then
440  tas_name = this%param_tasnames(n)
441  call mem_setptr(bound, idt%mf6varname, this%mf6_input%mempath)
442  ! set bound pointer
443  bndarrayptr => bound(:)
444  call this%tasmanager%MakeTasLink(this%mf6_input%subcomponent_name, &
445  bndarrayptr, &
446  this%ctx%iprpak, &
447  tas_name, idt%mf6varname, &
448  convertflux, nodelist, inunit)
449  end if
450  end if
451  end do
Here is the call graph for this function: