MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
BoundaryPackageExt.f90
Go to the documentation of this file.
1 !> @brief This module contains the extended boundary package
2 !!
3 !! This module contains the extended boundary type that itself
4 !! should be extended by model boundary packages that have been
5 !! updated to source static and dynamic input data from the
6 !! input context.
7 !!
8 !<
10 
11  use kindmodule, only: dp, lgp, i4b
13  use obsmodule, only: obs_cr
14  use simvariablesmodule, only: errmsg
16  use bndmodule, only: bndtype
17  use geomutilmodule, only: get_node, get_ijk
18 
19  implicit none
20 
21  private
22  public :: bndexttype
23 
24  !> @ brief BndExtType
25  !!
26  !! Generic extended boundary package type. This derived type can be
27  !! overridden to define concrete boundary package types that source
28  !! all input from the input context.
29  !<
30  type, extends(bndtype) :: bndexttype
31  ! -- characters
32  ! -- scalars
33  integer(I4B), pointer :: iper
34  logical(LGP), pointer :: readarraygrid
35  logical(LGP), pointer :: readasarrays
36  ! -- arrays
37  integer(I4B), dimension(:, :), pointer, contiguous :: cellid => null() !< input user cellid list
38  integer(I4B), dimension(:), pointer, contiguous :: nodeulist => null() !< input user nodelist
39  contains
40  procedure :: bnd_df => bndext_df
41  procedure :: bnd_rp => bndext_rp
42  procedure :: bnd_da => bndext_da
43  procedure :: allocate_scalars => bndext_allocate_scalars
44  procedure :: allocate_arrays => bndext_allocate_arrays
45  procedure :: source_options
46  procedure :: source_dimensions
47  procedure :: log_options
48  procedure :: cellid_to_nlist
49  procedure :: nodeu_to_nlist
50  procedure :: layarr_to_nlist
51  procedure :: default_nodelist
52  procedure :: check_cellid
53  procedure :: write_list
54  procedure :: bound_value
55  end type bndexttype
56 
57  !> @ brief BndExtFoundType
58  !!
59  !! This type is used to simplify the tracking of common parameters
60  !! that are sourced from the input context.
61  !<
63  logical :: naux = .false.
64  logical :: ipakcb = .false.
65  logical :: iprpak = .false.
66  logical :: iprflow = .false.
67  logical :: boundnames = .false.
68  logical :: auxmultname = .false.
69  logical :: inewton = .false.
70  logical :: auxiliary = .false.
71  logical :: maxbound = .false.
72  end type bndextfoundtype
73 
74 contains
75 
76  !> @ brief Define boundary package options and dimensions
77  !!
78  !! Define base boundary package options and dimensions for
79  !! a model boundary package.
80  !!
81  !<
82  subroutine bndext_df(this, neq, dis)
83  ! -- modules
84  use basedismodule, only: disbasetype
88  ! -- dummy variables
89  class(bndexttype), intent(inout) :: this !< BndExtType object
90  integer(I4B), intent(inout) :: neq !< number of equations
91  class(disbasetype), pointer :: dis !< discretization object
92  !
93  ! -- set pointer to dis object for the model
94  this%dis => dis
95  !
96  ! -- Create time series managers
97  ! -- Not in use by this type but BndType uses and deallocates
98  call tsmanager_cr(this%TsManager, this%iout)
99  call tasmanager_cr(this%TasManager, dis, this%name_model, this%iout)
100  !
101  ! -- create obs package
102  call obs_cr(this%obs, this%inobspkg)
103  !
104  ! -- Write information to model list file
105  write (this%iout, 1) trim(this%filtyp), trim(adjustl(this%text)), &
106  this%input_mempath
107 1 format(1x, /1x, a, ' -- ', a, ' PACKAGE, VERSION 8, 2/22/2014', &
108  ' INPUT READ FROM MEMPATH: ', a)
109  !
110  ! -- source options
111  call this%source_options()
112  !
113  ! -- Define time series managers
114  call this%tsmanager%tsmanager_df()
115  call this%tasmanager%tasmanager_df()
116  !
117  ! -- source dimensions
118  call this%source_dimensions()
119  !
120  ! -- update package moffset for packages that add rows
121  if (this%npakeq > 0) then
122  this%ioffset = neq - this%dis%nodes
123  end if
124  !
125  ! -- update neq
126  neq = neq + this%npakeq
127  !
128  ! -- Store information needed for observations
129  if (this%bnd_obs_supported()) then
130  call this%obs%obs_df(this%iout, this%packName, this%filtyp, this%dis)
131  call this%bnd_df_obs()
132  end if
133  !
134  ! -- Call define_listlabel to construct the list label that is written
135  ! when PRINT_INPUT option is used.
136  call this%define_listlabel()
137  end subroutine bndext_df
138 
139  subroutine bndext_rp(this)
140  ! -- modules
141  use tdismodule, only: kper
144  ! -- dummy variables
145  class(bndexttype), intent(inout) :: this !< BndExtType object
146  ! -- local variables
147  logical(LGP) :: found
148  integer(I4B) :: n
149  !
150  if (this%iper /= kper) return
151 
152  if (.not. this%readasarrays) then
153  ! -- copy nbound from input context
154  call mem_set_value(this%nbound, 'NBOUND', this%input_mempath, &
155  found)
156  end if
157 
158  if (this%readarraygrid) then
159  call this%nodeu_to_nlist()
160  else if (this%readasarrays) then
161  call this%layarr_to_nlist()
162  else
163  call this%cellid_to_nlist()
164  !
165  ! -- update boundname string list
166  if (this%inamedbound /= 0) then
167  do n = 1, size(this%boundname_cst)
168  this%boundname(n) = this%boundname_cst(n)
169  end do
170  end if
171  end if
172  end subroutine bndext_rp
173 
174  !> @ brief Deallocate package memory
175  !<
176  subroutine bndext_da(this)
177  ! -- modules
179  ! -- dummy variables
180  class(bndexttype) :: this !< BndExtType object
181  !
182  ! -- deallocate checkin paths
183  call mem_deallocate(this%cellid, 'CELLID', this%memoryPath)
184  call mem_deallocate(this%nodeulist, 'NODEULIST', this%memoryPath)
185  call mem_deallocate(this%boundname_cst, 'BOUNDNAME_IDM', this%memoryPath)
186  call mem_deallocate(this%auxvar, 'AUXVAR_IDM', this%memoryPath)
187  !
188  ! -- reassign pointers for base class _da
189  call mem_setptr(this%boundname_cst, 'BOUNDNAME_CST', this%memoryPath)
190  call mem_setptr(this%auxvar, 'AUXVAR', this%memoryPath)
191  !
192  ! -- scalars
193  deallocate (this%readarraygrid)
194  deallocate (this%readasarrays)
195  nullify (this%readarraygrid)
196  nullify (this%readasarrays)
197  nullify (this%iper)
198  !
199  ! -- deallocate
200  call this%BndType%bnd_da()
201  end subroutine bndext_da
202 
203  !> @ brief Allocate package scalars
204  !!
205  !! Allocate and initialize base boundary package scalars. This method
206  !! only needs to be overridden if additional scalars are defined
207  !! for a specific package.
208  !!
209  !<
210  subroutine bndext_allocate_scalars(this)
211  ! -- modules
213  ! -- dummy variables
214  class(bndexttype) :: this !< BndExtType object
215  ! -- local variables
216  !
217  ! -- allocate base BndType scalars
218  call this%BndType%allocate_scalars()
219  !
220  ! -- set IPER pointer
221  call mem_setptr(this%iper, 'IPER', this%input_mempath)
222 
223  ! -- allocate internal scalars
224  allocate (this%readarraygrid)
225  allocate (this%readasarrays)
226 
227  ! -- initialize internal scalars
228  this%readarraygrid = .false.
229  this%readasarrays = .false.
230  end subroutine bndext_allocate_scalars
231 
232  !> @ brief Allocate package arrays
233  !!
234  !! Allocate and initialize base boundary package arrays. This method
235  !! only needs to be overridden if additional arrays are defined
236  !! for a specific package.
237  !!
238  !<
239  subroutine bndext_allocate_arrays(this, nodelist, auxvar)
240  ! -- modules
242  ! -- dummy variables
243  class(bndexttype) :: this !< BndExtType object
244  ! -- local variables
245  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist !< package nodelist
246  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar !< package aux variable array
247  !
248  ! -- allocate base BndType arrays
249  call this%BndType%allocate_arrays(nodelist, auxvar)
250  !
251  ! -- set input context pointers
252  call mem_setptr(this%cellid, 'CELLID', this%input_mempath)
253  call mem_setptr(this%nodeulist, 'NODEULIST', this%input_mempath)
254  call mem_setptr(this%boundname_cst, 'BOUNDNAME', this%input_mempath)
255  !
256  ! -- checkin input context pointers
257  call mem_checkin(this%cellid, 'CELLID', this%memoryPath, &
258  'CELLID', this%input_mempath)
259  call mem_checkin(this%nodeulist, 'NODEULIST', this%memoryPath, &
260  'NODEULIST', this%input_mempath)
261  call mem_checkin(this%boundname_cst, lenboundname, 'BOUNDNAME_IDM', &
262  this%memoryPath, 'BOUNDNAME', this%input_mempath)
263  !
264  if (present(auxvar)) then
265  ! no-op
266  else
267  ! -- set auxvar input context pointer
268  call mem_setptr(this%auxvar, 'AUXVAR', this%input_mempath)
269  !
270  ! -- checkin auxvar input context pointer
271  call mem_checkin(this%auxvar, 'AUXVAR_IDM', this%memoryPath, &
272  'AUXVAR', this%input_mempath)
273  end if
274  end subroutine bndext_allocate_arrays
275 
276  !> @ brief Source package options from input context
277  !<
278  subroutine source_options(this)
279  ! -- modules
280  use memorymanagermodule, only: mem_reallocate, mem_setptr !, get_isize
285  ! -- dummy variables
286  class(bndexttype), intent(inout) :: this !< BndExtType object
287  ! -- local variables
288  type(bndextfoundtype) :: found
289  logical(LGP) :: found_readarr
290  character(len=LENAUXNAME) :: sfacauxname
291  integer(I4B) :: n
292  !
293  ! -- update defaults with idm sourced values
294  call mem_set_value(this%naux, 'NAUX', this%input_mempath, found%naux)
295  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
296  call mem_set_value(this%iprpak, 'IPRPAK', this%input_mempath, found%iprpak)
297  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
298  call mem_set_value(this%inamedbound, 'BOUNDNAMES', this%input_mempath, &
299  found%boundnames)
300  call mem_set_value(sfacauxname, 'AUXMULTNAME', this%input_mempath, &
301  found%auxmultname)
302  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
303  call mem_set_value(this%readarraygrid, 'READARRAYGRID', this%input_mempath, &
304  found_readarr)
305  call mem_set_value(this%readasarrays, 'READASARRAYS', this%input_mempath, &
306  found_readarr)
307  !
308  ! -- log found options
309  call this%log_options(found, sfacauxname)
310  !
311  ! -- reallocate aux arrays if aux variables provided
312  if (found%naux .and. this%naux > 0) then
313  call mem_reallocate(this%auxname, lenauxname, this%naux, &
314  'AUXNAME', this%memoryPath)
315  call mem_reallocate(this%auxname_cst, lenauxname, this%naux, &
316  'AUXNAME_CST', this%memoryPath)
317  call mem_set_value(this%auxname_cst, 'AUXILIARY', this%input_mempath, &
318  found%auxiliary)
319  !
320  do n = 1, this%naux
321  this%auxname(n) = this%auxname_cst(n)
322  end do
323  end if
324  !
325  ! -- save flows option active
326  if (found%ipakcb) this%ipakcb = -1
327  !
328  ! -- auxmultname provided
329  if (found%auxmultname) this%iauxmultcol = -1
330  !
331  !
332  ! -- enforce 0 or 1 OBS6_FILENAME entries in option block
333  if (filein_fname(this%obs%inputFilename, 'OBS6_FILENAME', &
334  this%input_mempath, this%input_fname)) then
335  this%obs%active = .true.
336  this%obs%inUnitObs = getunit()
337  call openfile(this%obs%inUnitObs, this%iout, this%obs%inputFilename, 'OBS')
338  end if
339  !
340  ! -- no newton specified
341  if (found%inewton) this%inewton = 0
342  !
343  ! -- AUXMULTNAME was specified, so find column of auxvar that will be multiplier
344  if (this%iauxmultcol < 0) then
345  !
346  ! -- Error if no aux variable specified
347  if (this%naux == 0) then
348  write (errmsg, '(a,2(1x,a))') &
349  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
350  'but no AUX variables specified.'
351  call store_error(errmsg)
352  end if
353  !
354  ! -- Assign mult column
355  this%iauxmultcol = 0
356  do n = 1, this%naux
357  if (sfacauxname == this%auxname(n)) then
358  this%iauxmultcol = n
359  exit
360  end if
361  end do
362  !
363  ! -- Error if aux variable cannot be found
364  if (this%iauxmultcol == 0) then
365  write (errmsg, '(a,2(1x,a))') &
366  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
367  'but no AUX variable found with this name.'
368  call store_error(errmsg)
369  end if
370  end if
371  !
372  if (this%readasarrays) then
373  if (.not. this%dis%supports_layers()) then
374  errmsg = 'READASARRAYS option is not compatible with selected'// &
375  ' discretization type.'
376  call store_error(errmsg)
377  call store_error_filename(this%input_fname)
378  end if
379  end if
380  !
381  ! -- terminate if errors were detected
382  if (count_errors() > 0) then
383  call store_error_filename(this%input_fname)
384  end if
385  end subroutine source_options
386 
387  !> @ brief Log package options
388  !<
389  subroutine log_options(this, found, sfacauxname)
390  ! -- modules
391  ! -- dummy variables
392  class(bndexttype), intent(inout) :: this !< BndExtType object
393  type(bndextfoundtype), intent(in) :: found
394  character(len=*), intent(in) :: sfacauxname
395  ! -- local variables
396  ! -- format
397  character(len=*), parameter :: fmtreadasarrays = &
398  &"(4x, 'PACKAGE INPUT WILL BE READ AS LAYER ARRAYS.')"
399  character(len=*), parameter :: fmtreadarraygrid = &
400  &"(4x, 'PACKAGE INPUT WILL BE READ AS GRID ARRAYS.')"
401  character(len=*), parameter :: fmtflow = &
402  &"(4x, 'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL')"
403  !
404  ! -- log found options
405  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
406  //' BASE OPTIONS'
407  !
408  if (this%readasarrays) then
409  write (this%iout, fmtreadasarrays)
410  end if
411  !
412  if (this%readarraygrid) then
413  write (this%iout, fmtreadarraygrid)
414  end if
415  !
416  if (found%ipakcb) then
417  write (this%iout, fmtflow)
418  end if
419  !
420  if (found%iprpak) then
421  write (this%iout, '(4x,a)') &
422  'LISTS OF '//trim(adjustl(this%text))//' CELLS WILL BE PRINTED.'
423  end if
424  !
425  if (found%iprflow) then
426  write (this%iout, '(4x,a)') trim(adjustl(this%text))// &
427  ' FLOWS WILL BE PRINTED TO LISTING FILE.'
428  end if
429  !
430  if (found%boundnames) then
431  write (this%iout, '(4x,a)') trim(adjustl(this%text))// &
432  ' BOUNDARIES HAVE NAMES IN LAST COLUMN.'
433  end if
434  !
435  if (found%auxmultname) then
436  write (this%iout, '(4x,a,a)') &
437  'AUXILIARY MULTIPLIER NAME: ', sfacauxname
438  end if
439  !
440  if (found%inewton) then
441  write (this%iout, '(4x,a)') &
442  'NEWTON-RAPHSON method disabled for unconfined cells'
443  end if
444  !
445  ! -- close logging block
446  write (this%iout, '(1x,a)') &
447  'END OF '//trim(adjustl(this%text))//' BASE OPTIONS'
448  end subroutine log_options
449 
450  !> @ brief Source package dimensions from input context
451  !<
452  subroutine source_dimensions(this)
454  ! -- dummy variables
455  class(bndexttype), intent(inout) :: this !< BndExtType object
456  ! -- local variables
457  type(bndextfoundtype) :: found
458  !
459  if (this%readasarrays) then
460  this%maxbound = this%dis%get_ncpl()
461  else
462  ! -- open dimensions logging block
463  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text))// &
464  ' BASE DIMENSIONS'
465 
466  ! -- update defaults with idm sourced values
467  call mem_set_value(this%maxbound, 'MAXBOUND', this%input_mempath, &
468  found%maxbound)
469 
470  write (this%iout, '(4x,a,i7)') 'MAXBOUND = ', this%maxbound
471 
472  ! -- close logging block
473  write (this%iout, '(1x,a)') &
474  'END OF '//trim(adjustl(this%text))//' BASE DIMENSIONS'
475  end if
476  !
477  ! -- verify dimensions were set
478  if (this%maxbound <= 0) then
479  write (errmsg, '(a)') 'MAXBOUND must be an integer greater than zero.'
480  call store_error(errmsg)
481  call store_error_filename(this%input_fname)
482  end if
483  end subroutine source_dimensions
484 
485  !> @ brief Update package nodelist
486  !!
487  !! Convert period updated cellids to node numbers.
488  !!
489  !<
490  subroutine cellid_to_nlist(this)
491  ! -- modules
492  use simvariablesmodule, only: errmsg
493  ! -- dummy
494  class(bndexttype) :: this !< BndExtType object
495  ! -- local
496  integer(I4B), dimension(:), pointer :: cellid
497  integer(I4B) :: n, nodeu, noder
498  character(len=LINELENGTH) :: nodestr
499  !
500  ! -- update nodelist
501  do n = 1, this%nbound
502  !
503  ! -- set cellid
504  cellid => this%cellid(:, n)
505  !
506  ! -- ensure cellid is valid, store an error otherwise
507  call this%check_cellid(n, cellid, this%dis%mshape, this%dis%ndim)
508  !
509  ! -- Determine user node number
510  if (this%dis%ndim == 1) then
511  nodeu = cellid(1)
512  elseif (this%dis%ndim == 2) then
513  nodeu = get_node(cellid(1), 1, cellid(2), &
514  this%dis%mshape(1), 1, &
515  this%dis%mshape(2))
516  else
517  nodeu = get_node(cellid(1), cellid(2), cellid(3), &
518  this%dis%mshape(1), &
519  this%dis%mshape(2), &
520  this%dis%mshape(3))
521  end if
522  !
523  ! -- update the nodelist
524  if (this%dis%nodes < this%dis%nodesuser) then
525  ! -- convert user to reduced node numbers
526  noder = this%dis%get_nodenumber(nodeu, 0)
527  if (noder <= 0) then
528  call this%dis%nodeu_to_string(nodeu, nodestr)
529  write (errmsg, *) &
530  ' Cell is outside active grid domain: '// &
531  trim(adjustl(nodestr))
532  call store_error(errmsg)
533  end if
534  this%nodelist(n) = noder
535  else
536  this%nodelist(n) = nodeu
537  end if
538  end do
539  !
540  ! -- exit if errors were found
541  if (count_errors() > 0) then
542  write (errmsg, *) count_errors(), ' errors encountered.'
543  call store_error(errmsg)
544  call store_error_filename(this%input_fname)
545  end if
546  end subroutine cellid_to_nlist
547 
548  !> @ brief Update package nodelist
549  !!
550  !! Convert period user nodes to reduced nodes
551  !!
552  !<
553  subroutine nodeu_to_nlist(this)
554  ! -- modules
556  ! -- dummy
557  class(bndexttype) :: this !< BndExtType object
558  integer(I4B) :: n, noder, nodeuser, ninactive
559 
560  ninactive = 0
561 
562  ! -- Set the nodelist
563  do n = 1, this%nbound
564  nodeuser = this%nodeulist(n)
565  noder = this%dis%get_nodenumber(nodeuser, 0)
566  if (noder > 0) then
567  this%nodelist(n) = noder
568  else
569  ninactive = ninactive + 1
570  end if
571  end do
572 
573  ! update nbound
574  this%nbound = this%nbound - ninactive
575  end subroutine nodeu_to_nlist
576 
577  !> @brief Update the nodelist based on layer number variable input
578  !!
579  !! This is a module scoped routine to check for I<filtyp>
580  !! input. If array input was provided, INI<filtyp> and I<filtyp>
581  !! will be allocated in the input context. If the read
582  !! state variable INI<filtyp> is set to 1 during this period
583  !! update, I<filtyp> input was read and is used here to update
584  !! the nodelist.
585  !!
586  !<
587  subroutine layarr_to_nlist(this)
588  ! -- modules
590  use constantsmodule, only: lenvarname
591  ! -- dummy
592  class(bndexttype) :: this !< BndExtType object
593  character(len=LENVARNAME) :: ilayname, inilayname
594  character(len=24) :: aname = ' LAYER OR NODE INDEX'
595  ! -- local
596  integer(I4B), dimension(:), contiguous, &
597  pointer :: ilay => null()
598  integer(I4B), pointer :: inilay => null()
599  !
600  ! set ilay and read state variable names
601  ilayname = 'I'//trim(this%filtyp)
602  inilayname = 'INI'//trim(this%filtyp)
603  !
604  ! -- set pointer to input context read state variable
605  call mem_setptr(inilay, inilayname, this%input_mempath)
606  !
607  ! -- check read state
608  if (inilay == 1) then
609  ! -- ilay variable was read this period
610  !
611  ! -- set pointer to input context layer index variable
612  call mem_setptr(ilay, ilayname, this%input_mempath)
613  !
614  ! -- update nodelist
615  call this%dis%nlarray_to_nodelist(ilay, this%nodelist, this%maxbound, &
616  this%nbound, aname)
617  end if
618  end subroutine layarr_to_nlist
619 
620  !> @brief Assign default nodelist when READASARRAYS is specified.
621  !!
622  !! Equivalent to reading layer number array as CONSTANT 1
623  !<
624  subroutine default_nodelist(this)
625  ! -- dummy
626  class(bndexttype) :: this
627  ! -- local
628  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nodeu, noder, ipos
629  !
630  if (this%readasarrays) then
631  !
632  ! -- set variables
633  if (this%dis%ndim == 3) then
634  nlay = this%dis%mshape(1)
635  nrow = this%dis%mshape(2)
636  ncol = this%dis%mshape(3)
637  elseif (this%dis%ndim == 2) then
638  nlay = this%dis%mshape(1)
639  nrow = 1
640  ncol = this%dis%mshape(2)
641  end if
642  !
643  ! -- Populate nodelist
644  ipos = 1
645  il = 1
646  do ir = 1, nrow
647  do ic = 1, ncol
648  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
649  noder = this%dis%get_nodenumber(nodeu, 0)
650  this%nodelist(ipos) = noder
651  ipos = ipos + 1
652  end do
653  end do
654  !
655  ! -- Assign nbound
656  this%nbound = ipos - 1
657  end if
658  end subroutine default_nodelist
659 
660  !> @ brief Check for valid cellid
661  !<
662  subroutine check_cellid(this, ii, cellid, mshape, ndim)
663  ! -- modules
664  use simvariablesmodule, only: errmsg
665  ! -- dummy
666  class(bndexttype) :: this !< BndExtType object
667  ! -- local
668  integer(I4B), intent(in) :: ii
669  integer(I4B), dimension(:), intent(in) :: cellid !< cellid
670  integer(I4B), dimension(:), intent(in) :: mshape !< model shape
671  integer(I4B), intent(in) :: ndim !< size of mshape
672  character(len=20) :: cellstr, mshstr
673  character(len=*), parameter :: fmterr = &
674  "('List entry ',i0,' contains cellid ',a,' but this cellid is invalid &
675  &for model with shape ', a)"
676  character(len=*), parameter :: fmtndim1 = &
677  "('(',i0,')')"
678  character(len=*), parameter :: fmtndim2 = &
679  "('(',i0,',',i0,')')"
680  character(len=*), parameter :: fmtndim3 = &
681  "('(',i0,',',i0,',',i0,')')"
682  select case (ndim)
683  case (1)
684  !
685  if (cellid(1) < 1 .or. cellid(1) > mshape(1)) then
686  write (cellstr, fmtndim1) cellid(1)
687  write (mshstr, fmtndim1) mshape(1)
688  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
689  call store_error(errmsg)
690  end if
691  !
692  case (2)
693  !
694  if (cellid(1) < 1 .or. cellid(1) > mshape(1) .or. &
695  cellid(2) < 1 .or. cellid(2) > mshape(2)) then
696  write (cellstr, fmtndim2) cellid(1), cellid(2)
697  write (mshstr, fmtndim2) mshape(1), mshape(2)
698  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
699  call store_error(errmsg)
700  end if
701  !
702  case (3)
703  !
704  if (cellid(1) < 1 .or. cellid(1) > mshape(1) .or. &
705  cellid(2) < 1 .or. cellid(2) > mshape(2) .or. &
706  cellid(3) < 1 .or. cellid(3) > mshape(3)) then
707  write (cellstr, fmtndim3) cellid(1), cellid(2), cellid(3)
708  write (mshstr, fmtndim3) mshape(1), mshape(2), mshape(3)
709  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
710  call store_error(errmsg)
711  end if
712  !
713  case default
714  end select
715  end subroutine check_cellid
716 
717  !> @ brief Log package list input
718  !!
719  !! Log period list based input. This routine requires a package specific
720  !! bound_value() routine to report accurate bound values.
721  !!
722  !<
723  subroutine write_list(this)
724  ! -- modules
727  use inputoutputmodule, only: ulstlb
728  use tablemodule, only: tabletype, table_cr
729  ! -- dummy
730  class(bndexttype) :: this !< BndExtType object
731  ! -- local
732  character(len=10) :: cpos
733  character(len=LINELENGTH) :: tag
734  character(len=LINELENGTH), allocatable, dimension(:) :: words
735  integer(I4B) :: ntabrows
736  integer(I4B) :: ntabcols
737  integer(I4B) :: ipos
738  integer(I4B) :: ii, jj, i, j, k, nod
739  integer(I4B) :: ldim
740  integer(I4B) :: naux
741  type(tabletype), pointer :: inputtab => null()
742  ! -- formats
743  character(len=LINELENGTH) :: fmtlstbn
744  !
745  ! -- Determine sizes
746  ldim = this%ncolbnd
747  naux = size(this%auxvar, 1)
748  !
749  ! -- dimension table
750  ntabrows = this%nbound
751  !
752  ! -- start building format statement to parse this%label, which
753  ! contains the column headers (except for boundname and auxnames)
754  ipos = index(this%listlabel, 'NO.')
755  if (ipos /= 0) then
756  write (cpos, '(i10)') ipos + 3
757  fmtlstbn = '(a'//trim(adjustl(cpos))
758  else
759  fmtlstbn = '(a7'
760  end if
761  ! -- sequence number, layer, row, and column.
762  if (size(this%dis%mshape) == 3) then
763  ntabcols = 4
764  fmtlstbn = trim(fmtlstbn)//',a7,a7,a7'
765  !
766  ! -- sequence number, layer, and cell2d.
767  else if (size(this%dis%mshape) == 2) then
768  ntabcols = 3
769  fmtlstbn = trim(fmtlstbn)//',a7,a7'
770  !
771  ! -- sequence number and node.
772  else
773  ntabcols = 2
774  fmtlstbn = trim(fmtlstbn)//',a7'
775  end if
776  !
777  ! -- Add fields for non-optional real values
778  ntabcols = ntabcols + ldim
779  do i = 1, ldim
780  fmtlstbn = trim(fmtlstbn)//',a16'
781  end do
782  !
783  ! -- Add field for boundary name
784  if (this%inamedbound == 1) then
785  ntabcols = ntabcols + 1
786  fmtlstbn = trim(fmtlstbn)//',a16'
787  end if
788  !
789  ! -- Add fields for auxiliary variables
790  ntabcols = ntabcols + naux
791  do i = 1, naux
792  fmtlstbn = trim(fmtlstbn)//',a16'
793  end do
794  fmtlstbn = trim(fmtlstbn)//')'
795  !
796  ! -- allocate words
797  allocate (words(ntabcols))
798  !
799  ! -- parse this%listlabel into words
800  read (this%listlabel, fmtlstbn) (words(i), i=1, ntabcols)
801  !
802  ! -- initialize the input table object
803  call table_cr(inputtab, ' ', ' ')
804  call inputtab%table_df(ntabrows, ntabcols, this%iout)
805  !
806  ! -- add the columns
807  ipos = 1
808  call inputtab%initialize_column(words(ipos), 10, alignment=tabcenter)
809  !
810  ! -- discretization
811  do i = 1, size(this%dis%mshape)
812  ipos = ipos + 1
813  call inputtab%initialize_column(words(ipos), 7, alignment=tabcenter)
814  end do
815  !
816  ! -- non-optional variables
817  do i = 1, ldim
818  ipos = ipos + 1
819  call inputtab%initialize_column(words(ipos), 16, alignment=tabcenter)
820  end do
821  !
822  ! -- boundname
823  if (this%inamedbound == 1) then
824  ipos = ipos + 1
825  tag = 'BOUNDNAME'
826  call inputtab%initialize_column(tag, lenboundname, alignment=tableft)
827  end if
828  !
829  ! -- aux variables
830  do i = 1, naux
831  call inputtab%initialize_column(this%auxname(i), 16, alignment=tabcenter)
832  end do
833  !
834  ! -- Write the table
835  do ii = 1, this%nbound
836  call inputtab%add_term(ii)
837  !
838  ! -- discretization
839  if (size(this%dis%mshape) == 3) then
840  nod = this%nodelist(ii)
841  call get_ijk(nod, this%dis%mshape(2), this%dis%mshape(3), &
842  this%dis%mshape(1), i, j, k)
843  call inputtab%add_term(k)
844  call inputtab%add_term(i)
845  call inputtab%add_term(j)
846  else if (size(this%dis%mshape) == 2) then
847  nod = this%nodelist(ii)
848  call get_ijk(nod, 1, this%dis%mshape(2), this%dis%mshape(1), i, j, k)
849  call inputtab%add_term(k)
850  call inputtab%add_term(j)
851  else
852  nod = this%nodelist(ii)
853  call inputtab%add_term(nod)
854  end if
855  !
856  ! -- non-optional variables
857  do jj = 1, ldim
858  call inputtab%add_term(this%bound_value(jj, ii))
859  end do
860  !
861  ! -- boundname
862  if (this%inamedbound == 1) then
863  call inputtab%add_term(this%boundname(ii))
864  end if
865  !
866  ! -- aux variables
867  do jj = 1, naux
868  call inputtab%add_term(this%auxvar(jj, ii))
869  end do
870  end do
871  !
872  ! -- deallocate the local variables
873  call inputtab%table_da()
874  deallocate (inputtab)
875  nullify (inputtab)
876  deallocate (words)
877  end subroutine write_list
878 
879  !> @ brief Return a bound value
880  !!
881  !! Return a bound value associated with an ncolbnd index
882  !! and row. This function should be overridden in the
883  !! derived package class.
884  !!
885  !<
886  function bound_value(this, col, row) result(bndval)
887  ! -- modules
888  use constantsmodule, only: dnodata
889  ! -- dummy variables
890  class(bndexttype), intent(inout) :: this !< BndExtType object
891  integer(I4B), intent(in) :: col
892  integer(I4B), intent(in) :: row
893  ! -- result
894  real(dp) :: bndval
895  !
896  ! -- override this return value by redefining this
897  ! routine in the derived package.
898  bndval = dnodata
899  end function bound_value
900 
901 end module bndextmodule
This module contains the extended boundary package.
subroutine bndext_rp(this)
subroutine write_list(this)
@ brief Log package list input
subroutine bndext_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate package arrays
subroutine bndext_df(this, neq, dis)
@ brief Define boundary package options and dimensions
subroutine bndext_da(this)
@ brief Deallocate package memory
subroutine default_nodelist(this)
Assign default nodelist when READASARRAYS is specified.
subroutine log_options(this, found, sfacauxname)
@ brief Log package options
subroutine cellid_to_nlist(this)
@ brief Update package nodelist
subroutine source_dimensions(this)
@ brief Source package dimensions from input context
subroutine check_cellid(this, ii, cellid, mshape, ndim)
@ brief Check for valid cellid
subroutine nodeu_to_nlist(this)
@ brief Update package nodelist
subroutine layarr_to_nlist(this)
Update the nodelist based on layer number variable input.
real(dp) function bound_value(this, col, row)
@ brief Return a bound value
subroutine source_options(this)
@ brief Source package options from input context
subroutine bndext_allocate_scalars(this)
@ brief Allocate package scalars
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:45
@ tabcenter
centered table column
Definition: Constants.f90:172
@ tableft
left justified table column
Definition: Constants.f90:171
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
integer(i4b), parameter lenauxname
maximum length of a aux variable
Definition: Constants.f90:35
integer(i4b), parameter lenboundname
maximum length of a bound name
Definition: Constants.f90:36
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
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:83
subroutine, public get_ijk(nodenumber, nrow, ncol, nlay, irow, icol, ilay)
Get row, column and layer indices from node number and grid dimensions. If nodenumber is invalid,...
Definition: GeomUtil.f90:100
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public ulstlb(iout, label, caux, ncaux, naux)
Print a label for a list.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module defines variable data types.
Definition: kind.f90:8
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public obs_cr(obs, inobs)
@ brief Create a new ObsType object
Definition: Obs.f90:225
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
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
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
subroutine, public table_cr(this, name, title)
Definition: Table.f90:87
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
subroutine, public tasmanager_cr(this, dis, modelname, iout)
Create the time-array series manager.
subroutine, public tsmanager_cr(this, iout, removeTsLinksOnCompletion, extendTsToEndOfSimulation)
Create the tsmanager.
@ 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