MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
DiscretizationBase.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b, lgp
6  dis, disv, disu, &
7  dis2d, disv2d, disu2d, &
9 
13  use simvariablesmodule, only: errmsg
14  use simmodule, only: count_errors, store_error, &
20  use tdismodule, only: kstp, kper, pertim, totim, delt
23 
24  implicit none
25 
26  private
27  public :: disbasetype
28  public :: dis_transform_xy
29 
30  type :: disbasetype
31  character(len=LENMEMPATH) :: memorypath !< path for memory allocation
32  character(len=LENMEMPATH) :: input_mempath = '' !< input context mempath
33  character(len=LENMODELNAME), pointer :: name_model => null() !< name of the model
34  character(len=LINELENGTH), pointer :: input_fname => null() !< input file name
35  character(len=LINELENGTH), pointer :: output_fname => null() !< output file name
36  integer(I4B), pointer :: inunit => null() !< unit number for input file
37  integer(I4B), pointer :: iout => null() !< unit number for output file
38  integer(I4B), pointer :: nodes => null() !< number of nodes in solution
39  integer(I4B), pointer :: nodesuser => null() !< number of user nodes (same as nodes for disu grid)
40  integer(I4B), pointer :: nja => null() !< number of connections plus number of nodes
41  integer(I4B), pointer :: njas => null() !< (nja-nodes)/2
42  integer(I4B), pointer :: lenuni => null() !< length unit
43  integer(I4B), pointer :: ndim => null() !< number of spatial model dimensions (1 for disu grid)
44  integer(I4B), pointer :: icondir => null() !< flag indicating if grid has enough info to calculate connection vectors
45  integer(I4B), pointer :: nogrb => null() !< don't write binary grid file
46  real(dp), dimension(:), pointer, contiguous :: xc => null() !< x-coordinate of the cell center
47  real(dp), dimension(:), pointer, contiguous :: yc => null() !< y-coordinate of the cell center
48  real(dp), pointer :: yorigin => null() !< y-position of the lower-left grid corner (default is 0.)
49  real(dp), pointer :: xorigin => null() !< x-position of the lower-left grid corner (default is 0.)
50  real(dp), pointer :: angrot => null() !< counter-clockwise rotation angle of the lower-left corner (default is 0.0)
51  integer(I4B), dimension(:), pointer, contiguous :: mshape => null() !< shape of the model; (nodes) for DisBaseType
52  real(dp), dimension(:), pointer, contiguous :: top => null() !< (size:nodes) cell top elevation
53  real(dp), dimension(:), pointer, contiguous :: bot => null() !< (size:nodes) cell bottom elevation
54  real(dp), dimension(:), pointer, contiguous :: area => null() !< (size:nodes) cell area, in plan view
55  type(connectionstype), pointer :: con => null() !< connections object
56  type(blockparsertype) :: parser !< object to read blocks
57  real(dp), dimension(:), pointer, contiguous :: dbuff => null() !< helper double array of size nodesuser
58  integer(I4B), dimension(:), pointer, contiguous :: ibuff => null() !< helper int array of size nodesuser
59  integer(I4B), dimension(:), pointer, contiguous :: nodereduced => null() !< (size:nodesuser)contains reduced nodenumber (size 0 if not reduced); -1 means vertical pass through, 0 is idomain = 0
60  integer(I4B), dimension(:), pointer, contiguous :: nodeuser => null() !< (size:nodes) given a reduced nodenumber, provide the user nodenumber (size 0 if not reduced)
61  contains
62  procedure :: dis_df
63  procedure :: dis_ac
64  procedure :: dis_mc
65  procedure :: dis_ar
66  procedure :: dis_da
67  ! -- helper functions
68  !
69  ! -- get_nodenumber is an overloaded integer function that will always
70  ! return the reduced nodenumber. For all grids, get_nodenumber can
71  ! be passed the user nodenumber. For some other grids, it can also
72  ! be passed an index. For dis3d the index is k, i, j, and for
73  ! disv the index is k, n.
74  generic :: get_nodenumber => get_nodenumber_idx1, &
77  procedure :: get_nodenumber_idx1
78  procedure :: get_nodenumber_idx2
79  procedure :: get_nodenumber_idx3
80  procedure :: get_nodeuser
81  procedure :: nodeu_to_string
82  procedure :: nodeu_to_array
83  procedure :: nodeu_from_string
84  procedure :: nodeu_from_cellid
85  procedure :: noder_from_string
86  procedure :: noder_from_cellid
87  procedure :: connection_normal
88  procedure :: connection_vector
89  procedure :: get_dis_type
90  procedure :: get_dis_enum
91  procedure :: supports_layers
92  procedure :: allocate_scalars
93  procedure :: allocate_arrays
94  procedure :: get_ncpl
95  procedure :: get_cell_volume
96  procedure :: get_polyverts
97  procedure :: write_grb
98  !
99  procedure :: read_int_array
100  procedure :: read_dbl_array
101  generic, public :: read_grid_array => read_int_array, read_dbl_array
102  procedure, public :: read_layer_array
103  procedure :: fill_int_array
104  procedure :: fill_dbl_array
105  generic, public :: fill_grid_array => fill_int_array, fill_dbl_array
106  procedure, public :: read_list
107  !
108  procedure, public :: record_array
109  procedure, public :: record_connection_array
110  procedure, public :: noder_to_string
111  procedure, public :: noder_to_array
112  procedure, public :: record_srcdst_list_header
113  procedure, private :: record_srcdst_list_entry
114  generic, public :: record_mf6_list_entry => record_srcdst_list_entry
115  procedure, public :: nlarray_to_nodelist
116  procedure, public :: highest_active
117  procedure, public :: highest_saturated
118  procedure, public :: get_area
119  procedure, public :: get_area_factor
120  procedure, public :: get_flow_width
121  procedure, public :: is_3d
122  procedure, public :: is_2d
123  procedure, public :: is_1d
124  end type disbasetype
125 
126 contains
127 
128  !> @brief Define the discretization
129  subroutine dis_df(this)
130  class(disbasetype) :: this
131  call store_error('Programmer error: dis_df must be overridden', &
132  terminate=.true.)
133  end subroutine dis_df
134 
135  !> @brief Add connections to sparse cell connectivity matrix
136  subroutine dis_ac(this, moffset, sparse)
137  ! -- modules
138  use sparsemodule, only: sparsematrix
139  ! -- dummy
140  class(disbasetype) :: this
141  integer(I4B), intent(in) :: moffset
142  type(sparsematrix), intent(inout) :: sparse
143  ! -- local
144  integer(I4B) :: i, j, ipos, iglo, jglo
145  !
146  do i = 1, this%nodes
147  do ipos = this%con%ia(i), this%con%ia(i + 1) - 1
148  j = this%con%ja(ipos)
149  iglo = i + moffset
150  jglo = j + moffset
151  call sparse%addconnection(iglo, jglo, 1)
152  end do
153  end do
154  end subroutine dis_ac
155 
156  !> @brief Map cell connections in the numerical solution coefficient matrix.
157  subroutine dis_mc(this, moffset, idxglo, matrix_sln)
158  ! -- dummy
159  class(disbasetype) :: this
160  integer(I4B), intent(in) :: moffset
161  integer(I4B), dimension(:), intent(inout) :: idxglo
162  class(matrixbasetype), pointer :: matrix_sln
163  ! -- local
164  integer(I4B) :: i, j, ipos, iglo, jglo
165  !
166  do i = 1, this%nodes
167  iglo = i + moffset
168  do ipos = this%con%ia(i), this%con%ia(i + 1) - 1
169  j = this%con%ja(ipos)
170  jglo = j + moffset
171  idxglo(ipos) = matrix_sln%get_position(iglo, jglo)
172  end do
173  end do
174  end subroutine dis_mc
175 
176  !> @brief Allocate and setup variables, and write binary grid file.
177  subroutine dis_ar(this, icelltype)
178  ! -- dummy
179  class(disbasetype) :: this
180  integer(I4B), dimension(:), intent(in) :: icelltype
181  ! -- local
182  integer(I4B), dimension(:), allocatable :: ict
183  integer(I4B) :: nu, nr
184  !
185  ! -- Expand icelltype to full grid; fill with 0 if cell is excluded
186  allocate (ict(this%nodesuser))
187  do nu = 1, this%nodesuser
188  nr = this%get_nodenumber(nu, 0)
189  if (nr > 0) then
190  ict(nu) = icelltype(nr)
191  else
192  ict(nu) = 0
193  end if
194  end do
195  !
196  if (this%nogrb == 0) call this%write_grb(ict)
197  end subroutine dis_ar
198 
199  !> @brief Write a binary grid file
200  subroutine write_grb(this, icelltype)
201  class(disbasetype) :: this
202  integer(I4B), dimension(:), intent(in) :: icelltype
203  call store_error('Programmer error: write_grb must be overridden', &
204  terminate=.true.)
205  end subroutine write_grb
206 
207  !> @brier Deallocate variables
208  subroutine dis_da(this)
209  ! -- modules
211  ! -- dummy
212  class(disbasetype) :: this
213  !
214  ! -- Strings
215  deallocate (this%name_model)
216  deallocate (this%input_fname)
217  deallocate (this%output_fname)
218  !
219  ! -- Scalars
220  call mem_deallocate(this%inunit)
221  call mem_deallocate(this%iout)
222  call mem_deallocate(this%nodes)
223  call mem_deallocate(this%nodesuser)
224  call mem_deallocate(this%ndim)
225  call mem_deallocate(this%icondir)
226  call mem_deallocate(this%nogrb)
227  call mem_deallocate(this%xorigin)
228  call mem_deallocate(this%yorigin)
229  call mem_deallocate(this%angrot)
230  call mem_deallocate(this%nja)
231  call mem_deallocate(this%njas)
232  call mem_deallocate(this%lenuni)
233  !
234  ! -- Arrays
235  call mem_deallocate(this%mshape)
236  call mem_deallocate(this%xc)
237  call mem_deallocate(this%yc)
238  call mem_deallocate(this%top)
239  call mem_deallocate(this%bot)
240  call mem_deallocate(this%area)
241  call mem_deallocate(this%dbuff)
242  call mem_deallocate(this%ibuff)
243  !
244  ! -- Connections
245  call this%con%con_da()
246  deallocate (this%con)
247  end subroutine dis_da
248 
249  !> @brief Convert a user nodenumber to a string (nodenumber), (k,j), or (k,i,j)
250  subroutine nodeu_to_string(this, nodeu, str)
251  class(disbasetype) :: this
252  integer(I4B), intent(in) :: nodeu
253  character(len=*), intent(inout) :: str
254 
255  call store_error('Programmer error: nodeu_to_string must be overridden', &
256  terminate=.true.)
257  end subroutine nodeu_to_string
258 
259  !> @brief Convert a user nodenumber to an array (nodenumber), (k,j), or (k,i,j)
260  subroutine nodeu_to_array(this, nodeu, arr)
261  class(disbasetype) :: this
262  integer(I4B), intent(in) :: nodeu
263  integer(I4B), dimension(:), intent(inout) :: arr
264 
265  call store_error('Programmer error: nodeu_to_array must be overridden', &
266  terminate=.true.)
267  end subroutine nodeu_to_array
268 
269  !> @brief Convert a reduced nodenumber to a user node number
270  function get_nodeuser(this, noder) result(nodenumber)
271  class(disbasetype) :: this
272  integer(I4B), intent(in) :: noder
273  integer(I4B) :: nodenumber
274 
275  if (this%nodes < this%nodesuser) then
276  nodenumber = this%nodeuser(noder)
277  else
278  nodenumber = noder
279  end if
280  end function get_nodeuser
281 
282  function get_nodenumber_idx1(this, nodeu, icheck) result(nodenumber)
283  class(disbasetype), intent(in) :: this
284  integer(I4B), intent(in) :: nodeu
285  integer(I4B), intent(in) :: icheck
286  integer(I4B) :: nodenumber
287 
288  nodenumber = 0
289  call store_error('Programmer error: get_nodenumber_idx1 must be overridden', &
290  terminate=.true.)
291  end function get_nodenumber_idx1
292 
293  function get_nodenumber_idx2(this, k, j, icheck) result(nodenumber)
294  class(disbasetype), intent(in) :: this
295  integer(I4B), intent(in) :: k, j
296  integer(I4B), intent(in) :: icheck
297  integer(I4B) :: nodenumber
298 
299  nodenumber = 0
300  call store_error('Programmer error: get_nodenumber_idx2 must be overridden', &
301  terminate=.true.)
302  end function get_nodenumber_idx2
303 
304  function get_nodenumber_idx3(this, k, i, j, icheck) result(nodenumber)
305  class(disbasetype), intent(in) :: this
306  integer(I4B), intent(in) :: k, i, j
307  integer(I4B), intent(in) :: icheck
308  integer(I4B) :: nodenumber
309 
310  nodenumber = 0
311  call store_error('Programmer error: get_nodenumber_idx3 must be overridden', &
312  terminate=.true.)
313  end function get_nodenumber_idx3
314 
315  !> @brief Get normal vector components between the cell and a given neighbor.
316  !! The normal points outward from the shared face between noden and nodem.
317  subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, &
318  ipos)
319  class(disbasetype) :: this
320  integer(I4B), intent(in) :: noden !< cell (reduced nn)
321  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
322  integer(I4B), intent(in) :: ihc !< horizontal connection flag
323  real(DP), intent(inout) :: xcomp
324  real(DP), intent(inout) :: ycomp
325  real(DP), intent(inout) :: zcomp
326  integer(I4B), intent(in) :: ipos
327 
328  call store_error('Programmer error: connection_normal must be overridden', &
329  terminate=.true.)
330  end subroutine connection_normal
331 
332  !> @brief Get unit vector components between the cell and a given neighbor.
333  !! Saturation must be provided to compute cell center vertical coordinates.
334  !! Also return the straight-line connection length.
335  subroutine connection_vector(this, noden, nodem, nozee, satn, satm, ihc, &
336  xcomp, ycomp, zcomp, conlen)
337  class(disbasetype) :: this
338  integer(I4B), intent(in) :: noden !< cell (reduced nn)
339  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
340  logical, intent(in) :: nozee
341  real(DP), intent(in) :: satn
342  real(DP), intent(in) :: satm
343  integer(I4B), intent(in) :: ihc !< horizontal connection flag
344  real(DP), intent(inout) :: xcomp
345  real(DP), intent(inout) :: ycomp
346  real(DP), intent(inout) :: zcomp
347  real(DP), intent(inout) :: conlen
348 
349  call store_error('Programmer error: connection_vector must be overridden', &
350  terminate=.true.)
351  end subroutine connection_vector
352 
353  !> @brief Get global (x, y) coordinates from cell-local coordinates.
354  subroutine dis_transform_xy(x, y, xorigin, yorigin, angrot, xglo, yglo)
355  real(dp), intent(in) :: x !< the cell-x coordinate to transform
356  real(dp), intent(in) :: y !< the cell-y coordinate to transform
357  real(dp), intent(in) :: xorigin !< the cell-y coordinate to transform
358  real(dp), intent(in) :: yorigin !< the cell-y coordinate to transform
359  real(dp), intent(in) :: angrot !< the cell-y coordinate to transform
360  real(dp), intent(out) :: xglo !< the global cell-x coordinate
361  real(dp), intent(out) :: yglo !< the global cell-y coordinate
362  ! local
363  real(dp) :: ang
364 
365  xglo = x
366  yglo = y
367 
368  ! first _rotate_ to 'real world'
369  ang = angrot * dpio180
370  if (ang /= dzero) then
371  xglo = x * cos(ang) - y * sin(ang)
372  yglo = x * sin(ang) + y * cos(ang)
373  end if
374 
375  ! then _translate_
376  xglo = xglo + xorigin
377  yglo = yglo + yorigin
378  end subroutine dis_transform_xy
379 
380  !> @brief Get the discretization type (DIS, DISV, or DISU)
381  subroutine get_dis_type(this, dis_type)
382  class(disbasetype), intent(in) :: this
383  character(len=*), intent(out) :: dis_type
384 
385  dis_type = "Not implemented"
386  call store_error('Programmer error: get_dis_type must be overridden', &
387  terminate=.true.)
388  end subroutine get_dis_type
389 
390  !> @brief Get the discretization type enumeration
391  function get_dis_enum(this) result(dis_enum)
392  use constantsmodule, only: disundef
393  class(disbasetype), intent(in) :: this
394  integer(I4B) :: dis_enum
395 
396  dis_enum = disundef
397  call store_error('Programmer error: get_dis_enum must be overridden', &
398  terminate=.true.)
399  end function get_dis_enum
400 
401  !> @brief Allocate and initialize scalar variables
402  subroutine allocate_scalars(this, name_model, input_mempath)
403  ! -- dummy
404  class(disbasetype) :: this
405  character(len=*), intent(in) :: name_model
406  character(len=*), intent(in) :: input_mempath
407  logical(LGP) :: found
408  !
409  ! -- Create memory path
410  this%memoryPath = create_mem_path(name_model, 'DIS')
411  !
412  ! -- Allocate
413  allocate (this%name_model)
414  allocate (this%input_fname)
415  allocate (this%output_fname)
416  !
417  call mem_allocate(this%inunit, 'INUNIT', this%memoryPath)
418  call mem_allocate(this%iout, 'IOUT', this%memoryPath)
419  call mem_allocate(this%nodes, 'NODES', this%memoryPath)
420  call mem_allocate(this%nodesuser, 'NODESUSER', this%memoryPath)
421  call mem_allocate(this%ndim, 'NDIM', this%memoryPath)
422  call mem_allocate(this%icondir, 'ICONDIR', this%memoryPath)
423  call mem_allocate(this%nogrb, 'NOGRB', this%memoryPath)
424  call mem_allocate(this%xorigin, 'XORIGIN', this%memoryPath)
425  call mem_allocate(this%yorigin, 'YORIGIN', this%memoryPath)
426  call mem_allocate(this%angrot, 'ANGROT', this%memoryPath)
427  call mem_allocate(this%nja, 'NJA', this%memoryPath)
428  call mem_allocate(this%njas, 'NJAS', this%memoryPath)
429  call mem_allocate(this%lenuni, 'LENUNI', this%memoryPath)
430  !
431  ! -- Initialize
432  this%name_model = name_model
433  this%input_mempath = input_mempath
434  this%input_fname = ''
435  this%output_fname = ''
436  this%inunit = 0
437  this%iout = 0
438  this%nodes = 0
439  this%nodesuser = 0
440  this%ndim = 1
441  this%icondir = 1
442  this%nogrb = 0
443  this%xorigin = dzero
444  this%yorigin = dzero
445  this%angrot = dzero
446  this%nja = 0
447  this%njas = 0
448  this%lenuni = 0
449  !
450  ! -- update input and output filenames
451  call mem_set_value(this%input_fname, 'INPUT_FNAME', &
452  this%input_mempath, found)
453  call mem_set_value(this%output_fname, 'GRB6_FILENAME', &
454  this%input_mempath, found)
455  if (.not. found) then
456  this%output_fname = trim(this%input_fname)//'.grb'
457  end if
458  end subroutine allocate_scalars
459 
460  !> @brief Allocate and initialize arrays
461  subroutine allocate_arrays(this)
462  class(disbasetype) :: this
463  integer :: isize
464  !
465  ! -- Allocate
466  call mem_allocate(this%mshape, this%ndim, 'MSHAPE', this%memoryPath)
467  call mem_allocate(this%xc, this%nodes, 'XC', this%memoryPath)
468  call mem_allocate(this%yc, this%nodes, 'YC', this%memoryPath)
469  call mem_allocate(this%top, this%nodes, 'TOP', this%memoryPath)
470  call mem_allocate(this%bot, this%nodes, 'BOT', this%memoryPath)
471  call mem_allocate(this%area, this%nodes, 'AREA', this%memoryPath)
472  !
473  ! -- Initialize
474  this%mshape(1) = this%nodes
475  !
476  ! -- Determine size of buff memory
477  if (this%nodes < this%nodesuser) then
478  isize = this%nodesuser
479  else
480  isize = this%nodes
481  end if
482  !
483  ! -- Allocate the arrays
484  call mem_allocate(this%dbuff, isize, 'DBUFF', this%name_model)
485  call mem_allocate(this%ibuff, isize, 'IBUFF', this%name_model)
486  end subroutine allocate_arrays
487 
488  !> @brief Convert a string to a user nodenumber.
489  !!
490  !! If DIS or DISV, read indices. If DISU, read user node number directly.
491  !! If flag_string is present and true, the first token may be
492  !! non-numeric (e.g. boundary name). In this case, return -2.
493  !<
494  function nodeu_from_string(this, lloc, istart, istop, in, iout, line, &
495  flag_string, allow_zero) result(nodeu)
496  ! -- dummy
497  class(disbasetype) :: this
498  integer(I4B), intent(inout) :: lloc
499  integer(I4B), intent(inout) :: istart
500  integer(I4B), intent(inout) :: istop
501  integer(I4B), intent(in) :: in
502  integer(I4B), intent(in) :: iout
503  character(len=*), intent(inout) :: line
504  logical, optional, intent(in) :: flag_string
505  logical, optional, intent(in) :: allow_zero
506  integer(I4B) :: nodeu
507 
508  nodeu = 0
509  call store_error('Programmer error: nodeu_from_string must be overridden', &
510  terminate=.true.)
511  end function nodeu_from_string
512 
513  !> @brief Convert a cellid string to a user nodenumber.
514  !!
515  !! If flag_string is present and true, the first token may be
516  !! non-numeric (e.g. boundary name). In this case, return -2.
517  !!
518  !! If allow_zero is present and true, and all indices are zero, the
519  !! result can be zero. If allow_zero is false, a zero in any index is an error.
520  !<
521  function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, &
522  allow_zero) result(nodeu)
523  ! -- dummy
524  class(disbasetype) :: this
525  character(len=*), intent(inout) :: cellid
526  integer(I4B), intent(in) :: inunit
527  integer(I4B), intent(in) :: iout
528  logical, optional, intent(in) :: flag_string
529  logical, optional, intent(in) :: allow_zero
530  integer(I4B) :: nodeu
531 
532  nodeu = 0
533  call store_error('Programmer error: nodeu_from_cellid must be overridden', &
534  terminate=.true.)
535  end function nodeu_from_cellid
536 
537  !> @brief Convert a string to a reduced nodenumber.
538  !!
539  !! If the model is unstructured; just read user nodenumber.
540  !! If flag_string argument is present and true, the first token in string
541  !! is allowed to be a string (e.g. boundary name). In this case, if a string
542  !! is encountered, return value as -2.
543  !<
544  function noder_from_string(this, lloc, istart, istop, in, iout, line, &
545  flag_string) result(noder)
546  ! -- dummy
547  class(disbasetype) :: this
548  integer(I4B), intent(inout) :: lloc
549  integer(I4B), intent(inout) :: istart
550  integer(I4B), intent(inout) :: istop
551  integer(I4B), intent(in) :: in
552  integer(I4B), intent(in) :: iout
553  character(len=*), intent(inout) :: line
554  logical, optional, intent(in) :: flag_string
555  integer(I4B) :: noder
556  ! -- local
557  integer(I4B) :: nodeu
558  character(len=LINELENGTH) :: nodestr
559  logical :: flag_string_local
560  !
561  if (present(flag_string)) then
562  flag_string_local = flag_string
563  else
564  flag_string_local = .false.
565  end if
566  nodeu = this%nodeu_from_string(lloc, istart, istop, in, iout, line, &
567  flag_string_local)
568  !
569  ! -- Convert user-based nodenumber to reduced node number
570  if (nodeu > 0) then
571  noder = this%get_nodenumber(nodeu, 0)
572  else
573  noder = nodeu
574  end if
575  if (noder <= 0 .and. .not. flag_string_local) then
576  call this%nodeu_to_string(nodeu, nodestr)
577  write (errmsg, *) &
578  ' Cell is outside active grid domain: '// &
579  trim(adjustl(nodestr))
580  call store_error(errmsg)
581  end if
582  end function noder_from_string
583 
584  !> @brief Convert cellid string to reduced nodenumber
585  !!
586  !! If flag_string argument is present and true, the first token in string
587  !! is allowed to be a string (e.g. boundary name). In this case, if a string
588  !! is encountered, return value as -2.
589  !! If allow_zero argument is present and true, if all indices equal zero, the
590  !! result can be zero. If allow_zero is false, a zero in any index is an error.
591  !<
592  function noder_from_cellid(this, cellid, inunit, iout, flag_string, &
593  allow_zero) result(noder)
594  ! -- return
595  integer(I4B) :: noder
596  ! -- dummy
597  class(disbasetype) :: this
598  character(len=*), intent(inout) :: cellid
599  integer(I4B), intent(in) :: inunit
600  integer(I4B), intent(in) :: iout
601  logical, optional, intent(in) :: flag_string
602  logical, optional, intent(in) :: allow_zero
603  ! -- local
604  integer(I4B) :: nodeu
605  logical :: allowzerolocal
606  character(len=LINELENGTH) :: nodestr
607  logical :: flag_string_local
608  !
609  if (present(flag_string)) then
610  flag_string_local = flag_string
611  else
612  flag_string_local = .false.
613  end if
614  if (present(allow_zero)) then
615  allowzerolocal = allow_zero
616  else
617  allowzerolocal = .false.
618  end if
619  !
620  nodeu = this%nodeu_from_cellid(cellid, inunit, iout, flag_string_local, &
621  allowzerolocal)
622  !
623  ! -- Convert user-based nodenumber to reduced node number
624  if (nodeu > 0) then
625  noder = this%get_nodenumber(nodeu, 0)
626  else
627  noder = nodeu
628  end if
629  if (noder <= 0 .and. .not. flag_string_local) then
630  call this%nodeu_to_string(nodeu, nodestr)
631  write (errmsg, *) &
632  ' Cell is outside active grid domain: '// &
633  trim(adjustl(nodestr))
634  call store_error(errmsg)
635  end if
636  end function noder_from_cellid
637 
638  !> @brief Indicates whether the grid discretization supports layers.
639  logical function supports_layers(this)
640  class(disbasetype) :: this
641  supports_layers = .false.
642  call store_error('Programmer error: supports_layers must be overridden', &
643  terminate=.true.)
644  end function supports_layers
645 
646  !> @brief Return number of cells per layer.
647  !! This is nodes for a DISU grid, as there are no layers.
648  function get_ncpl(this)
649  integer(I4B) :: get_ncpl
650  class(disbasetype) :: this
651  get_ncpl = 0
652  call store_error('Programmer error: get_ncpl must be overridden', &
653  terminate=.true.)
654  end function get_ncpl
655 
656  !> @brief Return volume of cell n based on x value passed.
657  function get_cell_volume(this, n, x)
658  ! -- return
659  real(dp) :: get_cell_volume
660  ! -- dummy
661  class(disbasetype) :: this
662  integer(I4B), intent(in) :: n
663  real(dp), intent(in) :: x
664  ! -- local
665  real(dp) :: tp
666  real(dp) :: bt
667  real(dp) :: sat
668  real(dp) :: thick
669 
671  tp = this%top(n)
672  bt = this%bot(n)
673  sat = squadraticsaturation(tp, bt, x)
674  thick = (tp - bt) * sat
675  get_cell_volume = this%area(n) * thick
676  end function get_cell_volume
677 
678  !> @brief Get a 2D array of polygon vertices, listed in
679  !! clockwise order beginning with the lower left corner.
680  subroutine get_polyverts(this, ic, polyverts, closed)
681  class(disbasetype), intent(inout) :: this
682  integer(I4B), intent(in) :: ic !< cell number (reduced)
683  real(DP), allocatable, intent(out) :: polyverts(:, :) !< polygon vertices (column-major indexing)
684  logical(LGP), intent(in), optional :: closed !< whether to close the polygon, duplicating a vertex
685 
686  errmsg = 'Programmer error: get_polyverts must be overridden'
687  call store_error(errmsg, terminate=.true.)
688  end subroutine
689 
690  !> @brief Read an integer array
691  subroutine read_int_array(this, line, lloc, istart, istop, iout, in, &
692  iarray, aname)
693  ! -- dummy
694  class(disbasetype), intent(inout) :: this
695  character(len=*), intent(inout) :: line
696  integer(I4B), intent(inout) :: lloc
697  integer(I4B), intent(inout) :: istart
698  integer(I4B), intent(inout) :: istop
699  integer(I4B), intent(in) :: in
700  integer(I4B), intent(in) :: iout
701  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: iarray
702  character(len=*), intent(in) :: aname
703 
704  errmsg = 'Programmer error: read_int_array must be overridden'
705  call store_error(errmsg, terminate=.true.)
706  end subroutine read_int_array
707 
708  !> @brief Read a double precision array
709  subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, &
710  darray, aname)
711  ! -- dummy
712  class(disbasetype), intent(inout) :: this
713  character(len=*), intent(inout) :: line
714  integer(I4B), intent(inout) :: lloc
715  integer(I4B), intent(inout) :: istart
716  integer(I4B), intent(inout) :: istop
717  integer(I4B), intent(in) :: in
718  integer(I4B), intent(in) :: iout
719  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray
720  character(len=*), intent(in) :: aname
721 
722  errmsg = 'Programmer error: read_dbl_array must be overridden'
723  call store_error(errmsg, terminate=.true.)
724  end subroutine read_dbl_array
725 
726  !> @brief Fill an integer array
727  subroutine fill_int_array(this, ibuff1, ibuff2)
728  ! -- dummy
729  class(disbasetype), intent(inout) :: this
730  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: ibuff1
731  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: ibuff2
732  ! -- local
733  integer(I4B) :: nodeu
734  integer(I4B) :: noder
735 
736  do nodeu = 1, this%nodesuser
737  noder = this%get_nodenumber(nodeu, 0)
738  if (noder <= 0) cycle
739  ibuff2(noder) = ibuff1(nodeu)
740  end do
741  end subroutine fill_int_array
742 
743  !> @brief Fill a double precision array
744  subroutine fill_dbl_array(this, buff1, buff2)
745  ! -- dummy
746  class(disbasetype), intent(inout) :: this
747  real(DP), dimension(:), pointer, contiguous, intent(in) :: buff1
748  real(DP), dimension(:), pointer, contiguous, intent(inout) :: buff2
749  ! -- local
750  integer(I4B) :: nodeu
751  integer(I4B) :: noder
752 
753  do nodeu = 1, this%nodesuser
754  noder = this%get_nodenumber(nodeu, 0)
755  if (noder <= 0) cycle
756  buff2(noder) = buff1(nodeu)
757  end do
758  end subroutine fill_dbl_array
759 
760  !> @brief Read a list using the list reader.
761  !!
762  !! Convert user node numbers to reduced numbers.
763  !! Terminate if any nodenumbers are within an inactive domain.
764  !! Set up time series and multiply by iauxmultcol if it exists.
765  !! Write the list to iout if iprpak is set.
766  !<
767  subroutine read_list(this, line_reader, in, iout, iprpak, nlist, &
768  inamedbound, iauxmultcol, nodelist, rlist, auxvar, &
769  auxname, boundname, label, pkgname, tsManager, iscloc, &
770  indxconvertflux)
771  ! -- modules
776  use inputoutputmodule, only: urword
779  ! -- dummy
780  class(disbasetype) :: this
781  type(longlinereadertype), intent(inout) :: line_reader
782  integer(I4B), intent(in) :: in
783  integer(I4B), intent(in) :: iout
784  integer(I4B), intent(in) :: iprpak
785  integer(I4B), intent(inout) :: nlist
786  integer(I4B), intent(in) :: inamedbound
787  integer(I4B), intent(in) :: iauxmultcol
788  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: nodelist
789  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: rlist
790  real(DP), dimension(:, :), pointer, contiguous, intent(inout) :: auxvar
791  character(len=LENAUXNAME), dimension(:), intent(inout) :: auxname
792  character(len=LENBOUNDNAME), dimension(:), pointer, contiguous, &
793  intent(inout) :: boundname
794  character(len=*), intent(in) :: label
795  character(len=*), intent(in) :: pkgName
796  type(timeseriesmanagertype) :: tsManager
797  integer(I4B), intent(in) :: iscloc
798  integer(I4B), intent(in), optional :: indxconvertflux
799  ! -- local
800  integer(I4B) :: l
801  integer(I4B) :: nodeu, noder
802  character(len=LINELENGTH) :: nodestr
803  integer(I4B) :: ii, jj
804  real(DP), pointer :: bndElem => null()
805  type(listreadertype) :: lstrdobj
806  type(timeserieslinktype), pointer :: tsLinkBnd => null()
807  type(timeserieslinktype), pointer :: tsLinkAux => null()
808  !
809  ! -- Read the list
810  call lstrdobj%read_list(line_reader, in, iout, nlist, inamedbound, &
811  this%mshape, nodelist, rlist, auxvar, auxname, &
812  boundname, label)
813  !
814  ! -- Go through all locations where a text string was found instead of
815  ! a double precision value and make time-series links to rlist
816  if (lstrdobj%ntxtrlist > 0) then
817  do l = 1, lstrdobj%ntxtrlist
818  ii = lstrdobj%idxtxtrow(l)
819  jj = lstrdobj%idxtxtcol(l)
820  tslinkbnd => null()
821  bndelem => rlist(jj, ii)
822  call read_value_or_time_series(lstrdobj%txtrlist(l), ii, jj, bndelem, &
823  pkgname, 'BND', tsmanager, iprpak, &
824  tslinkbnd)
825  if (associated(tslinkbnd)) then
826  !
827  ! -- If iauxmultcol is active and this column is the column
828  ! to be scaled, then assign tsLinkBnd%RMultiplier to auxvar
829  ! multiplier
830  if (iauxmultcol > 0 .and. jj == iscloc) then
831  tslinkbnd%RMultiplier => auxvar(iauxmultcol, ii)
832  end if
833  !
834  ! -- If boundaries are named, save the name in the link
835  if (lstrdobj%inamedbound == 1) then
836  tslinkbnd%BndName = lstrdobj%boundname(tslinkbnd%IRow)
837  end if
838  !
839  ! -- if the value is a flux and needs to be converted to a flow
840  ! then set the tsLinkBnd appropriately
841  if (present(indxconvertflux)) then
842  if (indxconvertflux == jj) then
843  tslinkbnd%convertflux = .true.
844  nodeu = nodelist(ii)
845  noder = this%get_nodenumber(nodeu, 0)
846  tslinkbnd%CellArea = this%get_area(noder)
847  end if
848  end if
849  !
850  end if
851  end do
852  end if
853  !
854  ! -- Make time-series substitutions for auxvar
855  if (lstrdobj%ntxtauxvar > 0) then
856  do l = 1, lstrdobj%ntxtauxvar
857  ii = lstrdobj%idxtxtauxrow(l)
858  jj = lstrdobj%idxtxtauxcol(l)
859  tslinkaux => null()
860  bndelem => auxvar(jj, ii)
861  call read_value_or_time_series(lstrdobj%txtauxvar(l), ii, jj, bndelem, &
862  pkgname, 'AUX', tsmanager, iprpak, &
863  tslinkaux)
864  if (lstrdobj%inamedbound == 1) then
865  if (associated(tslinkaux)) then
866  tslinkaux%BndName = lstrdobj%boundname(tslinkaux%IRow)
867  end if
868  end if
869  end do
870  end if
871  !
872  ! -- Multiply rlist by the multiplier column in auxvar
873  if (iauxmultcol > 0) then
874  do l = 1, nlist
875  rlist(iscloc, l) = rlist(iscloc, l) * auxvar(iauxmultcol, l)
876  end do
877  end if
878  !
879  ! -- Write the list to iout if requested
880  if (iprpak /= 0) then
881  call lstrdobj%write_list()
882  end if
883  !
884  ! -- Convert user nodenumbers to reduced nodenumbers, if necessary.
885  ! Conversion to reduced nodenumbers must be done last, after the
886  ! list is written so that correct indices are written to the list.
887  if (this%nodes < this%nodesuser) then
888  do l = 1, nlist
889  nodeu = nodelist(l)
890  noder = this%get_nodenumber(nodeu, 0)
891  if (noder <= 0) then
892  call this%nodeu_to_string(nodeu, nodestr)
893  write (errmsg, *) &
894  ' Cell is outside active grid domain: '// &
895  trim(adjustl(nodestr))
896  call store_error(errmsg)
897  end if
898  nodelist(l) = noder
899  end do
900  !
901  ! -- Check for errors and terminate if encountered
902  if (count_errors() > 0) then
903  write (errmsg, *) count_errors(), ' errors encountered.'
904  call store_error(errmsg)
905  call store_error_unit(in)
906  end if
907  end if
908  end subroutine read_list
909 
910  !> @brief Read a 2d double array into col icolbnd of darray.
911  !!
912  !! For cells that are outside of the active domain,
913  !! do not copy the array value into darray.
914  !<
915  subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, &
916  icolbnd, aname, inunit, iout)
917  ! -- dummy
918  class(disbasetype) :: this
919  integer(I4B), intent(in) :: ncolbnd
920  integer(I4B), intent(in) :: maxbnd
921  integer(I4B), dimension(maxbnd) :: nodelist
922  real(DP), dimension(ncolbnd, maxbnd), intent(inout) :: darray
923  integer(I4B), intent(in) :: icolbnd
924  character(len=*), intent(in) :: aname
925  integer(I4B), intent(in) :: inunit
926  integer(I4B), intent(in) :: iout
927 
928  errmsg = 'Programmer error: read_layer_array must be overridden'
929  call store_error(errmsg, terminate=.true.)
930  end subroutine read_layer_array
931 
932  !> @brief Record a double precision array.
933  !!
934  !! The array is written to a formatted or unformatted external file
935  !! depending on the arguments.
936  subroutine record_array(this, darray, iout, iprint, idataun, aname, &
937  cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
938  ! -- dummy
939  class(disbasetype), intent(inout) :: this
940  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray !< double precision array to record
941  integer(I4B), intent(in) :: iout !< ascii output unit number
942  integer(I4B), intent(in) :: iprint !< whether to print the array
943  integer(I4B), intent(in) :: idataun !< binary output unit number
944  character(len=*), intent(in) :: aname !< text descriptor
945  character(len=*), intent(in) :: cdatafmp !< write format
946  integer(I4B), intent(in) :: nvaluesp !< values per line
947  integer(I4B), intent(in) :: nwidthp !< number width
948  character(len=*), intent(in) :: editdesc !< format type (I, G, F, S, E)
949  real(DP), intent(in) :: dinact !< double precision value for cells excluded from model domain
950 
951  errmsg = 'Programmer error: record_array must be overridden'
952  call store_error(errmsg, terminate=.true.)
953  end subroutine record_array
954 
955  !> @brief Record a connection-based double precision array
956  subroutine record_connection_array(this, flowja, ibinun, iout)
957  ! -- dummy
958  class(disbasetype) :: this
959  real(DP), dimension(:), intent(in) :: flowja
960  integer(I4B), intent(in) :: ibinun
961  integer(I4B), intent(in) :: iout
962  ! -- local
963  character(len=16), dimension(1) :: text
964  ! -- data
965  data text(1)/' FLOW-JA-FACE'/
966 
967  ! -- write full ja array
968  call ubdsv1(kstp, kper, text(1), ibinun, flowja, size(flowja), 1, 1, &
969  iout, delt, pertim, totim)
970  end subroutine record_connection_array
971 
972  !> @brief Convert reduced node number to string (nodenumber), (k,j) or (k,i,j)
973  subroutine noder_to_string(this, noder, str)
974  ! -- dummy
975  class(disbasetype) :: this
976  integer(I4B), intent(in) :: noder
977  character(len=*), intent(inout) :: str
978  ! -- local
979  integer(I4B) :: nodeu
980 
981  nodeu = this%get_nodeuser(noder)
982  call this%nodeu_to_string(nodeu, str)
983  end subroutine noder_to_string
984 
985  !> @brief Convert reduced node number to array (nodenumber), (k,j) or (k,i,j)
986  subroutine noder_to_array(this, noder, arr)
987  ! -- dummy
988  class(disbasetype) :: this
989  integer(I4B), intent(in) :: noder
990  integer(I4B), dimension(:), intent(inout) :: arr
991  ! -- local
992  integer(I4B) :: nodeu
993 
994  nodeu = this%get_nodeuser(noder)
995  call this%nodeu_to_array(nodeu, arr)
996  end subroutine noder_to_array
997 
998  !> @brief Record list header for imeth=6
999  subroutine record_srcdst_list_header(this, text, textmodel, textpackage, &
1000  dstmodel, dstpackage, naux, auxtxt, &
1001  ibdchn, nlist, iout)
1002  class(disbasetype) :: this
1003  character(len=16), intent(in) :: text
1004  character(len=16), intent(in) :: textmodel
1005  character(len=16), intent(in) :: textpackage
1006  character(len=16), intent(in) :: dstmodel
1007  character(len=16), intent(in) :: dstpackage
1008  integer(I4B), intent(in) :: naux
1009  character(len=16), dimension(:), intent(in) :: auxtxt
1010  integer(I4B), intent(in) :: ibdchn
1011  integer(I4B), intent(in) :: nlist
1012  integer(I4B), intent(in) :: iout
1013 
1014  errmsg = 'Programmer error: record_srcdst_list_header must be overridden'
1015  call store_error(errmsg, terminate=.true.)
1016  end subroutine record_srcdst_list_header
1017 
1018  !> @brief Record list header
1019  subroutine record_srcdst_list_entry(this, ibdchn, noder, noder2, q, &
1020  naux, aux, olconv, olconv2)
1021  ! -- dummy
1022  class(disbasetype) :: this
1023  integer(I4B), intent(in) :: ibdchn
1024  integer(I4B), intent(in) :: noder
1025  integer(I4B), intent(in) :: noder2
1026  real(DP), intent(in) :: q
1027  integer(I4B), intent(in) :: naux
1028  real(DP), dimension(naux), intent(in) :: aux
1029  logical, optional, intent(in) :: olconv
1030  logical, optional, intent(in) :: olconv2
1031  ! -- local
1032  logical :: lconv
1033  logical :: lconv2
1034  integer(I4B) :: nodeu
1035  integer(I4B) :: nodeu2
1036  !
1037  ! -- Use ubdsvb to write list header
1038  if (present(olconv)) then
1039  lconv = olconv
1040  else
1041  lconv = .true.
1042  end if
1043  if (lconv) then
1044  nodeu = this%get_nodeuser(noder)
1045  else
1046  nodeu = noder
1047  end if
1048  if (present(olconv2)) then
1049  lconv2 = olconv2
1050  else
1051  lconv2 = .true.
1052  end if
1053  if (lconv2) then
1054  nodeu2 = this%get_nodeuser(noder2)
1055  else
1056  nodeu2 = noder2
1057  end if
1058  call ubdsvd(ibdchn, nodeu, nodeu2, q, naux, aux)
1059  end subroutine record_srcdst_list_entry
1060 
1061  !> @brief Convert an integer array to nodelist.
1062  !!
1063  !! For DIS/DISV, the array is layer number, for DISU it's node number.
1064  !<
1065  subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
1066  class(disbasetype) :: this
1067  integer(I4B), intent(in) :: maxbnd
1068  integer(I4B), dimension(:), pointer, contiguous :: darray
1069  integer(I4B), dimension(maxbnd), intent(inout) :: nodelist
1070  integer(I4B), intent(inout) :: nbound
1071  character(len=*), intent(in) :: aname
1072 
1073  errmsg = 'Programmer error: nlarray_to_nodelist must be overridden'
1074  call store_error(errmsg, terminate=.true.)
1075  end subroutine nlarray_to_nodelist
1076 
1077  !> @brief Find the first highest active cell beneath cell n
1078  subroutine highest_active(this, n, ibound)
1079  ! -- dummy
1080  class(disbasetype) :: this
1081  integer(I4B), intent(inout) :: n
1082  integer(I4B), dimension(:), intent(in) :: ibound
1083  ! -- locals
1084  integer(I4B) :: m, ii, iis
1085  logical done, bottomcell
1086  !
1087  ! -- Loop through connected cells until the highest active one (including a
1088  ! constant head cell) is found. Return that cell as n.
1089  done = .false.
1090  do while (.not. done)
1091  bottomcell = .true.
1092  cloop: do ii = this%con%ia(n) + 1, this%con%ia(n + 1) - 1
1093  m = this%con%ja(ii)
1094  iis = this%con%jas(ii)
1095  if (this%con%ihc(iis) == 0 .and. m > n) then
1096  !
1097  ! -- this cannot be a bottom cell
1098  bottomcell = .false.
1099  !
1100  ! -- vertical down
1101  if (ibound(m) /= 0) then
1102  n = m
1103  done = .true.
1104  exit cloop
1105  else
1106  n = m
1107  exit cloop
1108  end if
1109  end if
1110  end do cloop
1111  if (bottomcell) done = .true.
1112  end do
1113  end subroutine highest_active
1114 
1115  !> @brief Find the first saturated cell beneath cell n
1116  subroutine highest_saturated(this, n, sat)
1117  ! -- dummy
1118  class(disbasetype) :: this
1119  integer(I4B), intent(inout) :: n
1120  real(DP), dimension(:), intent(in) :: sat
1121  ! -- locals
1122  integer(I4B) :: m, ii, iis
1123  logical(LGP) :: is_done, bottomcell
1124  !
1125  ! -- Loop through connected cells until the highest saturated one (including a
1126  ! constant head cell) is found. Return that cell as n.
1127  is_done = .false.
1128  do while (.not. is_done)
1129  bottomcell = .true.
1130  cloop: do ii = this%con%ia(n) + 1, this%con%ia(n + 1) - 1
1131  m = this%con%ja(ii)
1132  iis = this%con%jas(ii)
1133  if (this%con%ihc(iis) == 0 .and. m > n) then
1134  !
1135  ! -- this cannot be a bottom cell
1136  bottomcell = .false.
1137  !
1138  ! -- vertical down
1139  if (sat(m) > dzero) then
1140  n = m
1141  is_done = .true.
1142  exit cloop
1143  else
1144  n = m
1145  exit cloop
1146  end if
1147  end if
1148  end do cloop
1149  if (bottomcell) is_done = .true.
1150  end do
1151  end subroutine highest_saturated
1152 
1153  !> @brief Return the cell area for the given node
1154  function get_area(this, node) result(area)
1155  class(disbasetype) :: this
1156  integer(I4B), intent(in) :: node !< reduced node number
1157  real(dp) :: area
1158 
1159  area = this%area(node)
1160  end function get_area
1161 
1162  !> @ brief Calculate the area factor for the cell connection
1163  !!
1164  !! Function calculates the area factor for the cell connection. The sum of
1165  !! all area factors for all cell connections to overlying or underlying
1166  !! cells cells will be 1.
1167  !!
1168  !! TODO: confirm that this works for cells that are only partially covered
1169  !! by overlying or underlying cells.
1170  !<
1171  function get_area_factor(this, node, idx_conn) result(area_factor)
1172  ! -- return
1173  real(dp) :: area_factor !< connection cell area factor
1174  ! -- dummy
1175  class(disbasetype) :: this
1176  integer(I4B), intent(in) :: node !< cell node number
1177  integer(I4B), intent(in) :: idx_conn !< connection index
1178  ! -- local
1179  real(dp) :: area_node
1180  real(dp) :: area_conn
1181  !
1182  ! -- calculate the cell area fraction
1183  area_node = this%area(node)
1184  area_conn = this%con%hwva(idx_conn)
1185  !
1186  ! -- return the cell area factor
1187  area_factor = area_conn / area_node
1188  end function get_area_factor
1189 
1190  !> @ brief Calculate the flow width between two cells
1191  !!
1192  !! This should only be called for connections with IHC > 0.
1193  !! Routine is needed, so it can be overridden by the linear
1194  !! network discretization, which allows for a separate flow
1195  !< width for each cell.
1196  !<
1197  subroutine get_flow_width(this, n, m, idx_conn, width_n, width_m)
1198  ! dummy
1199  class(disbasetype) :: this
1200  integer(I4B), intent(in) :: n !< cell node number
1201  integer(I4B), intent(in) :: m !< cell node number
1202  integer(I4B), intent(in) :: idx_conn !< connection index
1203  real(DP), intent(out) :: width_n !< flow width for cell n
1204  real(DP), intent(out) :: width_m !< flow width for cell m
1205  ! local
1206  integer(I4B) :: isympos
1207 
1208  ! For general case, width_n = width_m
1209  isympos = this%con%jas(idx_conn)
1210  width_n = this%con%hwva(isympos)
1211  width_m = width_n
1212 
1213  end subroutine get_flow_width
1214 
1215  !> @Brief return true if grid is three dimensional
1216  function is_3d(this) result(r)
1217  ! dummy
1218  class(disbasetype) :: this
1219  ! return
1220  logical(LGP) :: r
1221  r = .false.
1222  select case (this%get_dis_enum())
1223  case (dis, disv, disu)
1224  r = .true.
1225  end select
1226  end function is_3d
1227 
1228  !> @Brief return true if grid is two dimensional
1229  function is_2d(this) result(r)
1230  ! dummy
1231  class(disbasetype) :: this
1232  ! return
1233  logical(LGP) :: r
1234  r = .false.
1235  select case (this%get_dis_enum())
1236  case (dis2d, disv2d, disu2d)
1237  r = .true.
1238  end select
1239  end function is_2d
1240 
1241  !> @Brief return true if grid is one dimensional
1242  function is_1d(this) result(r)
1243  ! dummy
1244  class(disbasetype) :: this
1245  ! return
1246  logical(LGP) :: r
1247  r = .false.
1248  select case (this%get_dis_enum())
1249  case (dis1d, disv1d, disu1d)
1250  r = .true.
1251  end select
1252  end function is_1d
1253 
1254 end module basedismodule
real(dp) function get_area(this, node)
Return the cell area for the given node.
subroutine dis_df(this)
Define the discretization.
subroutine noder_to_string(this, noder, str)
Convert reduced node number to string (nodenumber), (k,j) or (k,i,j)
integer(i4b) function get_nodenumber_idx2(this, k, j, icheck)
real(dp) function get_area_factor(this, node, idx_conn)
@ brief Calculate the area factor for the cell connection
subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, icolbnd, aname, inunit, iout)
Read a 2d double array into col icolbnd of darray.
integer(i4b) function nodeu_from_string(this, lloc, istart, istop, in, iout, line, flag_string, allow_zero)
Convert a string to a user nodenumber.
subroutine allocate_scalars(this, name_model, input_mempath)
Allocate and initialize scalar variables.
subroutine dis_mc(this, moffset, idxglo, matrix_sln)
Map cell connections in the numerical solution coefficient matrix.
logical(lgp) function is_1d(this)
@Brief return true if grid is one dimensional
subroutine highest_active(this, n, ibound)
Find the first highest active cell beneath cell n.
subroutine record_array(this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
Record a double precision array.
subroutine record_srcdst_list_header(this, text, textmodel, textpackage, dstmodel, dstpackage, naux, auxtxt, ibdchn, nlist, iout)
Record list header for imeth=6.
integer(i4b) function get_nodenumber_idx3(this, k, i, j, icheck)
subroutine record_srcdst_list_entry(this, ibdchn, noder, noder2, q, naux, aux, olconv, olconv2)
Record list header.
subroutine write_grb(this, icelltype)
Write a binary grid file.
subroutine, public dis_transform_xy(x, y, xorigin, yorigin, angrot, xglo, yglo)
Get global (x, y) coordinates from cell-local coordinates.
subroutine read_int_array(this, line, lloc, istart, istop, iout, in, iarray, aname)
Read an integer array.
logical(lgp) function is_3d(this)
@Brief return true if grid is three dimensional
logical function supports_layers(this)
Indicates whether the grid discretization supports layers.
subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, ipos)
Get normal vector components between the cell and a given neighbor. The normal points outward from th...
subroutine nodeu_to_string(this, nodeu, str)
Convert a user nodenumber to a string (nodenumber), (k,j), or (k,i,j)
integer(i4b) function noder_from_string(this, lloc, istart, istop, in, iout, line, flag_string)
Convert a string to a reduced nodenumber.
subroutine fill_int_array(this, ibuff1, ibuff2)
Fill an integer array.
subroutine get_polyverts(this, ic, polyverts, closed)
Get a 2D array of polygon vertices, listed in clockwise order beginning with the lower left corner.
subroutine dis_da(this)
@brier Deallocate variables
subroutine fill_dbl_array(this, buff1, buff2)
Fill a double precision array.
subroutine dis_ac(this, moffset, sparse)
Add connections to sparse cell connectivity matrix.
subroutine record_connection_array(this, flowja, ibinun, iout)
Record a connection-based double precision array.
subroutine nodeu_to_array(this, nodeu, arr)
Convert a user nodenumber to an array (nodenumber), (k,j), or (k,i,j)
integer(i4b) function get_nodenumber_idx1(this, nodeu, icheck)
integer(i4b) function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, allow_zero)
Convert a cellid string to a user nodenumber.
subroutine highest_saturated(this, n, sat)
Find the first saturated cell beneath cell n.
logical(lgp) function is_2d(this)
@Brief return true if grid is two dimensional
subroutine noder_to_array(this, noder, arr)
Convert reduced node number to array (nodenumber), (k,j) or (k,i,j)
subroutine get_dis_type(this, dis_type)
Get the discretization type (DIS, DISV, or DISU)
integer(i4b) function get_nodeuser(this, noder)
Convert a reduced nodenumber to a user node number.
subroutine get_flow_width(this, n, m, idx_conn, width_n, width_m)
@ brief Calculate the flow width between two cells
subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, darray, aname)
Read a double precision array.
subroutine dis_ar(this, icelltype)
Allocate and setup variables, and write binary grid file.
subroutine allocate_arrays(this)
Allocate and initialize arrays.
integer(i4b) function get_dis_enum(this)
Get the discretization type enumeration.
integer(i4b) function noder_from_cellid(this, cellid, inunit, iout, flag_string, allow_zero)
Convert cellid string to reduced nodenumber.
subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
Convert an integer array to nodelist.
integer(i4b) function get_ncpl(this)
Return number of cells per layer. This is nodes for a DISU grid, as there are no layers.
real(dp) function get_cell_volume(this, n, x)
Return volume of cell n based on x value passed.
subroutine read_list(this, line_reader, in, iout, iprpak, nlist, inamedbound, iauxmultcol, nodelist, rlist, auxvar, auxname, boundname, label, pkgname, tsManager, iscloc, indxconvertflux)
Read a list using the list reader.
subroutine connection_vector(this, noden, nodem, nozee, satn, satm, ihc, xcomp, ycomp, zcomp, conlen)
Get unit vector components between the cell and a given neighbor. Saturation must be provided to comp...
This module contains block parser methods.
Definition: BlockParser.f90:7
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenmodelname
maximum length of the model name
Definition: Constants.f90:22
@ disu
DISV6 discretization.
Definition: Constants.f90:157
@ dis
DIS6 discretization.
Definition: Constants.f90:155
@ dis1d
DIS1D6 discretization.
Definition: Constants.f90:159
@ disv2d
DISV2D6 discretization.
Definition: Constants.f90:164
@ disv1d
DISV1D6 discretization.
Definition: Constants.f90:160
@ dis2d
DIS2D6 discretization.
Definition: Constants.f90:163
@ disv
DISU6 discretization.
Definition: Constants.f90:156
@ disu1d
DISU1D6 discretization.
Definition: Constants.f90:161
@ disundef
undefined discretization
Definition: Constants.f90:153
@ disu2d
DISU2D6 discretization.
Definition: Constants.f90:165
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 dpio180
real constant
Definition: Constants.f90:130
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
subroutine, public ubdsv1(kstp, kper, text, ibdchn, buff, ncol, nrow, nlay, iout, delt, pertim, totim)
Record cell-by-cell flow terms for one component of flow as a 3-D array with extra record to indicate...
subroutine, public ubdsvd(ibdchn, n, n2, q, naux, aux)
Write one value of cell-by-cell flow using a list structure.
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module defines variable data types.
Definition: kind.f90:8
Generic List Reader Module.
Definition: ListReader.f90:3
This module contains the LongLineReaderType.
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
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_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
real(dp) function squadraticsaturation(top, bot, x, eps)
@ brief sQuadraticSaturation
real(dp), pointer, public pertim
time relative to start of stress period
Definition: tdis.f90:30
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
subroutine, public read_value_or_time_series(textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, tsLink)
Call this subroutine if the time-series link is available or needed.