MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
gridarrayloadmodule Module Reference

This module contains the GridArrayLoadModule. More...

Data Types

type  gridarrayloadtype
 Ascii grid based 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 params_alloc (this)
 
subroutine param_load (this, parser, idt, mempath, layered, netcdf, iaux)
 

Detailed Description

This module contains the routines for reading period block array based input associated with the full grid, such as with the GHBA package.

Function/Subroutine Documentation

◆ ad()

subroutine gridarrayloadmodule::ad ( class(gridarrayloadtype), intent(inout)  this)
private

Definition at line 97 of file Mf6FileGridArray.f90.

98  class(GridArrayLoadType), intent(inout) :: this

◆ ainit()

subroutine gridarrayloadmodule::ainit ( class(gridarrayloadtype), 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 46 of file Mf6FileGridArray.f90.

52  class(GridArrayLoadType), intent(inout) :: this
53  type(ModflowInputType), intent(in) :: mf6_input
54  character(len=*), intent(in) :: component_name
55  character(len=*), intent(in) :: component_input_name
56  character(len=*), intent(in) :: input_name
57  integer(I4B), intent(in) :: iperblock
58  type(BlockParserType), pointer, intent(inout) :: parser
59  integer(I4B), intent(in) :: iout
60  type(LoadMf6FileType) :: loader
61  integer(I4B) :: isize
62  integer(I4B), pointer :: maxbound
63 
64  ! initialize base type
65  call this%DynamicPkgLoadType%init(mf6_input, component_name, &
66  component_input_name, &
67  input_name, iperblock, iout)
68  ! initialize
69  this%iout = iout
70 
71  ! load static input
72  call loader%load(parser, mf6_input, this%nc_vars, this%input_name, iout)
73 
74  ! maxbound is optional
75  call get_isize('MAXBOUND', mf6_input%mempath, isize)
76  if (isize < 0) then
77  ! set maxbound to grid nodes
78  call mem_allocate(maxbound, 'MAXBOUND', mf6_input%mempath)
79  maxbound = product(loader%mshape)
80  end if
81 
82  ! initialize input context memory
83  call this%ctx%init(mf6_input)
84 
85  ! allocate user nodelist
86  call mem_allocate(this%nodeulist, this%ctx%maxbound, &
87  'NODEULIST', mf6_input%mempath)
88 
89  ! allocate dfn params
90  call this%params_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 gridarrayloadmodule::destroy ( class(gridarrayloadtype), intent(inout)  this)

Definition at line 164 of file Mf6FileGridArray.f90.

166  class(GridArrayLoadType), intent(inout) :: this
167  call mem_deallocate(this%nodeulist)

◆ df()

subroutine gridarrayloadmodule::df ( class(gridarrayloadtype), intent(inout)  this)

Definition at line 93 of file Mf6FileGridArray.f90.

94  class(GridArrayLoadType), intent(inout) :: this

◆ param_load()

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

Definition at line 207 of file Mf6FileGridArray.f90.

208  use tdismodule, only: kper
209  use constantsmodule, only: dnodata
210  use arrayhandlersmodule, only: ifind
218  use idmloggermodule, only: idm_log_var
219  class(GridArrayLoadType), intent(inout) :: this
220  type(BlockParserType), intent(in) :: parser
221  type(InputParamDefinitionType), intent(in) :: idt
222  character(len=*), intent(in) :: mempath
223  logical(LGP), intent(in) :: layered
224  logical(LGP), intent(in) :: netcdf
225  real(DP), dimension(:), pointer, contiguous :: dbl1d, nodes
226  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
227  integer(I4B), dimension(:), allocatable :: layer_shape
228  integer(I4B) :: iaux, iparam, n, nlay, nnode
229 
230  nnode = 0
231 
232  select case (idt%datatype)
233  case ('DOUBLE1D')
234  call mem_setptr(dbl1d, idt%mf6varname, mempath)
235  allocate (nodes(this%ctx%nodes))
236  if (netcdf) then
237  call netcdf_read_array(nodes, this%ctx%mshape, idt, &
238  this%mf6_input, this%nc_vars, this%input_name, &
239  this%iout, kper)
240  else if (layered) then
241  call get_layered_shape(this%ctx%mshape, nlay, layer_shape)
242  call read_dbl1d_layered(parser, nodes, idt%mf6varname, nlay, layer_shape)
243  else
244  call read_dbl1d(parser, nodes, idt%mf6varname)
245  end if
246 
247  call idm_log_var(nodes, idt%tagname, mempath, this%iout)
248 
249  do n = 1, this%ctx%nodes
250  if (nodes(n) /= dnodata) then
251  nnode = nnode + 1
252  dbl1d(nnode) = nodes(n)
253  if (this%ctx%nbound == 0) then
254  this%nodeulist(nnode) = n
255  else if (this%nodeulist(nnode) /= n) then
256  write (errmsg, '(a,i0)') 'Grid input position mismatch param='// &
257  trim(idt%tagname)//', period=', kper
258  call store_error(errmsg)
259  call store_error_filename(this%input_name)
260  end if
261  end if
262  end do
263  deallocate (nodes)
264  case ('DOUBLE2D')
265  call mem_setptr(dbl2d, idt%mf6varname, mempath)
266  allocate (nodes(this%ctx%nodes))
267 
268  if (netcdf) then
269  call netcdf_read_array(nodes, this%ctx%mshape, idt, &
270  this%mf6_input, this%nc_vars, this%input_name, &
271  this%iout, kper, iaux)
272  else if (layered) then
273  call get_layered_shape(this%ctx%mshape, nlay, layer_shape)
274  call read_dbl1d_layered(parser, nodes, idt%mf6varname, nlay, layer_shape)
275  else
276  call read_dbl1d(parser, nodes, idt%mf6varname)
277  end if
278 
279  call idm_log_var(nodes, idt%tagname, mempath, this%iout)
280 
281  do n = 1, this%ctx%nodes
282  if (nodes(n) /= dnodata) then
283  nnode = nnode + 1
284  dbl2d(iaux, nnode) = nodes(n)
285  if (this%ctx%nbound == 0) then
286  this%nodeulist(nnode) = n
287  else if (this%nodeulist(nnode) /= n) then
288  write (errmsg, '(a,i0)') 'Grid input position mismatch param='// &
289  trim(idt%tagname)//', period=', kper
290  call store_error(errmsg)
291  call store_error_filename(this%input_name)
292  end if
293  end if
294  end do
295  deallocate (nodes)
296  case default
297  errmsg = 'IDM unimplemented. GridArrayLoad::param_load &
298  &datatype='//trim(idt%datatype)
299  call store_error(errmsg)
300  call store_error_filename(this%input_name)
301  end select
302 
303  ! set nbound
304  if (this%ctx%nbound == 0) this%ctx%nbound = nnode
305 
306  ! if param is tracked set read state
307  iparam = ifind(this%param_names, idt%tagname)
308  if (iparam > 0) then
309  this%param_reads(iparam)%invar = 1
310  end if
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
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_dbl1d_layered(parser, dbl1d, aname, nlay, layer_shape)
This module contains the LoadNCInputModule.
Definition: LoadNCInput.F90:7
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
subroutine, public get_layered_shape(mshape, nlay, layer_shape)
subroutine, public get_shape_from_string(shape_string, array_shape, memoryPath)
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ params_alloc()

subroutine gridarrayloadmodule::params_alloc ( class(gridarrayloadtype), intent(inout)  this)

Definition at line 183 of file Mf6FileGridArray.f90.

184  class(GridArrayLoadType), intent(inout) :: this
185  character(len=LENVARNAME) :: rs_varname
186  integer(I4B), pointer :: intvar
187  integer(I4B) :: iparam
188 
189  ! set in scope param names
190  call this%ctx%tags(this%param_names, this%nparam, this%input_name, &
191  create=.true.)
192  call this%ctx%allocate_arrays()
193 
194  ! allocate and set param_reads pointer array
195  allocate (this%param_reads(this%nparam))
196 
197  ! store read state variable pointers
198  do iparam = 1, this%nparam
199  ! allocate and store name of read state variable
200  rs_varname = this%ctx%rsv_alloc(this%param_names(iparam))
201  call mem_setptr(intvar, rs_varname, this%mf6_input%mempath)
202  this%param_reads(iparam)%invar => intvar
203  this%param_reads(iparam)%invar = 0
204  end do

◆ reset()

subroutine gridarrayloadmodule::reset ( class(gridarrayloadtype), intent(inout)  this)

Definition at line 170 of file Mf6FileGridArray.f90.

171  use constantsmodule, only: dnodata
172  class(GridArrayLoadType), intent(inout) :: this
173  integer(I4B) :: n
174 
175  this%ctx%nbound = 0
176 
177  do n = 1, this%nparam
178  ! reset read state
179  this%param_reads(n)%invar = 0
180  end do

◆ rp()

subroutine gridarrayloadmodule::rp ( class(gridarrayloadtype), intent(inout)  this,
type(blockparsertype), intent(inout), pointer  parser 
)
private

Definition at line 101 of file Mf6FileGridArray.f90.

105  use arrayhandlersmodule, only: ifind
108  class(GridArrayLoadType), intent(inout) :: this
109  type(BlockParserType), pointer, intent(inout) :: parser
110  logical(LGP) :: endOfBlock, netcdf, layered
111  character(len=LINELENGTH) :: keyword, param_tag
112  type(InputParamDefinitionType), pointer :: idt
113  integer(I4B) :: iaux
114 
115  ! reset for this period
116  call this%reset()
117 
118  ! log lst file header
119  call idm_log_header(this%mf6_input%component_name, &
120  this%mf6_input%subcomponent_name, this%iout)
121 
122  ! read array block
123  do
124  ! initialize
125  iaux = 0
126  netcdf = .false.
127  layered = .false.
128 
129  ! read next line
130  call parser%GetNextLine(endofblock)
131  if (endofblock) exit
132  ! read param_tag
133  call parser%GetStringCaps(param_tag)
134 
135  ! is param tag an auxvar?
136  iaux = ifind_charstr(this%ctx%auxname_cst, param_tag)
137 
138  ! any auvxar corresponds to the definition tag 'AUX'
139  if (iaux > 0) param_tag = 'AUX'
140 
141  ! set input definition
142  idt => get_param_definition_type(this%mf6_input%param_dfns, &
143  this%mf6_input%component_type, &
144  this%mf6_input%subcomponent_type, &
145  'PERIOD', param_tag, this%input_name)
146  ! look for Layered and NetCDF keywords
147  call parser%GetStringCaps(keyword)
148  if (keyword == 'LAYERED' .and. idt%layered) then
149  layered = .true.
150  else if (keyword == 'NETCDF') then
151  netcdf = .true.
152  end if
153 
154  ! read and load the parameter
155  call this%param_load(parser, idt, this%mf6_input%mempath, layered, &
156  netcdf, iaux)
157  end do
158 
159  ! log lst file header
160  call idm_log_close(this%mf6_input%component_name, &
161  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
integer(i4b) function, public ifind_charstr(array, str)
Here is the call graph for this function: