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