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)
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:23

◆ 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 901 of file BoundaryPackageExt.f90.

902  ! -- modules
903  use constantsmodule, only: dnodata
904  ! -- dummy variables
905  class(BndExtType), intent(inout) :: this !< BndExtType object
906  integer(I4B), intent(in) :: col
907  integer(I4B), intent(in) :: row
908  ! -- result
909  real(DP) :: bndval
910  !
911  ! -- override this return value by redefining this
912  ! routine in the derived package.
913  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 505 of file BoundaryPackageExt.f90.

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

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

640  ! -- dummy
641  class(BndExtType) :: this
642  ! -- local
643  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nodeu, noder, ipos
644  !
645  if (this%readasarrays) then
646  !
647  ! -- set variables
648  if (this%dis%ndim == 3) then
649  nlay = this%dis%mshape(1)
650  nrow = this%dis%mshape(2)
651  ncol = this%dis%mshape(3)
652  elseif (this%dis%ndim == 2) then
653  nlay = this%dis%mshape(1)
654  nrow = 1
655  ncol = this%dis%mshape(2)
656  end if
657  !
658  ! -- Populate nodelist
659  ipos = 1
660  il = 1
661  do ir = 1, nrow
662  do ic = 1, ncol
663  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
664  noder = this%dis%get_nodenumber(nodeu, 0)
665  this%nodelist(ipos) = noder
666  ipos = ipos + 1
667  end do
668  end do
669  !
670  ! -- Assign nbound
671  this%nbound = ipos - 1
672  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 602 of file BoundaryPackageExt.f90.

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

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

569  ! -- modules
571  ! -- dummy
572  class(BndExtType) :: this !< BndExtType object
573  integer(I4B) :: n, noder, nodeuser, ninactive
574 
575  ninactive = 0
576 
577  ! -- Set the nodelist
578  do n = 1, this%nbound
579  nodeuser = this%nodeulist(n)
580  noder = this%dis%get_nodenumber(nodeuser, 0)
581  if (noder > 0) then
582  this%nodelist(n) = noder
583  else
584  ninactive = ninactive + 1
585  end if
586  end do
587 
588  ! update nbound
589  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 467 of file BoundaryPackageExt.f90.

469  ! -- dummy variables
470  class(BndExtType), intent(inout) :: this !< BndExtType object
471  ! -- local variables
472  type(BndExtFoundType) :: found
473  !
474  if (this%readasarrays) then
475  this%maxbound = this%dis%get_ncpl()
476  else
477  ! -- open dimensions logging block
478  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text))// &
479  ' BASE DIMENSIONS'
480 
481  ! -- update defaults with idm sourced values
482  call mem_set_value(this%maxbound, 'MAXBOUND', this%input_mempath, &
483  found%maxbound)
484 
485  write (this%iout, '(4x,a,i7)') 'MAXBOUND = ', this%maxbound
486 
487  ! -- close logging block
488  write (this%iout, '(1x,a)') &
489  'END OF '//trim(adjustl(this%text))//' BASE DIMENSIONS'
490  end if
491  !
492  ! -- verify dimensions were set
493  if (this%maxbound <= 0) then
494  write (errmsg, '(a)') 'MAXBOUND must be an integer greater than zero.'
495  call store_error(errmsg)
496  call store_error_filename(this%input_fname)
497  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  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
311  call mem_set_value(this%iprpak, 'IPRPAK', this%input_mempath, found%iprpak)
312  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
313  call mem_set_value(this%inamedbound, 'BOUNDNAMES', this%input_mempath, &
314  found%boundnames)
315  call mem_set_value(sfacauxname, 'AUXMULTNAME', this%input_mempath, &
316  found%auxmultname)
317  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
318  call mem_set_value(this%readarraygrid, 'READARRAYGRID', this%input_mempath, &
319  found_readarr)
320  call mem_set_value(this%readasarrays, 'READASARRAYS', this%input_mempath, &
321  found_readarr)
322  !
323  ! -- log found options
324  call this%log_options(found, sfacauxname)
325  !
326  ! -- reallocate aux arrays if aux variables provided
327  if (found%naux .and. this%naux > 0) then
328  call mem_reallocate(this%auxname, lenauxname, this%naux, &
329  'AUXNAME', this%memoryPath)
330  call mem_reallocate(this%auxname_cst, lenauxname, this%naux, &
331  'AUXNAME_CST', this%memoryPath)
332  call mem_set_value(this%auxname_cst, 'AUXILIARY', this%input_mempath, &
333  found%auxiliary)
334  !
335  do n = 1, this%naux
336  this%auxname(n) = this%auxname_cst(n)
337  end do
338  end if
339  !
340  ! -- save flows option active
341  if (found%ipakcb) this%ipakcb = -1
342  !
343  ! -- auxmultname provided
344  if (found%auxmultname) this%iauxmultcol = -1
345  !
346  !
347  ! -- enforce 0 or 1 OBS6_FILENAME entries in option block
348  if (filein_fname(this%obs%inputFilename, 'OBS6_FILENAME', &
349  this%input_mempath, this%input_fname)) then
350  this%obs%active = .true.
351  this%obs%inUnitObs = getunit()
352  call openfile(this%obs%inUnitObs, this%iout, this%obs%inputFilename, 'OBS')
353  end if
354  !
355  ! -- no newton specified
356  if (found%inewton) this%inewton = 0
357  !
358  ! -- AUXMULTNAME was specified, so find column of auxvar that will be multiplier
359  if (this%iauxmultcol < 0) then
360  !
361  ! -- Error if no aux variable specified
362  if (this%naux == 0) then
363  write (errmsg, '(a,2(1x,a))') &
364  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
365  'but no AUX variables specified.'
366  call store_error(errmsg)
367  end if
368  !
369  ! -- Assign mult column
370  this%iauxmultcol = 0
371  do n = 1, this%naux
372  if (sfacauxname == this%auxname(n)) then
373  this%iauxmultcol = n
374  exit
375  end if
376  end do
377  !
378  ! -- Error if aux variable cannot be found
379  if (this%iauxmultcol == 0) then
380  write (errmsg, '(a,2(1x,a))') &
381  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
382  'but no AUX variable found with this name.'
383  call store_error(errmsg)
384  end if
385  end if
386  !
387  if (this%readasarrays) then
388  if (.not. this%dis%supports_layers()) then
389  errmsg = 'READASARRAYS option is not compatible with selected'// &
390  ' discretization type.'
391  call store_error(errmsg)
392  call store_error_filename(this%input_fname)
393  end if
394  end if
395  !
396  ! -- terminate if errors were detected
397  if (count_errors() > 0) then
398  call store_error_filename(this%input_fname)
399  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 738 of file BoundaryPackageExt.f90.

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