MODFLOW 6  version 6.7.0.dev2
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:47
Here is the call graph for this function:

◆ destroy()

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

Definition at line 221 of file Mf6FileLayerArray.f90.

222  class(LayerArrayLoadType), intent(inout) :: this
223  !
224  ! deallocate tasmanager
225  call this%tasmanager%da()
226  deallocate (this%tasmanager)
227  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 255 of file Mf6FileLayerArray.f90.

257  class(LayerArrayLoadType) :: this
258  character(len=*), intent(in) :: varname
259  character(len=*), intent(in) :: input_name
260  type(CharacterStringType), dimension(:), pointer, &
261  contiguous :: charstr1d
262  integer(I4B) :: n
263  call mem_setptr(charstr1d, varname, this%mf6_input%mempath)
264  do n = 1, size(charstr1d)
265  charstr1d(n) = ''
266  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 293 of file Mf6FileLayerArray.f90.

294  use tdismodule, only: kper
296  use arrayhandlersmodule, only: ifind
303  use idmloggermodule, only: idm_log_var
304  class(LayerArrayLoadType), intent(inout) :: this
305  type(BlockParserType), intent(in) :: parser
306  type(InputParamDefinitionType), intent(in) :: idt
307  character(len=*), intent(in) :: mempath
308  logical(LGP), intent(in) :: netcdf
309  integer(I4B), intent(in) :: iaux
310  integer(I4B), dimension(:), pointer, contiguous :: int1d
311  real(DP), dimension(:), pointer, contiguous :: dbl1d
312  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
313  integer(I4B) :: iparam, n
314 
315  select case (idt%datatype)
316  case ('INTEGER1D')
317  call mem_setptr(int1d, idt%mf6varname, mempath)
318  if (netcdf) then
319  call netcdf_read_array(int1d, this%ctx%mshape, idt, &
320  this%mf6_input, this%nc_vars, this%input_name, &
321  this%iout, kper)
322  else
323  call read_int1d(parser, int1d, idt%mf6varname)
324  end if
325  call idm_log_var(int1d, idt%tagname, mempath, this%iout)
326  case ('DOUBLE1D')
327  call mem_setptr(dbl1d, idt%mf6varname, mempath)
328  if (netcdf) then
329  call netcdf_read_array(dbl1d, this%ctx%mshape, idt, &
330  this%mf6_input, this%nc_vars, this%input_name, &
331  this%iout, kper)
332  else
333  call read_dbl1d(parser, dbl1d, idt%mf6varname)
334  end if
335  call idm_log_var(dbl1d, idt%tagname, mempath, this%iout)
336  case ('DOUBLE2D')
337  call mem_setptr(dbl2d, idt%mf6varname, mempath)
338  allocate (dbl1d(this%ctx%ncpl))
339  if (netcdf) then
340  call netcdf_read_array(dbl1d, this%ctx%mshape, idt, &
341  this%mf6_input, this%nc_vars, this%input_name, &
342  this%iout, kper, iaux)
343  else
344  call read_dbl1d(parser, dbl1d, idt%mf6varname)
345  end if
346  do n = 1, this%ctx%ncpl
347  dbl2d(iaux, n) = dbl1d(n)
348  end do
349  call idm_log_var(dbl1d, idt%tagname, mempath, this%iout)
350  deallocate (dbl1d)
351  case default
352  errmsg = 'IDM unimplemented. LayerArrayLoad::param_load &
353  &datatype='//trim(idt%datatype)
354  call store_error(errmsg)
355  call store_error_filename(this%input_name)
356  end select
357 
358  ! if param is tracked set read state
359  iparam = ifind(this%param_names, idt%tagname)
360  if (iparam > 0) then
361  this%param_reads(iparam)%invar = 1
362  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)
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
This module contains the InputDefinitionModule.
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
Here is the call graph for this function:

◆ params_alloc()

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

Definition at line 269 of file Mf6FileLayerArray.f90.

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

◆ reset()

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

Definition at line 230 of file Mf6FileLayerArray.f90.

231  class(LayerArrayLoadType), intent(inout) :: this
232  integer(I4B) :: n, m
233 
234  if (this%tas_active /= 0) then
235  ! reset tasmanager
236  call this%tasmanager%reset(this%mf6_input%subcomponent_name)
237  ! reinitialize tas name arrays
238  call this%init_charstr1d('AUXTASNAME', this%input_name)
239  call this%init_charstr1d('PARAMTASNAME', this%input_name)
240  end if
241 
242  do n = 1, this%nparam
243  ! reset read state
244  this%param_reads(n)%invar = 0
245  end do
246 
247  ! explicitly reset auxvar array each period
248  do m = 1, this%ctx%ncpl
249  do n = 1, this%ctx%naux
250  this%ctx%auxvar(n, m) = dzero
251  end do
252  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  ! cycle to next input param
187  cycle
188  else
189  ! TODO: throw error
190  end if
191  else if (keyword == 'NETCDF') then
192  netcdf = .true.
193  end if
194 
195  ! read and load the parameter
196  call this%param_load(parser, idt, this%mf6_input%mempath, netcdf, iaux)
197  end do
198 
199  ! check if layer index variable was read
200  ! TODO: assumes layer index variable is always in scope
201  if (this%param_reads(1)%invar == 0) then
202  ! set to default of 1 without updating invar
203  idt => get_param_definition_type(this%mf6_input%param_dfns, &
204  this%mf6_input%component_type, &
205  this%mf6_input%subcomponent_type, &
206  'PERIOD', this%param_names(1), &
207  this%input_name)
208  call mem_setptr(int1d, idt%mf6varname, this%mf6_input%mempath)
209  int1d = 1
210  end if
211 
212  if (this%tas_active /= 0) then
213  call this%tas_links_create(parser%iuactive)
214  end if
215 
216  ! log lst file header
217  call idm_log_close(this%mf6_input%component_name, &
218  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 365 of file Mf6FileLayerArray.f90.

367  class(LayerArrayLoadType), intent(inout) :: this
368 
369  ! count params other than AUX
370  if (this%tas_active /= 0) then
371  call mem_allocate(this%aux_tasnames, lentimeseriesname, &
372  this%ctx%naux, 'AUXTASNAME', &
373  this%mf6_input%mempath)
374  call mem_allocate(this%param_tasnames, lentimeseriesname, this%nparam, &
375  'PARAMTASNAME', this%mf6_input%mempath)
376  call this%init_charstr1d('AUXTASNAME', this%input_name)
377  call this%init_charstr1d('PARAMTASNAME', this%input_name)
378  else
379  call mem_allocate(this%aux_tasnames, lentimeseriesname, 0, &
380  'AUXTASNAME', this%mf6_input%mempath)
381  call mem_allocate(this%param_tasnames, lentimeseriesname, 0, &
382  'PARAMTASNAME', this%mf6_input%mempath)
383  end if

◆ tas_links_create()

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

Definition at line 387 of file Mf6FileLayerArray.f90.

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