MODFLOW 6  version 6.6.0.dev0
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  ! -- arrays
35  integer(I4B), dimension(:, :), pointer, contiguous :: cellid => null()
36  contains
37  procedure :: bnd_df => bndext_df
38  procedure :: bnd_rp => bndext_rp
39  procedure :: bnd_da => bndext_da
40  procedure :: allocate_scalars => bndext_allocate_scalars
41  procedure :: allocate_arrays => bndext_allocate_arrays
42  procedure :: source_options
43  procedure :: source_dimensions
44  procedure :: log_options
45  procedure :: nodelist_update
46  procedure :: check_cellid
47  procedure :: write_list
48  procedure :: bound_value
49  end type bndexttype
50 
51  !> @ brief BndExtFoundType
52  !!
53  !! This type is used to simplify the tracking of common parameters
54  !! that are sourced from the input context.
55  !<
57  logical :: naux = .false.
58  logical :: ipakcb = .false.
59  logical :: iprpak = .false.
60  logical :: iprflow = .false.
61  logical :: boundnames = .false.
62  logical :: auxmultname = .false.
63  logical :: inewton = .false.
64  logical :: auxiliary = .false.
65  logical :: maxbound = .false.
66  end type bndextfoundtype
67 
68 contains
69 
70  !> @ brief Define boundary package options and dimensions
71  !!
72  !! Define base boundary package options and dimensions for
73  !! a model boundary package.
74  !!
75  !<
76  subroutine bndext_df(this, neq, dis)
77  ! -- modules
78  use basedismodule, only: disbasetype
82  ! -- dummy variables
83  class(bndexttype), intent(inout) :: this !< BndExtType object
84  integer(I4B), intent(inout) :: neq !< number of equations
85  class(disbasetype), pointer :: dis !< discretization object
86  !
87  ! -- set pointer to dis object for the model
88  this%dis => dis
89  !
90  ! -- Create time series managers
91  ! -- Not in use by this type but BndType uses and deallocates
92  call tsmanager_cr(this%TsManager, this%iout)
93  call tasmanager_cr(this%TasManager, dis, this%name_model, this%iout)
94  !
95  ! -- create obs package
96  call obs_cr(this%obs, this%inobspkg)
97  !
98  ! -- Write information to model list file
99  write (this%iout, 1) this%filtyp, trim(adjustl(this%text)), this%input_mempath
100 1 format(1x, /1x, a, ' -- ', a, ' PACKAGE, VERSION 8, 2/22/2014', &
101  ' INPUT READ FROM MEMPATH: ', a)
102  !
103  ! -- source options
104  call this%source_options()
105  !
106  ! -- Define time series managers
107  call this%tsmanager%tsmanager_df()
108  call this%tasmanager%tasmanager_df()
109  !
110  ! -- source dimensions
111  call this%source_dimensions()
112  !
113  ! -- update package moffset for packages that add rows
114  if (this%npakeq > 0) then
115  this%ioffset = neq - this%dis%nodes
116  end if
117  !
118  ! -- update neq
119  neq = neq + this%npakeq
120  !
121  ! -- Store information needed for observations
122  if (this%bnd_obs_supported()) then
123  call this%obs%obs_df(this%iout, this%packName, this%filtyp, this%dis)
124  call this%bnd_df_obs()
125  end if
126  end subroutine bndext_df
127 
128  subroutine bndext_rp(this)
129  ! -- modules
130  use tdismodule, only: kper
133  ! -- dummy variables
134  class(bndexttype), intent(inout) :: this !< BndExtType object
135  ! -- local variables
136  logical(LGP) :: found
137  integer(I4B) :: n
138  !
139  if (this%iper /= kper) return
140  !
141  ! -- copy nbound from input context
142  call mem_set_value(this%nbound, 'NBOUND', this%input_mempath, &
143  found)
144  !
145  ! -- convert cellids to node numbers
146  call this%nodelist_update()
147  !
148  ! -- update boundname string list
149  if (this%inamedbound /= 0) then
150  do n = 1, size(this%boundname_cst)
151  this%boundname(n) = this%boundname_cst(n)
152  end do
153  end if
154  end subroutine bndext_rp
155 
156  !> @ brief Deallocate package memory
157  !<
158  subroutine bndext_da(this)
159  ! -- modules
161  ! -- dummy variables
162  class(bndexttype) :: this !< BndExtType object
163  !
164  ! -- deallocate checkin paths
165  call mem_deallocate(this%cellid, 'CELLID', this%memoryPath)
166  call mem_deallocate(this%boundname_cst, 'BOUNDNAME_IDM', this%memoryPath)
167  call mem_deallocate(this%auxvar, 'AUXVAR_IDM', this%memoryPath)
168  !
169  ! -- reassign pointers for base class _da
170  call mem_setptr(this%boundname_cst, 'BOUNDNAME_CST', this%memoryPath)
171  call mem_setptr(this%auxvar, 'AUXVAR', this%memoryPath)
172  !
173  ! -- scalars
174  nullify (this%iper)
175  !
176  ! -- deallocate
177  call this%BndType%bnd_da()
178  end subroutine bndext_da
179 
180  !> @ brief Allocate package scalars
181  !!
182  !! Allocate and initialize base boundary package scalars. This method
183  !! only needs to be overridden if additional scalars are defined
184  !! for a specific package.
185  !!
186  !<
187  subroutine bndext_allocate_scalars(this)
188  ! -- modules
193  ! -- dummy variables
194  class(bndexttype) :: this !< BndExtType object
195  ! -- local variables
196  character(len=LENMEMPATH) :: input_mempath
197  !
198  ! -- set memory path
199  input_mempath = create_mem_path(this%name_model, this%packName, idm_context)
200  !
201  ! -- allocate base BndType scalars
202  call this%BndType%allocate_scalars()
203  !
204  ! -- set pointers to period input data scalars
205  call mem_setptr(this%iper, 'IPER', input_mempath)
206  end subroutine bndext_allocate_scalars
207 
208  !> @ brief Allocate package arrays
209  !!
210  !! Allocate and initialize base boundary package arrays. This method
211  !! only needs to be overridden if additional arrays are defined
212  !! for a specific package.
213  !!
214  !<
215  subroutine bndext_allocate_arrays(this, nodelist, auxvar)
216  ! -- modules
218  ! -- dummy variables
219  class(bndexttype) :: this !< BndExtType object
220  ! -- local variables
221  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist !< package nodelist
222  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar !< package aux variable array
223  !
224  ! -- allocate base BndType arrays
225  call this%BndType%allocate_arrays(nodelist, auxvar)
226  !
227  ! -- set input context pointers
228  call mem_setptr(this%cellid, 'CELLID', this%input_mempath)
229  call mem_setptr(this%boundname_cst, 'BOUNDNAME', this%input_mempath)
230  !
231  ! -- checkin input context pointers
232  call mem_checkin(this%cellid, 'CELLID', this%memoryPath, &
233  'CELLID', this%input_mempath)
234  call mem_checkin(this%boundname_cst, lenboundname, 'BOUNDNAME_IDM', &
235  this%memoryPath, 'BOUNDNAME', this%input_mempath)
236  !
237  if (present(auxvar)) then
238  ! no-op
239  else
240  ! -- set auxvar input context pointer
241  call mem_setptr(this%auxvar, 'AUXVAR', this%input_mempath)
242  !
243  ! -- checkin auxvar input context pointer
244  call mem_checkin(this%auxvar, 'AUXVAR_IDM', this%memoryPath, &
245  'AUXVAR', this%input_mempath)
246  end if
247  end subroutine bndext_allocate_arrays
248 
249  !> @ brief Source package options from input context
250  !<
251  subroutine source_options(this)
252  ! -- modules
253  use memorymanagermodule, only: mem_reallocate, mem_setptr !, get_isize
258  ! -- dummy variables
259  class(bndexttype), intent(inout) :: this !< BndExtType object
260  ! -- local variables
261  type(bndextfoundtype) :: found
262  character(len=LENAUXNAME) :: sfacauxname
263  integer(I4B) :: n
264  !
265  ! -- update defaults with idm sourced values
266  call mem_set_value(this%naux, 'NAUX', this%input_mempath, found%naux)
267  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
268  call mem_set_value(this%iprpak, 'IPRPAK', this%input_mempath, found%iprpak)
269  call mem_set_value(this%iprflow, 'IPRFLOW', this%input_mempath, found%iprflow)
270  call mem_set_value(this%inamedbound, 'BOUNDNAMES', this%input_mempath, &
271  found%boundnames)
272  call mem_set_value(sfacauxname, 'AUXMULTNAME', this%input_mempath, &
273  found%auxmultname)
274  call mem_set_value(this%inewton, 'INEWTON', this%input_mempath, found%inewton)
275  !
276  ! -- log found options
277  call this%log_options(found, sfacauxname)
278  !
279  ! -- reallocate aux arrays if aux variables provided
280  if (found%naux .and. this%naux > 0) then
281  call mem_reallocate(this%auxname, lenauxname, this%naux, &
282  'AUXNAME', this%memoryPath)
283  call mem_reallocate(this%auxname_cst, lenauxname, this%naux, &
284  'AUXNAME_CST', this%memoryPath)
285  call mem_set_value(this%auxname_cst, 'AUXILIARY', this%input_mempath, &
286  found%auxiliary)
287  !
288  do n = 1, this%naux
289  this%auxname(n) = this%auxname_cst(n)
290  end do
291  end if
292  !
293  ! -- save flows option active
294  if (found%ipakcb) this%ipakcb = -1
295  !
296  ! -- auxmultname provided
297  if (found%auxmultname) this%iauxmultcol = -1
298  !
299  !
300  ! -- enforce 0 or 1 OBS6_FILENAME entries in option block
301  if (filein_fname(this%obs%inputFilename, 'OBS6_FILENAME', &
302  this%input_mempath, this%input_fname)) then
303  this%obs%active = .true.
304  this%obs%inUnitObs = getunit()
305  call openfile(this%obs%inUnitObs, this%iout, this%obs%inputFilename, 'OBS')
306  end if
307  !
308  ! -- no newton specified
309  if (found%inewton) this%inewton = 0
310  !
311  ! -- AUXMULTNAME was specified, so find column of auxvar that will be multiplier
312  if (this%iauxmultcol < 0) then
313  !
314  ! -- Error if no aux variable specified
315  if (this%naux == 0) then
316  write (errmsg, '(a,2(1x,a))') &
317  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
318  'but no AUX variables specified.'
319  call store_error(errmsg)
320  end if
321  !
322  ! -- Assign mult column
323  this%iauxmultcol = 0
324  do n = 1, this%naux
325  if (sfacauxname == this%auxname(n)) then
326  this%iauxmultcol = n
327  exit
328  end if
329  end do
330  !
331  ! -- Error if aux variable cannot be found
332  if (this%iauxmultcol == 0) then
333  write (errmsg, '(a,2(1x,a))') &
334  'AUXMULTNAME was specified as', trim(adjustl(sfacauxname)), &
335  'but no AUX variable found with this name.'
336  call store_error(errmsg)
337  end if
338  end if
339  !
340  ! -- terminate if errors were detected
341  if (count_errors() > 0) then
342  call store_error_filename(this%input_fname)
343  end if
344  end subroutine source_options
345 
346  !> @ brief Log package options
347  !<
348  subroutine log_options(this, found, sfacauxname)
349  ! -- modules
350  ! -- dummy variables
351  class(bndexttype), intent(inout) :: this !< BndExtType object
352  type(bndextfoundtype), intent(in) :: found
353  character(len=*), intent(in) :: sfacauxname
354  ! -- local variables
355  ! -- format
356  character(len=*), parameter :: fmtflow = &
357  &"(4x, 'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL')"
358  character(len=*), parameter :: fmttas = &
359  &"(4x, 'TIME-ARRAY SERIES DATA WILL BE READ FROM FILE: ', a)"
360  character(len=*), parameter :: fmtts = &
361  &"(4x, 'TIME-SERIES DATA WILL BE READ FROM FILE: ', a)"
362  character(len=*), parameter :: fmtnme = &
363  &"(a, i0, a)"
364  !
365  ! -- log found options
366  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
367  //' BASE OPTIONS'
368  !
369  if (found%ipakcb) then
370  write (this%iout, fmtflow)
371  end if
372  !
373  if (found%iprpak) then
374  write (this%iout, '(4x,a)') &
375  'LISTS OF '//trim(adjustl(this%text))//' CELLS WILL BE PRINTED.'
376  end if
377  !
378  if (found%iprflow) then
379  write (this%iout, '(4x,a)') trim(adjustl(this%text))// &
380  ' FLOWS WILL BE PRINTED TO LISTING FILE.'
381  end if
382  !
383  if (found%boundnames) then
384  write (this%iout, '(4x,a)') trim(adjustl(this%text))// &
385  ' BOUNDARIES HAVE NAMES IN LAST COLUMN.'
386  end if
387  !
388  if (found%auxmultname) then
389  write (this%iout, '(4x,a,a)') &
390  'AUXILIARY MULTIPLIER NAME: ', sfacauxname
391  end if
392  !
393  if (found%inewton) then
394  write (this%iout, '(4x,a)') &
395  'NEWTON-RAPHSON method disabled for unconfined cells'
396  end if
397  !
398  ! -- close logging block
399  write (this%iout, '(1x,a)') &
400  'END OF '//trim(adjustl(this%text))//' BASE OPTIONS'
401  end subroutine log_options
402 
403  !> @ brief Source package dimensions from input context
404  !<
405  subroutine source_dimensions(this)
407  ! -- dummy variables
408  class(bndexttype), intent(inout) :: this !< BndExtType object
409  ! -- local variables
410  type(bndextfoundtype) :: found
411  !
412  ! -- open dimensions logging block
413  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text))// &
414  ' BASE DIMENSIONS'
415  !
416  ! -- update defaults with idm sourced values
417  call mem_set_value(this%maxbound, 'MAXBOUND', this%input_mempath, &
418  found%maxbound)
419  !
420  write (this%iout, '(4x,a,i7)') 'MAXBOUND = ', this%maxbound
421  !
422  ! -- close logging block
423  write (this%iout, '(1x,a)') &
424  'END OF '//trim(adjustl(this%text))//' BASE DIMENSIONS'
425  !
426  ! -- verify dimensions were set
427  if (this%maxbound <= 0) then
428  write (errmsg, '(a)') 'MAXBOUND must be an integer greater than zero.'
429  call store_error(errmsg)
430  call store_error_filename(this%input_fname)
431  end if
432  !
433  ! -- Call define_listlabel to construct the list label that is written
434  ! when PRINT_INPUT option is used.
435  call this%define_listlabel()
436  end subroutine source_dimensions
437 
438  !> @ brief Update package nodelist
439  !!
440  !! Convert period updated cellids to node numbers.
441  !!
442  !<
443  subroutine nodelist_update(this)
444  ! -- modules
445  use simvariablesmodule, only: errmsg
446  ! -- dummy
447  class(bndexttype) :: this !< BndExtType object
448  ! -- local
449  integer(I4B), dimension(:), pointer :: cellid
450  integer(I4B) :: n, nodeu, noder
451  character(len=LINELENGTH) :: nodestr
452  !
453  ! -- update nodelist
454  do n = 1, this%nbound
455  !
456  ! -- set cellid
457  cellid => this%cellid(:, n)
458  !
459  ! -- ensure cellid is valid, store an error otherwise
460  call this%check_cellid(n, cellid, this%dis%mshape, this%dis%ndim)
461  !
462  ! -- Determine user node number
463  if (this%dis%ndim == 1) then
464  nodeu = cellid(1)
465  elseif (this%dis%ndim == 2) then
466  nodeu = get_node(cellid(1), 1, cellid(2), &
467  this%dis%mshape(1), 1, &
468  this%dis%mshape(2))
469  else
470  nodeu = get_node(cellid(1), cellid(2), cellid(3), &
471  this%dis%mshape(1), &
472  this%dis%mshape(2), &
473  this%dis%mshape(3))
474  end if
475  !
476  ! -- update the nodelist
477  if (this%dis%nodes < this%dis%nodesuser) then
478  ! -- convert user to reduced node numbers
479  noder = this%dis%get_nodenumber(nodeu, 0)
480  if (noder <= 0) then
481  call this%dis%nodeu_to_string(nodeu, nodestr)
482  write (errmsg, *) &
483  ' Cell is outside active grid domain: '// &
484  trim(adjustl(nodestr))
485  call store_error(errmsg)
486  end if
487  this%nodelist(n) = noder
488  else
489  this%nodelist(n) = nodeu
490  end if
491  end do
492  !
493  ! -- exit if errors were found
494  if (count_errors() > 0) then
495  write (errmsg, *) count_errors(), ' errors encountered.'
496  call store_error(errmsg)
497  call store_error_filename(this%input_fname)
498  end if
499  end subroutine nodelist_update
500 
501  !> @ brief Check for valid cellid
502  !<
503  subroutine check_cellid(this, ii, cellid, mshape, ndim)
504  ! -- modules
505  use simvariablesmodule, only: errmsg
506  ! -- dummy
507  class(bndexttype) :: this !< BndExtType object
508  ! -- local
509  integer(I4B), intent(in) :: ii
510  integer(I4B), dimension(:), intent(in) :: cellid !< cellid
511  integer(I4B), dimension(:), intent(in) :: mshape !< model shape
512  integer(I4B), intent(in) :: ndim !< size of mshape
513  character(len=20) :: cellstr, mshstr
514  character(len=*), parameter :: fmterr = &
515  "('List entry ',i0,' contains cellid ',a,' but this cellid is invalid &
516  &for model with shape ', a)"
517  character(len=*), parameter :: fmtndim1 = &
518  "('(',i0,')')"
519  character(len=*), parameter :: fmtndim2 = &
520  "('(',i0,',',i0,')')"
521  character(len=*), parameter :: fmtndim3 = &
522  "('(',i0,',',i0,',',i0,')')"
523  select case (ndim)
524  case (1)
525  !
526  if (cellid(1) < 1 .or. cellid(1) > mshape(1)) then
527  write (cellstr, fmtndim1) cellid(1)
528  write (mshstr, fmtndim1) mshape(1)
529  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
530  call store_error(errmsg)
531  end if
532  !
533  case (2)
534  !
535  if (cellid(1) < 1 .or. cellid(1) > mshape(1) .or. &
536  cellid(2) < 1 .or. cellid(2) > mshape(2)) then
537  write (cellstr, fmtndim2) cellid(1), cellid(2)
538  write (mshstr, fmtndim2) mshape(1), mshape(2)
539  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
540  call store_error(errmsg)
541  end if
542  !
543  case (3)
544  !
545  if (cellid(1) < 1 .or. cellid(1) > mshape(1) .or. &
546  cellid(2) < 1 .or. cellid(2) > mshape(2) .or. &
547  cellid(3) < 1 .or. cellid(3) > mshape(3)) then
548  write (cellstr, fmtndim3) cellid(1), cellid(2), cellid(3)
549  write (mshstr, fmtndim3) mshape(1), mshape(2), mshape(3)
550  write (errmsg, fmterr) ii, trim(adjustl(cellstr)), trim(adjustl(mshstr))
551  call store_error(errmsg)
552  end if
553  !
554  case default
555  end select
556  end subroutine check_cellid
557 
558  !> @ brief Log package list input
559  !!
560  !! Log period list based input. This routine requires a package specific
561  !! bound_value() routine to report accurate bound values.
562  !!
563  !<
564  subroutine write_list(this)
565  ! -- modules
568  use inputoutputmodule, only: ulstlb
569  use tablemodule, only: tabletype, table_cr
570  ! -- dummy
571  class(bndexttype) :: this !< BndExtType object
572  ! -- local
573  character(len=10) :: cpos
574  character(len=LINELENGTH) :: tag
575  character(len=LINELENGTH), allocatable, dimension(:) :: words
576  integer(I4B) :: ntabrows
577  integer(I4B) :: ntabcols
578  integer(I4B) :: ipos
579  integer(I4B) :: ii, jj, i, j, k, nod
580  integer(I4B) :: ldim
581  integer(I4B) :: naux
582  type(tabletype), pointer :: inputtab => null()
583  ! -- formats
584  character(len=LINELENGTH) :: fmtlstbn
585  !
586  ! -- Determine sizes
587  ldim = this%ncolbnd
588  naux = size(this%auxvar, 1)
589  !
590  ! -- dimension table
591  ntabrows = this%nbound
592  !
593  ! -- start building format statement to parse this%label, which
594  ! contains the column headers (except for boundname and auxnames)
595  ipos = index(this%listlabel, 'NO.')
596  if (ipos /= 0) then
597  write (cpos, '(i10)') ipos + 3
598  fmtlstbn = '(a'//trim(adjustl(cpos))
599  else
600  fmtlstbn = '(a7'
601  end if
602  ! -- sequence number, layer, row, and column.
603  if (size(this%dis%mshape) == 3) then
604  ntabcols = 4
605  fmtlstbn = trim(fmtlstbn)//',a7,a7,a7'
606  !
607  ! -- sequence number, layer, and cell2d.
608  else if (size(this%dis%mshape) == 2) then
609  ntabcols = 3
610  fmtlstbn = trim(fmtlstbn)//',a7,a7'
611  !
612  ! -- sequence number and node.
613  else
614  ntabcols = 2
615  fmtlstbn = trim(fmtlstbn)//',a7'
616  end if
617  !
618  ! -- Add fields for non-optional real values
619  ntabcols = ntabcols + ldim
620  do i = 1, ldim
621  fmtlstbn = trim(fmtlstbn)//',a16'
622  end do
623  !
624  ! -- Add field for boundary name
625  if (this%inamedbound == 1) then
626  ntabcols = ntabcols + 1
627  fmtlstbn = trim(fmtlstbn)//',a16'
628  end if
629  !
630  ! -- Add fields for auxiliary variables
631  ntabcols = ntabcols + naux
632  do i = 1, naux
633  fmtlstbn = trim(fmtlstbn)//',a16'
634  end do
635  fmtlstbn = trim(fmtlstbn)//')'
636  !
637  ! -- allocate words
638  allocate (words(ntabcols))
639  !
640  ! -- parse this%listlabel into words
641  read (this%listlabel, fmtlstbn) (words(i), i=1, ntabcols)
642  !
643  ! -- initialize the input table object
644  call table_cr(inputtab, ' ', ' ')
645  call inputtab%table_df(ntabrows, ntabcols, this%iout)
646  !
647  ! -- add the columns
648  ipos = 1
649  call inputtab%initialize_column(words(ipos), 10, alignment=tabcenter)
650  !
651  ! -- discretization
652  do i = 1, size(this%dis%mshape)
653  ipos = ipos + 1
654  call inputtab%initialize_column(words(ipos), 7, alignment=tabcenter)
655  end do
656  !
657  ! -- non-optional variables
658  do i = 1, ldim
659  ipos = ipos + 1
660  call inputtab%initialize_column(words(ipos), 16, alignment=tabcenter)
661  end do
662  !
663  ! -- boundname
664  if (this%inamedbound == 1) then
665  ipos = ipos + 1
666  tag = 'BOUNDNAME'
667  call inputtab%initialize_column(tag, lenboundname, alignment=tableft)
668  end if
669  !
670  ! -- aux variables
671  do i = 1, naux
672  call inputtab%initialize_column(this%auxname(i), 16, alignment=tabcenter)
673  end do
674  !
675  ! -- Write the table
676  do ii = 1, this%nbound
677  call inputtab%add_term(ii)
678  !
679  ! -- discretization
680  if (size(this%dis%mshape) == 3) then
681  nod = this%nodelist(ii)
682  call get_ijk(nod, this%dis%mshape(2), this%dis%mshape(3), &
683  this%dis%mshape(1), i, j, k)
684  call inputtab%add_term(k)
685  call inputtab%add_term(i)
686  call inputtab%add_term(j)
687  else if (size(this%dis%mshape) == 2) then
688  nod = this%nodelist(ii)
689  call get_ijk(nod, 1, this%dis%mshape(2), this%dis%mshape(1), i, j, k)
690  call inputtab%add_term(k)
691  call inputtab%add_term(j)
692  else
693  nod = this%nodelist(ii)
694  call inputtab%add_term(nod)
695  end if
696  !
697  ! -- non-optional variables
698  do jj = 1, ldim
699  call inputtab%add_term(this%bound_value(jj, ii))
700  end do
701  !
702  ! -- boundname
703  if (this%inamedbound == 1) then
704  call inputtab%add_term(this%boundname(ii))
705  end if
706  !
707  ! -- aux variables
708  do jj = 1, naux
709  call inputtab%add_term(this%auxvar(jj, ii))
710  end do
711  end do
712  !
713  ! -- deallocate the local variables
714  call inputtab%table_da()
715  deallocate (inputtab)
716  nullify (inputtab)
717  deallocate (words)
718  end subroutine write_list
719 
720  !> @ brief Return a bound value
721  !!
722  !! Return a bound value associated with an ncolbnd index
723  !! and row. This function should be overridden in the
724  !! derived package class.
725  !!
726  !<
727  function bound_value(this, col, row) result(bndval)
728  ! -- modules
729  use constantsmodule, only: dnodata
730  ! -- dummy variables
731  class(bndexttype), intent(inout) :: this !< BndExtType object
732  integer(I4B), intent(in) :: col
733  integer(I4B), intent(in) :: row
734  ! -- result
735  real(dp) :: bndval
736  !
737  ! -- override this return value by redefining this
738  ! routine in the derived package.
739  bndval = dnodata
740  end function bound_value
741 
742 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 nodelist_update(this)
@ brief Update package nodelist
subroutine bndext_df(this, neq, dis)
@ brief Define boundary package options and dimensions
subroutine bndext_da(this)
@ brief Deallocate package memory
subroutine log_options(this, found, sfacauxname)
@ brief Log package options
subroutine source_dimensions(this)
@ brief Source package dimensions from input context
subroutine check_cellid(this, ii, cellid, mshape, ndim)
@ brief Check for valid cellid
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 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
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
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
character(len=linelength) idm_context
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