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