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

This module contains the extended boundary package. More...

Data Types

type  bndexttype
 @ brief BndExtType More...
 
type  bndextfoundtype
 @ brief BndExtFoundType More...
 

Functions/Subroutines

subroutine bndext_df (this, neq, dis)
 @ brief Define boundary package options and dimensions More...
 
subroutine bndext_rp (this)
 
subroutine bndext_rp_log (this)
 Write the input list to the listing file if requested. More...
 
subroutine bndext_da (this)
 @ brief Deallocate package memory More...
 
subroutine bndext_allocate_scalars (this)
 @ brief Allocate package scalars More...
 
subroutine bndext_allocate_arrays (this, nodelist, auxvar)
 @ brief Allocate package arrays More...
 
subroutine source_options (this)
 @ brief Source package options from input context More...
 
subroutine log_options (this, found, sfacauxname)
 @ brief Log package options More...
 
subroutine source_dimensions (this)
 @ brief Source package dimensions from input context More...
 
subroutine cellid_to_nlist (this)
 @ brief Update package nodelist More...
 
subroutine nodeu_to_nlist (this)
 @ brief Update package nodelist More...
 
subroutine layarr_to_nlist (this)
 Update the nodelist based on layer number variable input. More...
 
subroutine default_nodelist (this)
 Assign default nodelist when READASARRAYS is specified. More...
 
subroutine check_cellid (this, ii, cellid, mshape, ndim)
 @ brief Check for valid cellid More...
 
subroutine write_lstfile (this)
 @ brief Log package stress period input More...
 
real(dp) function bound_value (this, col, row)
 @ brief Return a bound value More...
 

Detailed Description

This module contains the extended boundary type that itself should be extended by model boundary packages that have been updated to source static and dynamic input data from the input context.

Function/Subroutine Documentation

◆ bndext_allocate_arrays()

subroutine bndextmodule::bndext_allocate_arrays ( class(bndexttype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)

Allocate and initialize base boundary package arrays. This method only needs to be overridden if additional arrays are defined for a specific package.

Parameters
thisBndExtType object
nodelistpackage nodelist
auxvarpackage aux variable array

Definition at line 254 of file BoundaryPackageExt.f90.

255  ! -- modules
257  ! -- dummy variables
258  class(BndExtType) :: this !< BndExtType object
259  ! -- local variables
260  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist !< package nodelist
261  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar !< package aux variable array
262  !
263  ! -- allocate base BndType arrays
264  call this%BndType%allocate_arrays(nodelist, auxvar)
265  !
266  ! -- set input context pointers
267  call mem_setptr(this%cellid, 'CELLID', this%input_mempath)
268  call mem_setptr(this%nodeulist, 'NODEULIST', this%input_mempath)
269  call mem_setptr(this%boundname_cst, 'BOUNDNAME', this%input_mempath)
270  !
271  ! -- checkin input context pointers
272  call mem_checkin(this%cellid, 'CELLID', this%memoryPath, &
273  'CELLID', this%input_mempath)
274  call mem_checkin(this%nodeulist, 'NODEULIST', this%memoryPath, &
275  'NODEULIST', this%input_mempath)
276  call mem_checkin(this%boundname_cst, lenboundname, 'BOUNDNAME_IDM', &
277  this%memoryPath, 'BOUNDNAME', this%input_mempath)
278  !
279  if (present(auxvar)) then
280  ! no-op
281  else
282  ! -- set auxvar input context pointer
283  call mem_setptr(this%auxvar, 'AUXVAR', this%input_mempath)
284  !
285  ! -- checkin auxvar input context pointer
286  call mem_checkin(this%auxvar, 'AUXVAR_IDM', this%memoryPath, &
287  'AUXVAR', this%input_mempath)
288  end if

◆ bndext_allocate_scalars()

subroutine bndextmodule::bndext_allocate_scalars ( class(bndexttype this)

Allocate and initialize base boundary package scalars. This method only needs to be overridden if additional scalars are defined for a specific package.

Parameters
thisBndExtType object

Definition at line 225 of file BoundaryPackageExt.f90.

226  ! -- modules
228  ! -- dummy variables
229  class(BndExtType) :: this !< BndExtType object
230  ! -- local variables
231  !
232  ! -- allocate base BndType scalars
233  call this%BndType%allocate_scalars()
234  !
235  ! -- set IPER pointer
236  call mem_setptr(this%iper, 'IPER', this%input_mempath)
237 
238  ! -- allocate internal scalars
239  allocate (this%readarraygrid)
240  allocate (this%readasarrays)
241 
242  ! -- initialize internal scalars
243  this%readarraygrid = .false.
244  this%readasarrays = .false.

◆ bndext_da()

subroutine bndextmodule::bndext_da ( class(bndexttype this)
private
Parameters
thisBndExtType object

Definition at line 191 of file BoundaryPackageExt.f90.

192  ! -- modules
194  ! -- dummy variables
195  class(BndExtType) :: this !< BndExtType object
196  !
197  ! -- deallocate checkin paths
198  call mem_deallocate(this%cellid, 'CELLID', this%memoryPath)
199  call mem_deallocate(this%nodeulist, 'NODEULIST', this%memoryPath)
200  call mem_deallocate(this%boundname_cst, 'BOUNDNAME_IDM', this%memoryPath)
201  call mem_deallocate(this%auxvar, 'AUXVAR_IDM', this%memoryPath)
202  !
203  ! -- reassign pointers for base class _da
204  call mem_setptr(this%boundname_cst, 'BOUNDNAME_CST', this%memoryPath)
205  call mem_setptr(this%auxvar, 'AUXVAR', this%memoryPath)
206  !
207  ! -- scalars
208  deallocate (this%readarraygrid)
209  deallocate (this%readasarrays)
210  nullify (this%readarraygrid)
211  nullify (this%readasarrays)
212  nullify (this%iper)
213  !
214  ! -- deallocate
215  call this%BndType%bnd_da()

◆ bndext_df()

subroutine bndextmodule::bndext_df ( class(bndexttype), intent(inout)  this,
integer(i4b), intent(inout)  neq,
class(disbasetype), pointer  dis 
)
private

Define base boundary package options and dimensions for a model boundary package.

Parameters
[in,out]thisBndExtType object
[in,out]neqnumber of equations
disdiscretization object

Definition at line 83 of file BoundaryPackageExt.f90.

84  ! -- modules
85  use basedismodule, only: disbasetype
89  ! -- dummy variables
90  class(BndExtType), intent(inout) :: this !< BndExtType object
91  integer(I4B), intent(inout) :: neq !< number of equations
92  class(DisBaseType), pointer :: dis !< discretization object
93  !
94  ! -- set pointer to dis object for the model
95  this%dis => dis
96  !
97  ! -- Create time series managers
98  ! -- Not in use by this type but BndType uses and deallocates
99  call tsmanager_cr(this%TsManager, this%iout)
100  call tasmanager_cr(this%TasManager, dis, this%name_model, this%iout)
101  !
102  ! -- create obs package
103  call obs_cr(this%obs, this%inobspkg)
104  !
105  ! -- Write information to model list file
106  write (this%iout, 1) trim(this%filtyp), trim(adjustl(this%text)), &
107  this%input_mempath
108 1 format(1x, /1x, a, ' -- ', a, ' PACKAGE, VERSION 8, 2/22/2014', &
109  ' INPUT READ FROM MEMPATH: ', a)
110  !
111  ! -- source options
112  call this%source_options()
113  !
114  ! -- Define time series managers
115  call this%tsmanager%tsmanager_df()
116  call this%tasmanager%tasmanager_df()
117  !
118  ! -- source dimensions
119  call this%source_dimensions()
120  !
121  ! -- update package moffset for packages that add rows
122  if (this%npakeq > 0) then
123  this%ioffset = neq - this%dis%nodes
124  end if
125  !
126  ! -- update neq
127  neq = neq + this%npakeq
128  !
129  ! -- Store information needed for observations
130  if (this%bnd_obs_supported()) then
131  call this%obs%obs_df(this%iout, this%packName, this%filtyp, this%dis)
132  call this%bnd_df_obs()
133  end if
134  !
135  ! -- Call define_listlabel to construct the list label that is written
136  ! when PRINT_INPUT option is used.
137  call this%define_listlabel()
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.
Here is the call graph for this function:

◆ bndext_rp()

subroutine bndextmodule::bndext_rp ( class(bndexttype), intent(inout)  this)
Parameters
[in,out]thisBndExtType object

Definition at line 140 of file BoundaryPackageExt.f90.

141  ! -- modules
142  use tdismodule, only: kper
145  ! -- dummy variables
146  class(BndExtType), intent(inout) :: this !< BndExtType object
147  ! -- local variables
148  logical(LGP) :: found
149  integer(I4B) :: n
150  !
151  if (this%iper /= kper) return
152 
153  if (.not. this%readasarrays) then
154  ! -- copy nbound from input context
155  call mem_set_value(this%nbound, 'NBOUND', this%input_mempath, &
156  found, release=.false.)
157  end if
158 
159  if (this%readarraygrid) then
160  call this%nodeu_to_nlist()
161  else if (this%readasarrays) then
162  call this%layarr_to_nlist()
163  else
164  call this%cellid_to_nlist()
165  !
166  ! -- update boundname string list
167  if (this%inamedbound /= 0) then
168  do n = 1, size(this%boundname_cst)
169  this%boundname(n) = this%boundname_cst(n)
170  end do
171  end if
172  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:26

◆ bndext_rp_log()

subroutine bndextmodule::bndext_rp_log ( class(bndexttype), intent(inout)  this)

Called from model control files after bnd_rp(), which ensures bound_value() dispatches to the correct derived type.

Definition at line 180 of file BoundaryPackageExt.f90.

181  ! -- dummy
182  class(BndExtType), intent(inout) :: this
183  !
184  if (this%iprpak /= 0) then
185  call this%write_lstfile()
186  end if

◆ bound_value()

real(dp) function bndextmodule::bound_value ( class(bndexttype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)

Return a bound value associated with an ncolbnd index and row. This function should be overridden in the derived package class.

Parameters
[in,out]thisBndExtType object

Definition at line 902 of file BoundaryPackageExt.f90.

903  ! -- modules
904  use constantsmodule, only: dnodata
905  ! -- dummy variables
906  class(BndExtType), intent(inout) :: this !< BndExtType object
907  integer(I4B), intent(in) :: col
908  integer(I4B), intent(in) :: row
909  ! -- result
910  real(DP) :: bndval
911  !
912  ! -- override this return value by redefining this
913  ! routine in the derived package.
914  bndval = dnodata
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95

◆ cellid_to_nlist()

subroutine bndextmodule::cellid_to_nlist ( class(bndexttype this)

Convert period updated cellids to node numbers.

Parameters
thisBndExtType object

Definition at line 506 of file BoundaryPackageExt.f90.

507  ! -- modules
508  use simvariablesmodule, only: errmsg
509  ! -- dummy
510  class(BndExtType) :: this !< BndExtType object
511  ! -- local
512  integer(I4B), dimension(:), pointer :: cellid
513  integer(I4B) :: n, nodeu, noder
514  character(len=LINELENGTH) :: nodestr
515  !
516  ! -- update nodelist
517  do n = 1, this%nbound
518  !
519  ! -- set cellid
520  cellid => this%cellid(:, n)
521  !
522  ! -- ensure cellid is valid, store an error otherwise
523  call this%check_cellid(n, cellid, this%dis%mshape, this%dis%ndim)
524  !
525  ! -- Determine user node number
526  if (this%dis%ndim == 1) then
527  nodeu = cellid(1)
528  elseif (this%dis%ndim == 2) then
529  nodeu = get_node(cellid(1), 1, cellid(2), &
530  this%dis%mshape(1), 1, &
531  this%dis%mshape(2))
532  else
533  nodeu = get_node(cellid(1), cellid(2), cellid(3), &
534  this%dis%mshape(1), &
535  this%dis%mshape(2), &
536  this%dis%mshape(3))
537  end if
538  !
539  ! -- update the nodelist
540  if (this%dis%nodes < this%dis%nodesuser) then
541  ! -- convert user to reduced node numbers
542  noder = this%dis%get_nodenumber(nodeu, 0)
543  if (noder <= 0) then
544  call this%dis%nodeu_to_string(nodeu, nodestr)
545  write (errmsg, *) &
546  ' Cell is outside active grid domain: '// &
547  trim(adjustl(nodestr))
548  call store_error(errmsg)
549  end if
550  this%nodelist(n) = noder
551  else
552  this%nodelist(n) = nodeu
553  end if
554  end do
555  !
556  ! -- exit if errors were found
557  if (count_errors() > 0) then
558  write (errmsg, *) count_errors(), ' errors encountered.'
559  call store_error(errmsg)
560  call store_error_filename(this%input_fname)
561  end if
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
Here is the call graph for this function:

◆ check_cellid()

subroutine bndextmodule::check_cellid ( class(bndexttype this,
integer(i4b), intent(in)  ii,
integer(i4b), dimension(:), intent(in)  cellid,
integer(i4b), dimension(:), intent(in)  mshape,
integer(i4b), intent(in)  ndim 
)
private
Parameters
thisBndExtType object
[in]mshapemodel shape
[in]ndimsize of mshape

Definition at line 678 of file BoundaryPackageExt.f90.

679  ! -- modules
680  use simvariablesmodule, only: errmsg
681  ! -- dummy
682  class(BndExtType) :: this !< BndExtType object
683  ! -- local
684  integer(I4B), intent(in) :: ii
685  integer(I4B), dimension(:), intent(in) :: cellid !< cellid
686  integer(I4B), dimension(:), intent(in) :: mshape !< model shape
687  integer(I4B), intent(in) :: ndim !< size of mshape
688  character(len=20) :: cellstr, mshstr
689  character(len=*), parameter :: fmterr = &
690  "('List entry ',i0,' contains cellid ',a,' but this cellid is invalid &
691  &for model with shape ', a)"
692  character(len=*), parameter :: fmtndim1 = &
693  "('(',i0,')')"
694  character(len=*), parameter :: fmtndim2 = &
695  "('(',i0,',',i0,')')"
696  character(len=*), parameter :: fmtndim3 = &
697  "('(',i0,',',i0,',',i0,')')"
698  select case (ndim)
699  case (1)
700  !
701  if (cellid(1) < 1 .or. cellid(1) > mshape(1)) then
702  write (cellstr, fmtndim1) cellid(1)
703  write (mshstr, fmtndim1) mshape(1)
704  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
705  call store_error(errmsg)
706  end if
707  !
708  case (2)
709  !
710  if (cellid(1) < 1 .or. cellid(1) > mshape(1) .or. &
711  cellid(2) < 1 .or. cellid(2) > mshape(2)) then
712  write (cellstr, fmtndim2) cellid(1), cellid(2)
713  write (mshstr, fmtndim2) mshape(1), mshape(2)
714  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
715  call store_error(errmsg)
716  end if
717  !
718  case (3)
719  !
720  if (cellid(1) < 1 .or. cellid(1) > mshape(1) .or. &
721  cellid(2) < 1 .or. cellid(2) > mshape(2) .or. &
722  cellid(3) < 1 .or. cellid(3) > mshape(3)) then
723  write (cellstr, fmtndim3) cellid(1), cellid(2), cellid(3)
724  write (mshstr, fmtndim3) mshape(1), mshape(2), mshape(3)
725  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
726  call store_error(errmsg)
727  end if
728  !
729  case default
730  end select
Here is the call graph for this function:

◆ default_nodelist()

subroutine bndextmodule::default_nodelist ( class(bndexttype this)

Equivalent to reading layer number array as CONSTANT 1

Definition at line 640 of file BoundaryPackageExt.f90.

641  ! -- dummy
642  class(BndExtType) :: this
643  ! -- local
644  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nodeu, noder, ipos
645  !
646  if (this%readasarrays) then
647  !
648  ! -- set variables
649  if (this%dis%ndim == 3) then
650  nlay = this%dis%mshape(1)
651  nrow = this%dis%mshape(2)
652  ncol = this%dis%mshape(3)
653  elseif (this%dis%ndim == 2) then
654  nlay = this%dis%mshape(1)
655  nrow = 1
656  ncol = this%dis%mshape(2)
657  end if
658  !
659  ! -- Populate nodelist
660  ipos = 1
661  il = 1
662  do ir = 1, nrow
663  do ic = 1, ncol
664  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
665  noder = this%dis%get_nodenumber(nodeu, 0)
666  this%nodelist(ipos) = noder
667  ipos = ipos + 1
668  end do
669  end do
670  !
671  ! -- Assign nbound
672  this%nbound = ipos - 1
673  end if
Here is the call graph for this function:

◆ layarr_to_nlist()

subroutine bndextmodule::layarr_to_nlist ( class(bndexttype this)

This is a module scoped routine to check for I<filtyp> input. If array input was provided, INI<filtyp> and I<filtyp> will be allocated in the input context. If the read state variable INI<filtyp> is set to 1 during this period update, I<filtyp> input was read and is used here to update the nodelist.

Parameters
thisBndExtType object

Definition at line 603 of file BoundaryPackageExt.f90.

604  ! -- modules
606  use constantsmodule, only: lenvarname
607  ! -- dummy
608  class(BndExtType) :: this !< BndExtType object
609  character(len=LENVARNAME) :: ilayname, inilayname
610  character(len=24) :: aname = ' LAYER OR NODE INDEX'
611  ! -- local
612  integer(I4B), dimension(:), contiguous, &
613  pointer :: ilay => null()
614  integer(I4B), pointer :: inilay => null()
615  !
616  ! set ilay and read state variable names
617  ilayname = 'I'//trim(this%filtyp)
618  inilayname = 'INI'//trim(this%filtyp)
619  !
620  ! -- set pointer to input context read state variable
621  call mem_setptr(inilay, inilayname, this%input_mempath)
622  !
623  ! -- check read state
624  if (inilay == 1) then
625  ! -- ilay variable was read this period
626  !
627  ! -- set pointer to input context layer index variable
628  call mem_setptr(ilay, ilayname, this%input_mempath)
629  !
630  ! -- update nodelist
631  call this%dis%nlarray_to_nodelist(ilay, this%nodelist, this%maxbound, &
632  this%nbound, aname)
633  end if
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17

◆ log_options()

subroutine bndextmodule::log_options ( class(bndexttype), intent(inout)  this,
type(bndextfoundtype), intent(in)  found,
character(len=*), intent(in)  sfacauxname 
)
Parameters
[in,out]thisBndExtType object

Definition at line 405 of file BoundaryPackageExt.f90.

406  ! -- modules
407  ! -- dummy variables
408  class(BndExtType), intent(inout) :: this !< BndExtType object
409  type(BndExtFoundType), intent(in) :: found
410  character(len=*), intent(in) :: sfacauxname
411  ! -- local variables
412  ! -- format
413  character(len=*), parameter :: fmtreadasarrays = &
414  &"(4x, 'PACKAGE INPUT WILL BE READ AS LAYER ARRAYS.')"
415  character(len=*), parameter :: fmtreadarraygrid = &
416  &"(4x, 'PACKAGE INPUT WILL BE READ AS GRID ARRAYS.')"
417  character(len=*), parameter :: fmtflow = &
418  &"(4x, 'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL')"
419  !
420  ! -- log found options
421  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
422  //' BASE OPTIONS'
423  !
424  if (this%readasarrays) then
425  write (this%iout, fmtreadasarrays)
426  end if
427  !
428  if (this%readarraygrid) then
429  write (this%iout, fmtreadarraygrid)
430  end if
431  !
432  if (found%ipakcb) then
433  write (this%iout, fmtflow)
434  end if
435  !
436  if (found%iprpak) then
437  write (this%iout, '(4x,a)') &
438  'LISTS OF '//trim(adjustl(this%text))//' CELLS WILL BE PRINTED.'
439  end if
440  !
441  if (found%iprflow) then
442  write (this%iout, '(4x,a)') trim(adjustl(this%text))// &
443  ' FLOWS WILL BE PRINTED TO LISTING FILE.'
444  end if
445  !
446  if (found%boundnames) then
447  write (this%iout, '(4x,a)') trim(adjustl(this%text))// &
448  ' BOUNDARIES HAVE NAMES IN LAST COLUMN.'
449  end if
450  !
451  if (found%auxmultname) then
452  write (this%iout, '(4x,a,a)') &
453  'AUXILIARY MULTIPLIER NAME: ', sfacauxname
454  end if
455  !
456  if (found%inewton) then
457  write (this%iout, '(4x,a)') &
458  'NEWTON-RAPHSON method disabled for unconfined cells'
459  end if
460  !
461  ! -- close logging block
462  write (this%iout, '(1x,a)') &
463  'END OF '//trim(adjustl(this%text))//' BASE OPTIONS'

◆ nodeu_to_nlist()

subroutine bndextmodule::nodeu_to_nlist ( class(bndexttype this)

Convert period user nodes to reduced nodes

Parameters
thisBndExtType object

Definition at line 569 of file BoundaryPackageExt.f90.

570  ! -- modules
572  ! -- dummy
573  class(BndExtType) :: this !< BndExtType object
574  integer(I4B) :: n, noder, nodeuser, ninactive
575 
576  ninactive = 0
577 
578  ! -- Set the nodelist
579  do n = 1, this%nbound
580  nodeuser = this%nodeulist(n)
581  noder = this%dis%get_nodenumber(nodeuser, 0)
582  if (noder > 0) then
583  this%nodelist(n) = noder
584  else
585  ninactive = ninactive + 1
586  end if
587  end do
588 
589  ! update nbound
590  this%nbound = this%nbound - ninactive

◆ source_dimensions()

subroutine bndextmodule::source_dimensions ( class(bndexttype), intent(inout)  this)
private
Parameters
[in,out]thisBndExtType object

Definition at line 468 of file BoundaryPackageExt.f90.

470  ! -- dummy variables
471  class(BndExtType), intent(inout) :: this !< BndExtType object
472  ! -- local variables
473  type(BndExtFoundType) :: found
474  !
475  if (this%readasarrays) then
476  this%maxbound = this%dis%get_ncpl()
477  else
478  ! -- open dimensions logging block
479  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text))// &
480  ' BASE DIMENSIONS'
481 
482  ! -- update defaults with idm sourced values
483  call mem_set_value(this%maxbound, 'MAXBOUND', this%input_mempath, &
484  found%maxbound, release=.false.)
485 
486  write (this%iout, '(4x,a,i7)') 'MAXBOUND = ', this%maxbound
487 
488  ! -- close logging block
489  write (this%iout, '(1x,a)') &
490  'END OF '//trim(adjustl(this%text))//' BASE DIMENSIONS'
491  end if
492  !
493  ! -- verify dimensions were set
494  if (this%maxbound <= 0) then
495  write (errmsg, '(a)') 'MAXBOUND must be an integer greater than zero.'
496  call store_error(errmsg)
497  call store_error_filename(this%input_fname)
498  end if
Here is the call graph for this function:

◆ source_options()

subroutine bndextmodule::source_options ( class(bndexttype), intent(inout)  this)
Parameters
[in,out]thisBndExtType object

Definition at line 293 of file BoundaryPackageExt.f90.

294  ! -- modules
295  use memorymanagermodule, only: mem_reallocate, mem_setptr !, get_isize
300  ! -- dummy variables
301  class(BndExtType), intent(inout) :: this !< BndExtType object
302  ! -- local variables
303  type(BndExtFoundType) :: found
304  logical(LGP) :: found_readarr
305  character(len=LENAUXNAME) :: sfacauxname
306  integer(I4B) :: n
307  !
308  ! -- update defaults with idm sourced values
309  call mem_set_value(this%naux, 'NAUX', this%input_mempath, found%naux, &
310  release=.false.)
311  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
312  call mem_set_value(this%iprpak, 'IPRPAK', this%input_mempath, found%iprpak)
313  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
314  call mem_set_value(this%inamedbound, 'BOUNDNAMES', this%input_mempath, &
315  found%boundnames, release=.false.)
316  call mem_set_value(sfacauxname, 'AUXMULTNAME', this%input_mempath, &
317  found%auxmultname)
318  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
319  call mem_set_value(this%readarraygrid, 'READARRAYGRID', this%input_mempath, &
320  found_readarr)
321  call mem_set_value(this%readasarrays, 'READASARRAYS', this%input_mempath, &
322  found_readarr)
323  !
324  ! -- log found options
325  call this%log_options(found, sfacauxname)
326  !
327  ! -- reallocate aux arrays if aux variables provided
328  if (found%naux .and. this%naux > 0) then
329  call mem_reallocate(this%auxname, lenauxname, this%naux, &
330  'AUXNAME', this%memoryPath)
331  call mem_reallocate(this%auxname_cst, lenauxname, this%naux, &
332  'AUXNAME_CST', this%memoryPath)
333  call mem_set_value(this%auxname_cst, 'AUXILIARY', this%input_mempath, &
334  found%auxiliary, release=.false.)
335  !
336  do n = 1, this%naux
337  this%auxname(n) = this%auxname_cst(n)
338  end do
339  end if
340  !
341  ! -- save flows option active
342  if (found%ipakcb) this%ipakcb = -1
343  !
344  ! -- auxmultname provided
345  if (found%auxmultname) this%iauxmultcol = -1
346  !
347  !
348  ! -- enforce 0 or 1 OBS6_FILENAME entries in option block
349  if (filein_fname(this%obs%inputFilename, 'OBS6_FILENAME', &
350  this%input_mempath, this%input_fname)) then
351  this%obs%active = .true.
352  this%obs%inUnitObs = getunit()
353  call openfile(this%obs%inUnitObs, this%iout, this%obs%inputFilename, 'OBS')
354  end if
355  !
356  ! -- no newton specified
357  if (found%inewton) this%inewton = 0
358  !
359  ! -- AUXMULTNAME was specified, so find column of auxvar that will be multiplier
360  if (this%iauxmultcol < 0) then
361  !
362  ! -- Error if no aux variable specified
363  if (this%naux == 0) then
364  write (errmsg, '(a,2(1x,a))') &
365  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
366  'but no AUX variables specified.'
367  call store_error(errmsg)
368  end if
369  !
370  ! -- Assign mult column
371  this%iauxmultcol = 0
372  do n = 1, this%naux
373  if (sfacauxname == this%auxname(n)) then
374  this%iauxmultcol = n
375  exit
376  end if
377  end do
378  !
379  ! -- Error if aux variable cannot be found
380  if (this%iauxmultcol == 0) then
381  write (errmsg, '(a,2(1x,a))') &
382  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
383  'but no AUX variable found with this name.'
384  call store_error(errmsg)
385  end if
386  end if
387  !
388  if (this%readasarrays) then
389  if (.not. this%dis%supports_layers()) then
390  errmsg = 'READASARRAYS option is not compatible with selected'// &
391  ' discretization type.'
392  call store_error(errmsg)
393  call store_error_filename(this%input_fname)
394  end if
395  end if
396  !
397  ! -- terminate if errors were detected
398  if (count_errors() > 0) then
399  call store_error_filename(this%input_fname)
400  end if
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
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
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function:

◆ write_lstfile()

subroutine bndextmodule::write_lstfile ( class(bndexttype this)

Log period based input. This routine requires a package specific bound_value() routine to report accurate bound values.

Parameters
thisBndExtType object

Definition at line 739 of file BoundaryPackageExt.f90.

740  ! -- modules
743  use inputoutputmodule, only: ulstlb
744  use tablemodule, only: tabletype, table_cr
745  ! -- dummy
746  class(BndExtType) :: this !< BndExtType object
747  ! -- local
748  character(len=10) :: cpos
749  character(len=LINELENGTH) :: tag
750  character(len=LINELENGTH), allocatable, dimension(:) :: words
751  integer(I4B) :: ntabrows
752  integer(I4B) :: ntabcols
753  integer(I4B) :: ipos
754  integer(I4B) :: ii, jj, i, j, k, nod
755  integer(I4B) :: ldim
756  integer(I4B) :: naux
757  type(TableType), pointer :: inputtab => null()
758  ! -- formats
759  character(len=LINELENGTH) :: fmtlstbn
760  !
761  ! -- Determine sizes
762  ldim = this%ncolbnd
763  naux = size(this%auxvar, 1)
764  !
765  ! -- dimension table
766  ntabrows = this%nbound
767  !
768  ! -- start building format statement to parse this%label, which
769  ! contains the column headers (except for boundname and auxnames)
770  ipos = index(this%listlabel, 'NO.')
771  if (ipos /= 0) then
772  write (cpos, '(i10)') ipos + 3
773  fmtlstbn = '(a'//trim(adjustl(cpos))
774  else
775  fmtlstbn = '(a7'
776  end if
777  ! -- sequence number, layer, row, and column.
778  if (size(this%dis%mshape) == 3) then
779  ntabcols = 4
780  fmtlstbn = trim(fmtlstbn)//',a7,a7,a7'
781  !
782  ! -- sequence number, layer, and cell2d.
783  else if (size(this%dis%mshape) == 2) then
784  ntabcols = 3
785  fmtlstbn = trim(fmtlstbn)//',a7,a7'
786  !
787  ! -- sequence number and node.
788  else
789  ntabcols = 2
790  fmtlstbn = trim(fmtlstbn)//',a7'
791  end if
792  !
793  ! -- Add fields for non-optional real values
794  ntabcols = ntabcols + ldim
795  do i = 1, ldim
796  fmtlstbn = trim(fmtlstbn)//',a16'
797  end do
798  !
799  ! -- Add field for boundary name
800  if (this%inamedbound == 1) then
801  ntabcols = ntabcols + 1
802  fmtlstbn = trim(fmtlstbn)//',a16'
803  end if
804  !
805  ! -- Add fields for auxiliary variables
806  ntabcols = ntabcols + naux
807  do i = 1, naux
808  fmtlstbn = trim(fmtlstbn)//',a16'
809  end do
810  fmtlstbn = trim(fmtlstbn)//')'
811  !
812  ! -- allocate words
813  allocate (words(ntabcols))
814  !
815  ! -- parse this%listlabel into words
816  read (this%listlabel, fmtlstbn) (words(i), i=1, ntabcols)
817  !
818  ! -- initialize the input table object
819  call table_cr(inputtab, ' ', ' ')
820  call inputtab%table_df(ntabrows, ntabcols, this%iout)
821  !
822  ! -- add the columns
823  ipos = 1
824  call inputtab%initialize_column(words(ipos), 10, alignment=tabcenter)
825  !
826  ! -- discretization
827  do i = 1, size(this%dis%mshape)
828  ipos = ipos + 1
829  call inputtab%initialize_column(words(ipos), 7, alignment=tabcenter)
830  end do
831  !
832  ! -- non-optional variables
833  do i = 1, ldim
834  ipos = ipos + 1
835  call inputtab%initialize_column(words(ipos), 16, alignment=tabcenter)
836  end do
837  !
838  ! -- boundname
839  if (this%inamedbound == 1) then
840  ipos = ipos + 1
841  tag = 'BOUNDNAME'
842  call inputtab%initialize_column(tag, lenboundname, alignment=tableft)
843  end if
844  !
845  ! -- aux variables
846  do i = 1, naux
847  call inputtab%initialize_column(this%auxname(i), 16, alignment=tabcenter)
848  end do
849  !
850  ! -- Write the table
851  do ii = 1, this%nbound
852  call inputtab%add_term(ii)
853  !
854  ! -- discretization
855  if (size(this%dis%mshape) == 3) then
856  nod = this%nodelist(ii)
857  call get_ijk(nod, this%dis%mshape(2), this%dis%mshape(3), &
858  this%dis%mshape(1), i, j, k)
859  call inputtab%add_term(k)
860  call inputtab%add_term(i)
861  call inputtab%add_term(j)
862  else if (size(this%dis%mshape) == 2) then
863  nod = this%nodelist(ii)
864  call get_ijk(nod, 1, this%dis%mshape(2), this%dis%mshape(1), i, j, k)
865  call inputtab%add_term(k)
866  call inputtab%add_term(j)
867  else
868  nod = this%nodelist(ii)
869  call inputtab%add_term(nod)
870  end if
871  !
872  ! -- non-optional variables
873  do jj = 1, ldim
874  call inputtab%add_term(this%bound_value(jj, ii))
875  end do
876  !
877  ! -- boundname
878  if (this%inamedbound == 1) then
879  call inputtab%add_term(this%boundname(ii))
880  end if
881  !
882  ! -- aux variables
883  do jj = 1, naux
884  call inputtab%add_term(this%auxvar(jj, ii))
885  end do
886  end do
887  !
888  ! -- deallocate the local variables
889  call inputtab%table_da()
890  deallocate (inputtab)
891  nullify (inputtab)
892  deallocate (words)
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
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
subroutine, public ulstlb(iout, label, caux, ncaux, naux)
Print a label for a list.
subroutine, public table_cr(this, name, title)
Definition: Table.f90:87
Here is the call graph for this function: