MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
Dis.f90
Go to the documentation of this file.
1 module dismodule
2 
4  use kindmodule, only: dp, i4b, lgp
5  use constantsmodule, only: linelength, dhalf, done, dzero, &
7  use basedismodule, only: disbasetype
16  use tdismodule, only: kstp, kper, pertim, totim, delt
17 
18  implicit none
19  private
20  public dis_cr, distype
21 
22  !> @brief Structured grid discretization
23  type, extends(disbasetype) :: distype
24  integer(I4B), pointer :: nlay => null() !< number of layers
25  integer(I4B), pointer :: nrow => null() !< number of rows
26  integer(I4B), pointer :: ncol => null() !< number of columns
27  real(dp), dimension(:), pointer, contiguous :: delr => null() !< spacing along a row
28  real(dp), dimension(:), pointer, contiguous :: delc => null() !< spacing along a column
29  real(dp), dimension(:, :), pointer, contiguous :: top2d => null() !< top elevations for each cell at top of model (ncol, nrow)
30  real(dp), dimension(:, :, :), pointer, contiguous :: bot3d => null() !< bottom elevations for each cell (ncol, nrow, nlay)
31  integer(I4B), dimension(:, :, :), pointer, contiguous :: idomain => null() !< idomain (ncol, nrow, nlay)
32  real(dp), dimension(:, :, :), pointer :: botm => null() !< top and bottom elevations for each cell (ncol, nrow, nlay)
33  real(dp), dimension(:), pointer, contiguous :: cellx => null() !< cell center x coordinate for column j
34  real(dp), dimension(:), pointer, contiguous :: celly => null() !< cell center y coordinate for row i
35 
36  contains
37 
38  procedure :: dis_df => dis3d_df
39  procedure :: dis_da => dis3d_da
40  procedure :: get_dis_type => get_dis_type
41  procedure :: get_dis_enum => get_dis_enum
42  procedure, public :: record_array
43  procedure, public :: read_layer_array
44  procedure, public :: record_srcdst_list_header
45  procedure, public :: nlarray_to_nodelist
46  ! -- helper functions
47  procedure :: get_nodenumber_idx1
48  procedure :: get_nodenumber_idx3
49  procedure :: nodeu_to_string
50  procedure :: nodeu_to_array
51  procedure :: nodeu_from_string
52  procedure :: nodeu_from_cellid
53  procedure :: supports_layers
54  procedure :: get_ncpl
55  procedure :: get_polyverts
56  procedure :: connection_vector
57  procedure :: connection_normal
58  ! -- private
59  procedure :: source_options
60  procedure :: source_dimensions
61  procedure :: source_griddata
62  procedure :: log_options
63  procedure :: log_dimensions
64  procedure :: log_griddata
65  procedure :: grid_finalize
66  procedure :: write_grb
67  procedure :: allocate_scalars
68  procedure :: allocate_arrays
69  !
70  ! -- Read a node-sized model array (reduced or not)
71  procedure :: read_int_array
72  procedure :: read_dbl_array
73  end type distype
74 
75  !> @brief Simplifies tracking parameters sourced from the input context.
77  logical :: length_units = .false.
78  logical :: nogrb = .false.
79  logical :: xorigin = .false.
80  logical :: yorigin = .false.
81  logical :: angrot = .false.
82  logical :: nlay = .false.
83  logical :: nrow = .false.
84  logical :: ncol = .false.
85  logical :: delr = .false.
86  logical :: delc = .false.
87  logical :: top = .false.
88  logical :: botm = .false.
89  logical :: idomain = .false.
90  end type disfoundtype
91 
92 contains
93 
94  !> @brief Create a new structured discretization object
95  !<
96  subroutine dis_cr(dis, name_model, input_mempath, inunit, iout)
97  ! -- dummy
98  class(disbasetype), pointer :: dis
99  character(len=*), intent(in) :: name_model
100  character(len=*), intent(in) :: input_mempath
101  integer(I4B), intent(in) :: inunit
102  integer(I4B), intent(in) :: iout
103  ! -- locals
104  type(distype), pointer :: disnew
105  ! -- formats
106  character(len=*), parameter :: fmtheader = &
107  "(1X, /1X, 'DIS -- STRUCTURED GRID DISCRETIZATION PACKAGE,', &
108  &' VERSION 2 : 3/27/2014 - INPUT READ FROM MEMPATH: ', A, /)"
109  !
110  allocate (disnew)
111  dis => disnew
112  call disnew%allocate_scalars(name_model, input_mempath)
113  dis%inunit = inunit
114  dis%iout = iout
115  !
116  ! -- If dis enabled
117  if (inunit > 0) then
118  !
119  ! -- Identify package
120  if (iout > 0) then
121  write (iout, fmtheader) dis%input_mempath
122  end if
123  end if
124  !
125  end subroutine dis_cr
126 
127  !> @brief Define the discretization
128  !<
129  subroutine dis3d_df(this)
130  ! -- dummy
131  class(distype) :: this
132  !
133  ! -- Transfer the data from the memory manager into this package object
134  if (this%inunit /= 0) then
135  !
136  ! -- source input options
137  call this%source_options()
138  !
139  ! -- source input dimensions
140  call this%source_dimensions()
141  !
142  ! -- source input griddata
143  call this%source_griddata()
144  end if
145  !
146  ! -- Final grid initialization
147  call this%grid_finalize()
148  !
149  end subroutine dis3d_df
150 
151  !> @brief Deallocate variables
152  !<
153  subroutine dis3d_da(this)
154  ! -- dummy
155  class(distype) :: this
156  !
157  ! -- Deallocate idm memory
158  call memorylist_remove(this%name_model, 'DIS', idm_context)
159  !
160  ! -- DisBaseType deallocate
161  call this%DisBaseType%dis_da()
162  !
163  ! -- Deallocate scalars
164  call mem_deallocate(this%nlay)
165  call mem_deallocate(this%nrow)
166  call mem_deallocate(this%ncol)
167  call mem_deallocate(this%delr)
168  call mem_deallocate(this%delc)
169  call mem_deallocate(this%cellx)
170  call mem_deallocate(this%celly)
171  !
172  ! -- Deallocate Arrays
173  call mem_deallocate(this%nodereduced)
174  call mem_deallocate(this%nodeuser)
175  call mem_deallocate(this%top2d)
176  call mem_deallocate(this%bot3d)
177  call mem_deallocate(this%idomain)
178  !
179  end subroutine dis3d_da
180 
181  !> @brief Copy options from IDM into package
182  !<
183  subroutine source_options(this)
184  ! -- dummy
185  class(distype) :: this
186  ! -- locals
187  character(len=LENVARNAME), dimension(3) :: lenunits = &
188  &[character(len=LENVARNAME) :: 'FEET', 'METERS', 'CENTIMETERS']
189  type(disfoundtype) :: found
190  !
191  ! -- update defaults with idm sourced values
192  call mem_set_value(this%lenuni, 'LENGTH_UNITS', this%input_mempath, &
193  lenunits, found%length_units)
194  call mem_set_value(this%nogrb, 'NOGRB', this%input_mempath, found%nogrb)
195  call mem_set_value(this%xorigin, 'XORIGIN', this%input_mempath, found%xorigin)
196  call mem_set_value(this%yorigin, 'YORIGIN', this%input_mempath, found%yorigin)
197  call mem_set_value(this%angrot, 'ANGROT', this%input_mempath, found%angrot)
198  !
199  ! -- log values to list file
200  if (this%iout > 0) then
201  call this%log_options(found)
202  end if
203  !
204  end subroutine source_options
205 
206  !> @brief Write user options to list file
207  !<
208  subroutine log_options(this, found)
209  ! -- dummy
210  class(distype) :: this
211  type(disfoundtype), intent(in) :: found
212  !
213  write (this%iout, '(1x,a)') 'Setting Discretization Options'
214  !
215  if (found%length_units) then
216  write (this%iout, '(4x,a,i0)') 'Model length unit [0=UND, 1=FEET, &
217  &2=METERS, 3=CENTIMETERS] set as ', this%lenuni
218  end if
219  !
220  if (found%nogrb) then
221  write (this%iout, '(4x,a,i0)') 'Binary grid file [0=GRB, 1=NOGRB] &
222  &set as ', this%nogrb
223  end if
224  !
225  if (found%xorigin) then
226  write (this%iout, '(4x,a,G0)') 'XORIGIN = ', this%xorigin
227  end if
228  !
229  if (found%yorigin) then
230  write (this%iout, '(4x,a,G0)') 'YORIGIN = ', this%yorigin
231  end if
232  !
233  if (found%angrot) then
234  write (this%iout, '(4x,a,G0)') 'ANGROT = ', this%angrot
235  end if
236  !
237  write (this%iout, '(1x,a,/)') 'End Setting Discretization Options'
238  !
239  end subroutine log_options
240 
241  !> @brief Copy dimensions from IDM into package
242  !<
243  subroutine source_dimensions(this)
244  ! -- dummy
245  class(distype) :: this
246  ! -- locals
247  integer(I4B) :: i, j, k
248  type(disfoundtype) :: found
249  !
250  ! -- update defaults with idm sourced values
251  call mem_set_value(this%nlay, 'NLAY', this%input_mempath, found%nlay)
252  call mem_set_value(this%nrow, 'NROW', this%input_mempath, found%nrow)
253  call mem_set_value(this%ncol, 'NCOL', this%input_mempath, found%ncol)
254  !
255  ! -- log simulation values
256  if (this%iout > 0) then
257  call this%log_dimensions(found)
258  end if
259  !
260  ! -- verify dimensions were set
261  if (this%nlay < 1) then
262  call store_error( &
263  'NLAY was not specified or was specified incorrectly.')
264  call store_error_filename(this%input_fname)
265  end if
266  if (this%nrow < 1) then
267  call store_error( &
268  'NROW was not specified or was specified incorrectly.')
269  call store_error_filename(this%input_fname)
270  end if
271  if (this%ncol < 1) then
272  call store_error( &
273  'NCOL was not specified or was specified incorrectly.')
274  call store_error_filename(this%input_fname)
275  end if
276  !
277  ! -- calculate nodesuser
278  this%nodesuser = this%nlay * this%nrow * this%ncol
279  !
280  ! -- Allocate delr, delc, and non-reduced vectors for dis
281  call mem_allocate(this%delr, this%ncol, 'DELR', this%memoryPath)
282  call mem_allocate(this%delc, this%nrow, 'DELC', this%memoryPath)
283  call mem_allocate(this%idomain, this%ncol, this%nrow, this%nlay, 'IDOMAIN', &
284  this%memoryPath)
285  call mem_allocate(this%top2d, this%ncol, this%nrow, 'TOP2D', this%memoryPath)
286  call mem_allocate(this%bot3d, this%ncol, this%nrow, this%nlay, 'BOT3D', &
287  this%memoryPath)
288  call mem_allocate(this%cellx, this%ncol, 'CELLX', this%memoryPath)
289  call mem_allocate(this%celly, this%nrow, 'CELLY', this%memoryPath)
290  !
291  ! -- initialize all cells to be active (idomain = 1)
292  do k = 1, this%nlay
293  do i = 1, this%nrow
294  do j = 1, this%ncol
295  this%idomain(j, i, k) = 1
296  end do
297  end do
298  end do
299  !
300  end subroutine source_dimensions
301 
302  !> @brief Write dimensions to list file
303  !<
304  subroutine log_dimensions(this, found)
305  ! -- dummy
306  class(distype) :: this
307  type(disfoundtype), intent(in) :: found
308  !
309  write (this%iout, '(1x,a)') 'Setting Discretization Dimensions'
310  !
311  if (found%nlay) then
312  write (this%iout, '(4x,a,i0)') 'NLAY = ', this%nlay
313  end if
314  !
315  if (found%nrow) then
316  write (this%iout, '(4x,a,i0)') 'NROW = ', this%nrow
317  end if
318  !
319  if (found%ncol) then
320  write (this%iout, '(4x,a,i0)') 'NCOL = ', this%ncol
321  end if
322  !
323  write (this%iout, '(1x,a,/)') 'End Setting Discretization Dimensions'
324  !
325  end subroutine log_dimensions
326 
327  !> @brief Copy grid data from IDM into package
328  !<
329  subroutine source_griddata(this)
330  ! -- dummy
331  class(distype) :: this
332  type(disfoundtype) :: found
333  !
334  ! -- update defaults with idm sourced values
335  call mem_set_value(this%delr, 'DELR', this%input_mempath, found%delr)
336  call mem_set_value(this%delc, 'DELC', this%input_mempath, found%delc)
337  call mem_set_value(this%top2d, 'TOP', this%input_mempath, found%top)
338  call mem_set_value(this%bot3d, 'BOTM', this%input_mempath, found%botm)
339  call mem_set_value(this%idomain, 'IDOMAIN', this%input_mempath, found%idomain)
340  !
341  ! -- log simulation values
342  if (this%iout > 0) then
343  call this%log_griddata(found)
344  end if
345  !
346  end subroutine source_griddata
347 
348  !> @brief Write dimensions to list file
349  !<
350  subroutine log_griddata(this, found)
351  ! -- dummy
352  class(distype) :: this
353  type(disfoundtype), intent(in) :: found
354  !
355  write (this%iout, '(1x,a)') 'Setting Discretization Griddata'
356  !
357  if (found%delr) then
358  write (this%iout, '(4x,a)') 'DELR set from input file'
359  end if
360  !
361  if (found%delc) then
362  write (this%iout, '(4x,a)') 'DELC set from input file'
363  end if
364  !
365  if (found%top) then
366  write (this%iout, '(4x,a)') 'TOP set from input file'
367  end if
368  !
369  if (found%botm) then
370  write (this%iout, '(4x,a)') 'BOTM set from input file'
371  end if
372  !
373  if (found%idomain) then
374  write (this%iout, '(4x,a)') 'IDOMAIN set from input file'
375  end if
376  !
377  write (this%iout, '(1x,a,/)') 'End Setting Discretization Griddata'
378  !
379  end subroutine log_griddata
380 
381  !> @brief Finalize grid (check properties, allocate arrays, compute connections)
382  !<
383  subroutine grid_finalize(this)
384  ! -- modules
386  ! -- dummy
387  class(distype) :: this
388  ! -- locals
389  integer(I4B) :: n, i, j, k
390  integer(I4B) :: node
391  integer(I4B) :: noder
392  integer(I4B) :: nrsize
393  real(DP) :: top
394  real(DP) :: dz
395  ! -- formats
396  character(len=*), parameter :: fmtdz = &
397  "('CELL (',i0,',',i0,',',i0,') THICKNESS <= 0. ', &
398  &'TOP, BOT: ',2(1pg24.15))"
399  character(len=*), parameter :: fmtnr = &
400  "(/1x, 'The specified IDOMAIN results in a reduced number of cells.',&
401  &/1x, 'Number of user nodes: ',I0,&
402  &/1X, 'Number of nodes in solution: ', I0, //)"
403  !
404  ! -- count active cells
405  this%nodes = 0
406  do k = 1, this%nlay
407  do i = 1, this%nrow
408  do j = 1, this%ncol
409  if (this%idomain(j, i, k) > 0) this%nodes = this%nodes + 1
410  end do
411  end do
412  end do
413  !
414  ! -- Check to make sure nodes is a valid number
415  if (this%nodes == 0) then
416  call store_error('Model does not have any active nodes. &
417  &Ensure IDOMAIN array has some values greater &
418  &than zero.')
419  call store_error_filename(this%input_fname)
420  end if
421  !
422  ! -- Check cell thicknesses
423  n = 0
424  do k = 1, this%nlay
425  do i = 1, this%nrow
426  do j = 1, this%ncol
427  if (this%idomain(j, i, k) < 1) cycle
428  if (k > 1) then
429  top = this%bot3d(j, i, k - 1)
430  else
431  top = this%top2d(j, i)
432  end if
433  dz = top - this%bot3d(j, i, k)
434  if (dz <= dzero) then
435  n = n + 1
436  write (errmsg, fmt=fmtdz) k, i, j, top, this%bot3d(j, i, k)
437  call store_error(errmsg)
438  end if
439  end do
440  end do
441  end do
442  if (count_errors() > 0) then
443  call store_error_filename(this%input_fname)
444  end if
445  !
446  ! -- Write message if reduced grid
447  if (this%nodes < this%nodesuser) then
448  write (this%iout, fmtnr) this%nodesuser, this%nodes
449  end if
450  !
451  ! -- Array size is now known, so allocate
452  call this%allocate_arrays()
453  !
454  ! -- Fill the nodereduced array with the reduced nodenumber, or
455  ! a negative number to indicate it is a pass-through cell, or
456  ! a zero to indicate that the cell is excluded from the
457  ! solution.
458  if (this%nodes < this%nodesuser) then
459  node = 1
460  noder = 1
461  do k = 1, this%nlay
462  do i = 1, this%nrow
463  do j = 1, this%ncol
464  if (this%idomain(j, i, k) > 0) then
465  this%nodereduced(node) = noder
466  noder = noder + 1
467  elseif (this%idomain(j, i, k) < 0) then
468  this%nodereduced(node) = -1
469  else
470  this%nodereduced(node) = 0
471  end if
472  node = node + 1
473  end do
474  end do
475  end do
476  end if
477  !
478  ! -- allocate and fill nodeuser if a reduced grid
479  if (this%nodes < this%nodesuser) then
480  node = 1
481  noder = 1
482  do k = 1, this%nlay
483  do i = 1, this%nrow
484  do j = 1, this%ncol
485  if (this%idomain(j, i, k) > 0) then
486  this%nodeuser(noder) = node
487  noder = noder + 1
488  end if
489  node = node + 1
490  end do
491  end do
492  end do
493  end if
494  !
495  ! -- fill x,y coordinate arrays
496  this%cellx(1) = dhalf * this%delr(1)
497  this%celly(this%nrow) = dhalf * this%delc(this%nrow)
498  do j = 2, this%ncol
499  this%cellx(j) = this%cellx(j - 1) + dhalf * this%delr(j - 1) + &
500  dhalf * this%delr(j)
501  end do
502  ! -- row number increases in negative y direction:
503  do i = this%nrow - 1, 1, -1
504  this%celly(i) = this%celly(i + 1) + dhalf * this%delc(i + 1) + &
505  dhalf * this%delc(i)
506  end do
507  !
508  ! -- Move top2d and botm3d into top and bot, and calculate area
509  node = 0
510  do k = 1, this%nlay
511  do i = 1, this%nrow
512  do j = 1, this%ncol
513  node = node + 1
514  noder = node
515  if (this%nodes < this%nodesuser) noder = this%nodereduced(node)
516  if (noder <= 0) cycle
517  if (k > 1) then
518  top = this%bot3d(j, i, k - 1)
519  else
520  top = this%top2d(j, i)
521  end if
522  this%top(noder) = top
523  this%bot(noder) = this%bot3d(j, i, k)
524  this%area(noder) = this%delr(j) * this%delc(i)
525  this%xc(noder) = this%cellx(j)
526  this%yc(noder) = this%celly(i)
527  end do
528  end do
529  end do
530  !
531  ! -- create and fill the connections object
532  nrsize = 0
533  if (this%nodes < this%nodesuser) nrsize = this%nodes
534  allocate (this%con)
535  call this%con%disconnections(this%name_model, this%nodes, &
536  this%ncol, this%nrow, this%nlay, &
537  nrsize, this%delr, this%delc, &
538  this%top, this%bot, this%nodereduced, &
539  this%nodeuser)
540  this%nja = this%con%nja
541  this%njas = this%con%njas
542  !
543  end subroutine grid_finalize
544 
545  !> @brief Write a binary grid file
546  !<
547  subroutine write_grb(this, icelltype)
548  ! -- modules
549  use openspecmodule, only: access, form
550  ! -- dummy
551  class(distype) :: this
552  integer(I4B), dimension(:), intent(in) :: icelltype
553  ! -- local
554  integer(I4B) :: iunit, ntxt, ncpl
555  integer(I4B), parameter :: lentxt = 100
556  character(len=50) :: txthdr
557  character(len=lentxt) :: txt
558  character(len=LINELENGTH) :: fname
559  character(len=*), parameter :: fmtgrdsave = &
560  "(4X,'BINARY GRID INFORMATION WILL BE WRITTEN TO:', &
561  &/,6X,'UNIT NUMBER: ', I0,/,6X, 'FILE NAME: ', A)"
562  !
563  ! -- Initialize
564  ntxt = 16
565  ncpl = this%nrow * this%ncol
566  !
567  ! -- Open the file
568  fname = trim(this%input_fname)//'.grb'
569  iunit = getunit()
570  write (this%iout, fmtgrdsave) iunit, trim(adjustl(fname))
571  call openfile(iunit, this%iout, trim(adjustl(fname)), 'DATA(BINARY)', &
572  form, access, 'REPLACE')
573  !
574  ! -- write header information
575  write (txthdr, '(a)') 'GRID DIS'
576  txthdr(50:50) = new_line('a')
577  write (iunit) txthdr
578  write (txthdr, '(a)') 'VERSION 1'
579  txthdr(50:50) = new_line('a')
580  write (iunit) txthdr
581  write (txthdr, '(a, i0)') 'NTXT ', ntxt
582  txthdr(50:50) = new_line('a')
583  write (iunit) txthdr
584  write (txthdr, '(a, i0)') 'LENTXT ', lentxt
585  txthdr(50:50) = new_line('a')
586  write (iunit) txthdr
587  !
588  ! -- write variable definitions
589  write (txt, '(3a, i0)') 'NCELLS ', 'INTEGER ', 'NDIM 0 # ', this%nodesuser
590  txt(lentxt:lentxt) = new_line('a')
591  write (iunit) txt
592  write (txt, '(3a, i0)') 'NLAY ', 'INTEGER ', 'NDIM 0 # ', this%nlay
593  txt(lentxt:lentxt) = new_line('a')
594  write (iunit) txt
595  write (txt, '(3a, i0)') 'NROW ', 'INTEGER ', 'NDIM 0 # ', this%nrow
596  txt(lentxt:lentxt) = new_line('a')
597  write (iunit) txt
598  write (txt, '(3a, i0)') 'NCOL ', 'INTEGER ', 'NDIM 0 # ', this%ncol
599  txt(lentxt:lentxt) = new_line('a')
600  write (iunit) txt
601  write (txt, '(3a, i0)') 'NJA ', 'INTEGER ', 'NDIM 0 # ', this%nja
602  txt(lentxt:lentxt) = new_line('a')
603  write (iunit) txt
604  write (txt, '(3a, 1pg24.15)') 'XORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%xorigin
605  txt(lentxt:lentxt) = new_line('a')
606  write (iunit) txt
607  write (txt, '(3a, 1pg24.15)') 'YORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%yorigin
608  txt(lentxt:lentxt) = new_line('a')
609  write (iunit) txt
610  write (txt, '(3a, 1pg24.15)') 'ANGROT ', 'DOUBLE ', 'NDIM 0 # ', this%angrot
611  txt(lentxt:lentxt) = new_line('a')
612  write (iunit) txt
613  write (txt, '(3a, i0)') 'DELR ', 'DOUBLE ', 'NDIM 1 ', this%ncol
614  txt(lentxt:lentxt) = new_line('a')
615  write (iunit) txt
616  write (txt, '(3a, i0)') 'DELC ', 'DOUBLE ', 'NDIM 1 ', this%nrow
617  txt(lentxt:lentxt) = new_line('a')
618  write (iunit) txt
619  write (txt, '(3a, i0)') 'TOP ', 'DOUBLE ', 'NDIM 1 ', ncpl
620  txt(lentxt:lentxt) = new_line('a')
621  write (iunit) txt
622  write (txt, '(3a, i0)') 'BOTM ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
623  txt(lentxt:lentxt) = new_line('a')
624  write (iunit) txt
625  write (txt, '(3a, i0)') 'IA ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
626  txt(lentxt:lentxt) = new_line('a')
627  write (iunit) txt
628  write (txt, '(3a, i0)') 'JA ', 'INTEGER ', 'NDIM 1 ', size(this%con%jausr)
629  txt(lentxt:lentxt) = new_line('a')
630  write (iunit) txt
631  write (txt, '(3a, i0)') 'IDOMAIN ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
632  txt(lentxt:lentxt) = new_line('a')
633  write (iunit) txt
634  write (txt, '(3a, i0)') 'ICELLTYPE ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
635  txt(lentxt:lentxt) = new_line('a')
636  write (iunit) txt
637  !
638  ! -- write data
639  write (iunit) this%nodesuser ! ncells
640  write (iunit) this%nlay ! nlay
641  write (iunit) this%nrow ! nrow
642  write (iunit) this%ncol ! ncol
643  write (iunit) this%nja ! nja
644  write (iunit) this%xorigin ! xorigin
645  write (iunit) this%yorigin ! yorigin
646  write (iunit) this%angrot ! angrot
647  write (iunit) this%delr ! delr
648  write (iunit) this%delc ! delc
649  write (iunit) this%top2d ! top2d
650  write (iunit) this%bot3d ! bot3d
651  write (iunit) this%con%iausr ! iausr
652  write (iunit) this%con%jausr ! jausr
653  write (iunit) this%idomain ! idomain
654  write (iunit) icelltype ! icelltype
655  !
656  ! -- Close the file
657  close (iunit)
658  !
659  end subroutine write_grb
660 
661  !> @brief Convert a user nodenumber to a string (nodenumber) or (k,i,j)
662  !<
663  subroutine nodeu_to_string(this, nodeu, str)
664  ! -- dummy
665  class(distype) :: this
666  integer(I4B), intent(in) :: nodeu
667  character(len=*), intent(inout) :: str
668  ! -- local
669  integer(I4B) :: i, j, k
670  character(len=10) :: kstr, istr, jstr
671  !
672  call get_ijk(nodeu, this%nrow, this%ncol, this%nlay, i, j, k)
673  write (kstr, '(i10)') k
674  write (istr, '(i10)') i
675  write (jstr, '(i10)') j
676  str = '('//trim(adjustl(kstr))//','// &
677  trim(adjustl(istr))//','// &
678  trim(adjustl(jstr))//')'
679  !
680  end subroutine nodeu_to_string
681 
682  !> @brief Convert a user nodenumber to an array (nodenumber) or (k,i,j)
683  !<
684  subroutine nodeu_to_array(this, nodeu, arr)
685  ! -- dummy
686  class(distype) :: this
687  integer(I4B), intent(in) :: nodeu
688  integer(I4B), dimension(:), intent(inout) :: arr
689  ! -- local
690  integer(I4B) :: isize
691  integer(I4B) :: i, j, k
692  !
693  ! -- check the size of arr
694  isize = size(arr)
695  if (isize /= this%ndim) then
696  write (errmsg, '(a,i0,a,i0,a)') &
697  'Program error: nodeu_to_array size of array (', isize, &
698  ') is not equal to the discretization dimension (', this%ndim, ')'
699  call store_error(errmsg, terminate=.true.)
700  end if
701  !
702  ! -- get k, i, j
703  call get_ijk(nodeu, this%nrow, this%ncol, this%nlay, i, j, k)
704  !
705  ! -- fill array
706  arr(1) = k
707  arr(2) = i
708  arr(3) = j
709  !
710  end subroutine nodeu_to_array
711 
712  !> @brief Get reduced node number from user node number
713  !<
714  function get_nodenumber_idx1(this, nodeu, icheck) result(nodenumber)
715  ! -- return
716  integer(I4B) :: nodenumber
717  ! -- dummy
718  class(distype), intent(in) :: this
719  integer(I4B), intent(in) :: nodeu
720  integer(I4B), intent(in) :: icheck
721  !
722  ! -- check the node number if requested
723  if (icheck /= 0) then
724  !
725  ! -- If within valid range, convert to reduced nodenumber
726  if (nodeu < 1 .or. nodeu > this%nodesuser) then
727  write (errmsg, '(a,i0,a)') &
728  'Node number (', nodeu, &
729  ') less than 1 or greater than the number of nodes.'
730  call store_error(errmsg)
731  nodenumber = 0
732  else
733  nodenumber = nodeu
734  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
735  end if
736  else
737  nodenumber = nodeu
738  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
739  end if
740  !
741  end function get_nodenumber_idx1
742 
743  !> @brief Get reduced node number from layer, row and column indices
744  !<
745  function get_nodenumber_idx3(this, k, i, j, icheck) result(nodenumber)
746  ! -- return
747  integer(I4B) :: nodenumber
748  ! -- dummy
749  class(distype), intent(in) :: this
750  integer(I4B), intent(in) :: k, i, j
751  integer(I4B), intent(in) :: icheck
752  ! -- local
753  integer(I4B) :: nodeu
754  ! formats
755  character(len=*), parameter :: fmterr = &
756  "('Error in structured-grid cell indices: layer = ',i0,', &
757  &row = ',i0,', column = ',i0)"
758  !
759  nodeu = get_node(k, i, j, this%nlay, this%nrow, this%ncol)
760  if (nodeu < 1) then
761  write (errmsg, fmterr) k, i, j
762  call store_error(errmsg, terminate=.true.)
763  end if
764  nodenumber = nodeu
765  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
766  !
767  ! -- check the node number if requested
768  if (icheck /= 0) then
769  !
770  if (k < 1 .or. k > this%nlay) &
771  call store_error('Layer less than one or greater than nlay')
772  if (i < 1 .or. i > this%nrow) &
773  call store_error('Row less than one or greater than nrow')
774  if (j < 1 .or. j > this%ncol) &
775  call store_error('Column less than one or greater than ncol')
776  !
777  ! -- Error if outside of range
778  if (nodeu < 1 .or. nodeu > this%nodesuser) then
779  write (errmsg, '(a,i0,a)') &
780  'Node number (', nodeu, ')less than 1 or greater than nodes.'
781  call store_error(errmsg)
782  end if
783  end if
784  !
785  end function get_nodenumber_idx3
786 
787  !> @brief Allocate and initialize scalar variables
788  !<
789  subroutine allocate_scalars(this, name_model, input_mempath)
790  ! -- dummy
791  class(distype) :: this
792  character(len=*), intent(in) :: name_model
793  character(len=*), intent(in) :: input_mempath
794  !
795  ! -- Allocate parent scalars
796  call this%DisBaseType%allocate_scalars(name_model, input_mempath)
797  !
798  ! -- Allocate
799  call mem_allocate(this%nlay, 'NLAY', this%memoryPath)
800  call mem_allocate(this%nrow, 'NROW', this%memoryPath)
801  call mem_allocate(this%ncol, 'NCOL', this%memoryPath)
802  !
803  ! -- Initialize
804  this%nlay = 0
805  this%nrow = 0
806  this%ncol = 0
807  this%ndim = 3
808  !
809  end subroutine allocate_scalars
810 
811  !> @brief Allocate and initialize arrays
812  !<
813  subroutine allocate_arrays(this)
814  ! -- dummy
815  class(distype) :: this
816  !
817  ! -- Allocate arrays in DisBaseType (mshape, top, bot, area)
818  call this%DisBaseType%allocate_arrays()
819  !
820  ! -- Allocate arrays for DisType
821  if (this%nodes < this%nodesuser) then
822  call mem_allocate(this%nodeuser, this%nodes, 'NODEUSER', this%memoryPath)
823  call mem_allocate(this%nodereduced, this%nodesuser, 'NODEREDUCED', &
824  this%memoryPath)
825  else
826  call mem_allocate(this%nodeuser, 1, 'NODEUSER', this%memoryPath)
827  call mem_allocate(this%nodereduced, 1, 'NODEREDUCED', this%memoryPath)
828  end if
829  !
830  ! -- Initialize
831  this%mshape(1) = this%nlay
832  this%mshape(2) = this%nrow
833  this%mshape(3) = this%ncol
834  !
835  end subroutine allocate_arrays
836 
837  !> @brief Convert a string to a user nodenumber
838  !!
839  !! Parse layer, row and column and return user nodenumber.
840  !! If flag_string is present and true, the first token may be
841  !! non-numeric (e.g. boundary name). In this case, return -2.
842  !<
843  function nodeu_from_string(this, lloc, istart, istop, in, iout, line, &
844  flag_string, allow_zero) result(nodeu)
845  ! -- dummy
846  class(distype) :: this
847  integer(I4B), intent(inout) :: lloc
848  integer(I4B), intent(inout) :: istart
849  integer(I4B), intent(inout) :: istop
850  integer(I4B), intent(in) :: in
851  integer(I4B), intent(in) :: iout
852  character(len=*), intent(inout) :: line
853  logical, optional, intent(in) :: flag_string
854  logical, optional, intent(in) :: allow_zero
855  integer(I4B) :: nodeu
856  ! -- local
857  integer(I4B) :: k, i, j, nlay, nrow, ncol
858  integer(I4B) :: lloclocal, ndum, istat, n
859  real(dp) :: r
860  !
861  if (present(flag_string)) then
862  if (flag_string) then
863  ! Check to see if first token in line can be read as an integer.
864  lloclocal = lloc
865  call urword(line, lloclocal, istart, istop, 1, ndum, r, iout, in)
866  read (line(istart:istop), *, iostat=istat) n
867  if (istat /= 0) then
868  ! First token in line is not an integer; return flag to this effect.
869  nodeu = -2
870  return
871  end if
872  end if
873  end if
874  !
875  nlay = this%mshape(1)
876  nrow = this%mshape(2)
877  ncol = this%mshape(3)
878  !
879  call urword(line, lloc, istart, istop, 2, k, r, iout, in)
880  call urword(line, lloc, istart, istop, 2, i, r, iout, in)
881  call urword(line, lloc, istart, istop, 2, j, r, iout, in)
882  !
883  if (k == 0 .and. i == 0 .and. j == 0) then
884  if (present(allow_zero)) then
885  if (allow_zero) then
886  nodeu = 0
887  return
888  end if
889  end if
890  end if
891  !
892  errmsg = ""
893  !
894  if (k < 1 .or. k > nlay) then
895  write (errmsg, '(a,i0,a)') &
896  'Layer number in list (', k, ') is outside of the grid.'
897  end if
898  if (i < 1 .or. i > nrow) then
899  write (errmsg, '(a,1x,a,i0,a)') &
900  trim(adjustl(errmsg)), 'Row number in list (', i, &
901  ') is outside of the grid.'
902  end if
903  if (j < 1 .or. j > ncol) then
904  write (errmsg, '(a,1x,a,i0,a)') &
905  trim(adjustl(errmsg)), 'Column number in list (', j, &
906  ') is outside of the grid.'
907  end if
908  !
909  nodeu = get_node(k, i, j, nlay, nrow, ncol)
910  !
911  if (nodeu < 1 .or. nodeu > this%nodesuser) then
912  write (errmsg, '(a,1x,a,i0,a)') &
913  trim(adjustl(errmsg)), &
914  "Node number in list (", nodeu, ") is outside of the grid. "// &
915  "Cell number cannot be determined in line '"// &
916  trim(adjustl(line))//"'."
917  end if
918  !
919  if (len_trim(adjustl(errmsg)) > 0) then
920  call store_error(errmsg)
921  call store_error_unit(in)
922  end if
923  !
924  end function nodeu_from_string
925 
926  !> @brief Convert a cellid string to a user nodenumber
927  !!
928  !! If flag_string is present and true, the first token may be
929  !! non-numeric (e.g. boundary name). In this case, return -2.
930  !!
931  !! If allow_zero is present and true, and all indices are zero, the
932  !! result can be zero. If allow_zero is false, a zero in any index is an error.
933  !<
934  function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, &
935  allow_zero) result(nodeu)
936  ! -- return
937  integer(I4B) :: nodeu
938  ! -- dummy
939  class(distype) :: this
940  character(len=*), intent(inout) :: cellid
941  integer(I4B), intent(in) :: inunit
942  integer(I4B), intent(in) :: iout
943  logical, optional, intent(in) :: flag_string
944  logical, optional, intent(in) :: allow_zero
945  ! -- local
946  integer(I4B) :: lloclocal, istart, istop, ndum, n
947  integer(I4B) :: k, i, j, nlay, nrow, ncol
948  integer(I4B) :: istat
949  real(dp) :: r
950  !
951  if (present(flag_string)) then
952  if (flag_string) then
953  ! Check to see if first token in cellid can be read as an integer.
954  lloclocal = 1
955  call urword(cellid, lloclocal, istart, istop, 1, ndum, r, iout, inunit)
956  read (cellid(istart:istop), *, iostat=istat) n
957  if (istat /= 0) then
958  ! First token in cellid is not an integer; return flag to this effect.
959  nodeu = -2
960  return
961  end if
962  end if
963  end if
964  !
965  nlay = this%mshape(1)
966  nrow = this%mshape(2)
967  ncol = this%mshape(3)
968  !
969  lloclocal = 1
970  call urword(cellid, lloclocal, istart, istop, 2, k, r, iout, inunit)
971  call urword(cellid, lloclocal, istart, istop, 2, i, r, iout, inunit)
972  call urword(cellid, lloclocal, istart, istop, 2, j, r, iout, inunit)
973  !
974  if (k == 0 .and. i == 0 .and. j == 0) then
975  if (present(allow_zero)) then
976  if (allow_zero) then
977  nodeu = 0
978  return
979  end if
980  end if
981  end if
982  !
983  errmsg = ""
984  !
985  if (k < 1 .or. k > nlay) then
986  write (errmsg, '(a,i0,a)') &
987  'Layer number in list (', k, ') is outside of the grid.'
988  end if
989  if (i < 1 .or. i > nrow) then
990  write (errmsg, '(a,1x,a,i0,a)') &
991  trim(adjustl(errmsg)), 'Row number in list (', i, &
992  ') is outside of the grid.'
993  end if
994  if (j < 1 .or. j > ncol) then
995  write (errmsg, '(a,1x,a,i0,a)') &
996  trim(adjustl(errmsg)), 'Column number in list (', j, &
997  ') is outside of the grid.'
998  end if
999  !
1000  nodeu = get_node(k, i, j, nlay, nrow, ncol)
1001  !
1002  if (nodeu < 1 .or. nodeu > this%nodesuser) then
1003  write (errmsg, '(a,1x,a,i0,a)') &
1004  trim(adjustl(errmsg)), &
1005  "Cell number cannot be determined for cellid ("// &
1006  trim(adjustl(cellid))//") and results in a user "// &
1007  "node number (", nodeu, ") that is outside of the grid."
1008  end if
1009  !
1010  if (len_trim(adjustl(errmsg)) > 0) then
1011  call store_error(errmsg)
1012  call store_error_unit(inunit)
1013  end if
1014  !
1015  end function nodeu_from_cellid
1016 
1017  !> @brief Indicates whether the grid discretization supports layers
1018  !<
1019  logical function supports_layers(this)
1020  ! -- dummy
1021  class(distype) :: this
1022  !
1023  supports_layers = .true.
1024  !
1025  end function supports_layers
1026 
1027  !> @brief Return number of cells per layer (nrow * ncol)
1028  !<
1029  function get_ncpl(this)
1030  integer(I4B) :: get_ncpl
1031  class(distype) :: this
1032  get_ncpl = this%nrow * this%ncol
1033  end function get_ncpl
1034 
1035  !> @brief Get normal vector components between the cell and a given neighbor
1036  !!
1037  !! The normal points outward from the shared face between noden and nodem.
1038  !<
1039  subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, &
1040  ipos)
1041  ! -- dummy
1042  class(distype) :: this
1043  integer(I4B), intent(in) :: noden !< cell (reduced nn)
1044  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
1045  integer(I4B), intent(in) :: ihc !< horizontal connection flag
1046  real(DP), intent(inout) :: xcomp
1047  real(DP), intent(inout) :: ycomp
1048  real(DP), intent(inout) :: zcomp
1049  integer(I4B), intent(in) :: ipos
1050  ! -- local
1051  integer(I4B) :: nodeu1, i1, j1, k1
1052  integer(I4B) :: nodeu2, i2, j2, k2
1053  !
1054  ! -- Set vector components based on ihc
1055  if (ihc == 0) then
1056  xcomp = dzero
1057  ycomp = dzero
1058  if (nodem < noden) then
1059  !
1060  ! -- nodem must be above noden, so upward connection
1061  zcomp = done
1062  else
1063  !
1064  ! -- nodem must be below noden, so downward connection
1065  zcomp = -done
1066  end if
1067  else
1068  xcomp = dzero
1069  ycomp = dzero
1070  zcomp = dzero
1071  nodeu1 = this%get_nodeuser(noden)
1072  nodeu2 = this%get_nodeuser(nodem)
1073  call get_ijk(nodeu1, this%nrow, this%ncol, this%nlay, i1, j1, k1)
1074  call get_ijk(nodeu2, this%nrow, this%ncol, this%nlay, i2, j2, k2)
1075  if (i2 < i1) then ! back
1076  ycomp = done
1077  elseif (j2 < j1) then ! left
1078  xcomp = -done
1079  elseif (j2 > j1) then ! right
1080  xcomp = done
1081  else ! front
1082  ycomp = -done
1083  end if
1084  !
1085  end if
1086  !
1087  end subroutine connection_normal
1088 
1089  !> @brief Get unit vector components between the cell and a given neighbor
1090  !!
1091  !! Saturation must be provided to compute cell center vertical coordinates.
1092  !! Also return the straight-line connection length.
1093  !<
1094  subroutine connection_vector(this, noden, nodem, nozee, satn, satm, ihc, &
1095  xcomp, ycomp, zcomp, conlen)
1096  ! -- modules
1097  use disvgeom, only: line_unit_vector
1098  ! -- dummy
1099  class(distype) :: this
1100  integer(I4B), intent(in) :: noden !< cell (reduced nn)
1101  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
1102  logical, intent(in) :: nozee
1103  real(DP), intent(in) :: satn
1104  real(DP), intent(in) :: satm
1105  integer(I4B), intent(in) :: ihc !< horizontal connection flag
1106  real(DP), intent(inout) :: xcomp
1107  real(DP), intent(inout) :: ycomp
1108  real(DP), intent(inout) :: zcomp
1109  real(DP), intent(inout) :: conlen
1110  ! -- local
1111  real(DP) :: z1, z2
1112  real(DP) :: x1, y1, x2, y2
1113  real(DP) :: ds
1114  integer(I4B) :: i1, i2, j1, j2, k1, k2
1115  integer(I4B) :: nodeu1, nodeu2, ipos
1116  !
1117  ! -- Set vector components based on ihc
1118  if (ihc == 0) then
1119  !
1120  ! -- vertical connection; set zcomp positive upward
1121  xcomp = dzero
1122  ycomp = dzero
1123  if (nodem < noden) then
1124  zcomp = done
1125  else
1126  zcomp = -done
1127  end if
1128  z1 = this%bot(noden) + dhalf * (this%top(noden) - this%bot(noden))
1129  z2 = this%bot(nodem) + dhalf * (this%top(nodem) - this%bot(nodem))
1130  conlen = abs(z2 - z1)
1131  else
1132  !
1133  if (nozee) then
1134  z1 = dzero
1135  z2 = dzero
1136  else
1137  z1 = this%bot(noden) + dhalf * satn * (this%top(noden) - this%bot(noden))
1138  z2 = this%bot(nodem) + dhalf * satm * (this%top(nodem) - this%bot(nodem))
1139  end if
1140  ipos = this%con%getjaindex(noden, nodem)
1141  ds = this%con%cl1(this%con%jas(ipos)) + this%con%cl2(this%con%jas(ipos))
1142  nodeu1 = this%get_nodeuser(noden)
1143  nodeu2 = this%get_nodeuser(nodem)
1144  call get_ijk(nodeu1, this%nrow, this%ncol, this%nlay, i1, j1, k1)
1145  call get_ijk(nodeu2, this%nrow, this%ncol, this%nlay, i2, j2, k2)
1146  x1 = dzero
1147  x2 = dzero
1148  y1 = dzero
1149  y2 = dzero
1150  if (i2 < i1) then ! back
1151  y2 = ds
1152  elseif (j2 < j1) then ! left
1153  x2 = -ds
1154  elseif (j2 > j1) then ! right
1155  x2 = ds
1156  else ! front
1157  y2 = -ds
1158  end if
1159  call line_unit_vector(x1, y1, z1, x2, y2, z2, xcomp, ycomp, zcomp, conlen)
1160  end if
1161  !
1162  end subroutine
1163 
1164  !> @brief Get the discretization type
1165  !<
1166  subroutine get_dis_type(this, dis_type)
1167  ! -- dummy
1168  class(distype), intent(in) :: this
1169  character(len=*), intent(out) :: dis_type
1170  !
1171  dis_type = "DIS"
1172  !
1173  end subroutine get_dis_type
1174 
1175  !> @brief Get the discretization type enumeration
1176  function get_dis_enum(this) result(dis_enum)
1177  use constantsmodule, only: dis
1178  class(distype), intent(in) :: this
1179  integer(I4B) :: dis_enum
1180  dis_enum = dis
1181  end function get_dis_enum
1182 
1183  !> @brief Get a 2D array of polygon vertices, listed in
1184  !!
1185  !! clockwise order beginning with the lower left corner
1186  !<
1187  subroutine get_polyverts(this, ic, polyverts, closed)
1188  ! -- dummy
1189  class(distype), intent(inout) :: this
1190  integer(I4B), intent(in) :: ic !< cell number (reduced)
1191  real(DP), allocatable, intent(out) :: polyverts(:, :) !< polygon vertices (column-major indexing)
1192  logical(LGP), intent(in), optional :: closed !< whether to close the polygon, duplicating a vertex
1193  ! -- local
1194  integer(I4B) :: icu, nverts, irow, jcol, klay
1195  real(DP) :: cellx, celly, dxhalf, dyhalf
1196  logical(LGP) :: lclosed
1197  !
1198  nverts = 4
1199  !
1200  ! check closed option
1201  if (.not. (present(closed))) then
1202  lclosed = .false.
1203  else
1204  lclosed = closed
1205  end if
1206  !
1207  ! allocate vertices array
1208  if (lclosed) then
1209  allocate (polyverts(2, nverts + 1))
1210  else
1211  allocate (polyverts(2, nverts))
1212  end if
1213  !
1214  ! set vertices
1215  icu = this%get_nodeuser(ic)
1216  call get_ijk(icu, this%nrow, this%ncol, this%nlay, irow, jcol, klay)
1217  cellx = this%cellx(jcol)
1218  celly = this%celly(irow)
1219  dxhalf = dhalf * this%delr(jcol)
1220  dyhalf = dhalf * this%delc(irow)
1221  polyverts(:, 1) = (/cellx - dxhalf, celly - dyhalf/) ! SW
1222  polyverts(:, 2) = (/cellx - dxhalf, celly + dyhalf/) ! NW
1223  polyverts(:, 3) = (/cellx + dxhalf, celly + dyhalf/) ! NE
1224  polyverts(:, 4) = (/cellx + dxhalf, celly - dyhalf/) ! SE
1225  !
1226  ! close if enabled
1227  if (lclosed) &
1228  polyverts(:, nverts + 1) = polyverts(:, 1)
1229  !
1230  end subroutine
1231 
1232  !> @brief Read an integer array
1233  !<
1234  subroutine read_int_array(this, line, lloc, istart, istop, iout, in, &
1235  iarray, aname)
1236  ! -- dummy
1237  class(distype), intent(inout) :: this
1238  character(len=*), intent(inout) :: line
1239  integer(I4B), intent(inout) :: lloc
1240  integer(I4B), intent(inout) :: istart
1241  integer(I4B), intent(inout) :: istop
1242  integer(I4B), intent(in) :: in
1243  integer(I4B), intent(in) :: iout
1244  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: iarray
1245  character(len=*), intent(in) :: aname
1246  ! -- local
1247  integer(I4B) :: ival
1248  real(DP) :: rval
1249  integer(I4B) :: nlay
1250  integer(I4B) :: nrow
1251  integer(I4B) :: ncol
1252  integer(I4B) :: nval
1253  integer(I4B), dimension(:), pointer, contiguous :: itemp
1254  !
1255  ! -- Point the temporary pointer array, which is passed to the reading
1256  ! subroutine. The temporary array will point to ibuff if it is a
1257  ! reduced structured system, or to iarray if it is an unstructured
1258  ! model.
1259  nlay = this%mshape(1)
1260  nrow = this%mshape(2)
1261  ncol = this%mshape(3)
1262  !
1263  if (this%nodes < this%nodesuser) then
1264  nval = this%nodesuser
1265  itemp => this%ibuff
1266  else
1267  nval = this%nodes
1268  itemp => iarray
1269  end if
1270  !
1271  ! -- Read the array
1272  call urword(line, lloc, istart, istop, 1, ival, rval, iout, in)
1273  if (line(istart:istop) .EQ. 'LAYERED') then
1274  !
1275  ! -- Read layered input
1276  call readarray(in, itemp, aname, this%ndim, ncol, nrow, nlay, nval, &
1277  iout, 1, nlay)
1278  else
1279  !
1280  ! -- Read unstructured input
1281  call readarray(in, itemp, aname, this%ndim, nval, iout, 0)
1282  end if
1283  !
1284  ! -- If reduced model, then need to copy from itemp(=>ibuff) to iarray
1285  if (this%nodes < this%nodesuser) then
1286  call this%fill_grid_array(itemp, iarray)
1287  end if
1288  !
1289  end subroutine read_int_array
1290 
1291  !> @brief Read a double precision array
1292  !<
1293  subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, &
1294  darray, aname)
1295  ! -- dummy
1296  class(distype), intent(inout) :: this
1297  character(len=*), intent(inout) :: line
1298  integer(I4B), intent(inout) :: lloc
1299  integer(I4B), intent(inout) :: istart
1300  integer(I4B), intent(inout) :: istop
1301  integer(I4B), intent(in) :: in
1302  integer(I4B), intent(in) :: iout
1303  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray
1304  character(len=*), intent(in) :: aname
1305  ! -- local
1306  integer(I4B) :: ival
1307  real(DP) :: rval
1308  integer(I4B) :: nlay
1309  integer(I4B) :: nrow
1310  integer(I4B) :: ncol
1311  integer(I4B) :: nval
1312  real(DP), dimension(:), pointer, contiguous :: dtemp
1313  !
1314  ! -- Point the temporary pointer array, which is passed to the reading
1315  ! subroutine. The temporary array will point to dbuff if it is a
1316  ! reduced structured system, or to darray if it is an unstructured
1317  ! model.
1318  nlay = this%mshape(1)
1319  nrow = this%mshape(2)
1320  ncol = this%mshape(3)
1321  !
1322  if (this%nodes < this%nodesuser) then
1323  nval = this%nodesuser
1324  dtemp => this%dbuff
1325  else
1326  nval = this%nodes
1327  dtemp => darray
1328  end if
1329  !
1330  ! -- Read the array
1331  call urword(line, lloc, istart, istop, 1, ival, rval, iout, in)
1332  if (line(istart:istop) .EQ. 'LAYERED') then
1333  !
1334  ! -- Read structured input
1335  call readarray(in, dtemp, aname, this%ndim, ncol, nrow, nlay, nval, &
1336  iout, 1, nlay)
1337  else
1338  !
1339  ! -- Read unstructured input
1340  call readarray(in, dtemp, aname, this%ndim, nval, iout, 0)
1341  end if
1342  !
1343  ! -- If reduced model, then need to copy from dtemp(=>dbuff) to darray
1344  if (this%nodes < this%nodesuser) then
1345  call this%fill_grid_array(dtemp, darray)
1346  end if
1347  !
1348  end subroutine read_dbl_array
1349 
1350  !> @brief Read a 2d double array into col icolbnd of darray
1351  !!
1352  !! For cells that are outside of the active domain,
1353  !! do not copy the array value into darray.
1354  !<
1355  subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, &
1356  icolbnd, aname, inunit, iout)
1357  ! -- dummy
1358  class(distype) :: this
1359  integer(I4B), intent(in) :: maxbnd
1360  integer(I4B), dimension(maxbnd) :: nodelist
1361  integer(I4B), intent(in) :: ncolbnd
1362  real(DP), dimension(ncolbnd, maxbnd), intent(inout) :: darray
1363  integer(I4B), intent(in) :: icolbnd
1364  character(len=*), intent(in) :: aname
1365  integer(I4B), intent(in) :: inunit
1366  integer(I4B), intent(in) :: iout
1367  ! -- local
1368  integer(I4B) :: ir, ic, ncol, nrow, nlay, nval, ipos, nodeu
1369  !
1370  ! -- set variables
1371  nlay = this%mshape(1)
1372  nrow = this%mshape(2)
1373  ncol = this%mshape(3)
1374  !
1375  ! -- Read the array
1376  nval = ncol * nrow
1377  call readarray(inunit, this%dbuff, aname, this%ndim, ncol, nrow, nlay, &
1378  nval, iout, 0, 0)
1379  !
1380  ! -- Copy array into bound. Note that this routine was substantially
1381  ! changed on 9/21/2021 to support changes to READASARRAYS input
1382  ! for recharge and evapotranspiration. nodelist and bound are of
1383  ! size nrow * ncol and correspond directly to dbuff.
1384  ipos = 1
1385  do ir = 1, nrow
1386  do ic = 1, ncol
1387  nodeu = get_node(1, ir, ic, nlay, nrow, ncol)
1388  darray(icolbnd, ipos) = this%dbuff(nodeu)
1389  ipos = ipos + 1
1390  !
1391  end do
1392  end do
1393  !
1394  end subroutine read_layer_array
1395 
1396  !> @brief Record a double precision array.
1397  !!
1398  !! The array is written to a formatted or unformatted external file
1399  !! depending on the arguments.
1400  !<
1401  subroutine record_array(this, darray, iout, iprint, idataun, aname, &
1402  cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
1403  ! -- dummy
1404  class(distype), intent(inout) :: this
1405  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray !< double precision array to record
1406  integer(I4B), intent(in) :: iout !< ascii output unit number
1407  integer(I4B), intent(in) :: iprint !< whether to print the array
1408  integer(I4B), intent(in) :: idataun !< binary output unit number, if negative don't write by layers, write entire array
1409  character(len=*), intent(in) :: aname !< text descriptor
1410  character(len=*), intent(in) :: cdatafmp !< write format
1411  integer(I4B), intent(in) :: nvaluesp !< values per line
1412  integer(I4B), intent(in) :: nwidthp !< number width
1413  character(len=*), intent(in) :: editdesc !< format type (I, G, F, S, E)
1414  real(DP), intent(in) :: dinact !< double precision value for cells excluded from model domain
1415  ! -- local
1416  integer(I4B) :: k, ifirst
1417  integer(I4B) :: nlay
1418  integer(I4B) :: nrow
1419  integer(I4B) :: ncol
1420  integer(I4B) :: nval
1421  integer(I4B) :: nodeu, noder
1422  integer(I4B) :: istart, istop
1423  real(DP), dimension(:), pointer, contiguous :: dtemp
1424  ! -- formats
1425  character(len=*), parameter :: fmthsv = &
1426  "(1X,/1X,a,' WILL BE SAVED ON UNIT ',I4, &
1427  &' AT END OF TIME STEP',I5,', STRESS PERIOD ',I4)"
1428  !
1429  ! -- set variables
1430  nlay = this%mshape(1)
1431  nrow = this%mshape(2)
1432  ncol = this%mshape(3)
1433  !
1434  ! -- If this is a reduced model, then copy the values from darray into
1435  ! dtemp.
1436  if (this%nodes < this%nodesuser) then
1437  nval = this%nodes
1438  dtemp => this%dbuff
1439  do nodeu = 1, this%nodesuser
1440  noder = this%get_nodenumber(nodeu, 0)
1441  if (noder <= 0) then
1442  dtemp(nodeu) = dinact
1443  cycle
1444  end if
1445  dtemp(nodeu) = darray(noder)
1446  end do
1447  else
1448  nval = this%nodes
1449  dtemp => darray
1450  end if
1451  !
1452  ! -- Print to iout if iprint /= 0
1453  if (iprint /= 0) then
1454  istart = 1
1455  do k = 1, nlay
1456  istop = istart + nrow * ncol - 1
1457  call ulaprufw(ncol, nrow, kstp, kper, k, iout, dtemp(istart:istop), &
1458  aname, cdatafmp, nvaluesp, nwidthp, editdesc)
1459  istart = istop + 1
1460  end do
1461  end if
1462  !
1463  ! -- Save array to an external file.
1464  if (idataun > 0) then
1465  ! -- write to binary file by layer
1466  ifirst = 1
1467  istart = 1
1468  do k = 1, nlay
1469  istop = istart + nrow * ncol - 1
1470  if (ifirst == 1) write (iout, fmthsv) &
1471  trim(adjustl(aname)), idataun, &
1472  kstp, kper
1473  ifirst = 0
1474  call ulasav(dtemp(istart:istop), aname, kstp, kper, &
1475  pertim, totim, ncol, nrow, k, idataun)
1476  istart = istop + 1
1477  end do
1478  elseif (idataun < 0) then
1479  !
1480  ! -- write entire array as one record
1481  call ubdsv1(kstp, kper, aname, -idataun, dtemp, ncol, nrow, nlay, &
1482  iout, delt, pertim, totim)
1483  end if
1484  !
1485  end subroutine record_array
1486 
1487  !> @brief Record list header for imeth=6
1488  !<
1489  subroutine record_srcdst_list_header(this, text, textmodel, textpackage, &
1490  dstmodel, dstpackage, naux, auxtxt, &
1491  ibdchn, nlist, iout)
1492  ! -- dummy
1493  class(distype) :: this
1494  character(len=16), intent(in) :: text
1495  character(len=16), intent(in) :: textmodel
1496  character(len=16), intent(in) :: textpackage
1497  character(len=16), intent(in) :: dstmodel
1498  character(len=16), intent(in) :: dstpackage
1499  integer(I4B), intent(in) :: naux
1500  character(len=16), dimension(:), intent(in) :: auxtxt
1501  integer(I4B), intent(in) :: ibdchn
1502  integer(I4B), intent(in) :: nlist
1503  integer(I4B), intent(in) :: iout
1504  ! -- local
1505  integer(I4B) :: nlay, nrow, ncol
1506  !
1507  nlay = this%mshape(1)
1508  nrow = this%mshape(2)
1509  ncol = this%mshape(3)
1510  !
1511  ! -- Use ubdsv06 to write list header
1512  call ubdsv06(kstp, kper, text, textmodel, textpackage, dstmodel, dstpackage, &
1513  ibdchn, naux, auxtxt, ncol, nrow, nlay, &
1514  nlist, iout, delt, pertim, totim)
1515  !
1516  end subroutine record_srcdst_list_header
1517 
1518  !> @brief Convert an integer array (layer numbers) to nodelist
1519  !<
1520  subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
1521  ! -- dummy
1522  class(distype) :: this
1523  integer(I4B), intent(in) :: maxbnd
1524  integer(I4B), dimension(:), pointer, contiguous :: darray
1525  integer(I4B), dimension(maxbnd), intent(inout) :: nodelist
1526  integer(I4B), intent(inout) :: nbound
1527  character(len=*), intent(in) :: aname
1528  ! -- local
1529  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nval, nodeu, noder, ipos, ierr
1530  !
1531  ! -- set variables
1532  nlay = this%mshape(1)
1533  nrow = this%mshape(2)
1534  ncol = this%mshape(3)
1535  !
1536  if (this%ndim > 1) then
1537  !
1538  nval = ncol * nrow
1539  !
1540  ! -- Copy array into nodelist
1541  ipos = 1
1542  ierr = 0
1543  do ir = 1, nrow
1544  do ic = 1, ncol
1545  nodeu = get_node(1, ir, ic, nlay, nrow, ncol)
1546  il = darray(nodeu)
1547  if (il < 1 .or. il > nlay) then
1548  write (errmsg, '(a,1x,i0)') 'Invalid layer number:', il
1549  call store_error(errmsg, terminate=.true.)
1550  end if
1551  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
1552  noder = this%get_nodenumber(nodeu, 0)
1553  if (ipos > maxbnd) then
1554  ierr = ipos
1555  else
1556  nodelist(ipos) = noder
1557  end if
1558  ipos = ipos + 1
1559  end do
1560  end do
1561  !
1562  ! -- Check for errors
1563  nbound = ipos - 1
1564  if (ierr > 0) then
1565  write (errmsg, '(a,1x,i0)') &
1566  'MAXBOUND dimension is too small.'// &
1567  'INCREASE MAXBOUND TO:', ierr
1568  call store_error(errmsg, terminate=.true.)
1569  end if
1570  !
1571  ! -- If nbound < maxbnd, then initialize nodelist to zero in this range
1572  if (nbound < maxbnd) then
1573  do ipos = nbound + 1, maxbnd
1574  nodelist(ipos) = 0
1575  end do
1576  end if
1577  !
1578  else
1579  !
1580  ! -- For unstructured, read nodelist directly, then check node numbers
1581  nodelist = darray
1582  do noder = 1, maxbnd
1583  if (noder < 1 .or. noder > this%nodes) then
1584  write (errmsg, '(a,1x,i0)') 'Invalid node number:', noder
1585  call store_error(errmsg, terminate=.true.)
1586  end if
1587  end do
1588  nbound = maxbnd
1589  !
1590  end if
1591  !
1592  end subroutine nlarray_to_nodelist
1593 
1594 end module dismodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
@ dis
DIS6 discretization.
Definition: Constants.f90:154
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:67
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
Definition: Dis.f90:1
subroutine dis3d_df(this)
Define the discretization.
Definition: Dis.f90:130
integer(i4b) function get_ncpl(this)
Return number of cells per layer (nrow * ncol)
Definition: Dis.f90:1030
subroutine nodeu_to_array(this, nodeu, arr)
Convert a user nodenumber to an array (nodenumber) or (k,i,j)
Definition: Dis.f90:685
subroutine dis3d_da(this)
Deallocate variables.
Definition: Dis.f90:154
integer(i4b) function get_dis_enum(this)
Get the discretization type enumeration.
Definition: Dis.f90:1177
subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
Convert an integer array (layer numbers) to nodelist.
Definition: Dis.f90:1521
subroutine record_array(this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
Record a double precision array.
Definition: Dis.f90:1403
subroutine log_options(this, found)
Write user options to list file.
Definition: Dis.f90:209
integer(i4b) function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, allow_zero)
Convert a cellid string to a user nodenumber.
Definition: Dis.f90:936
subroutine allocate_scalars(this, name_model, input_mempath)
Allocate and initialize scalar variables.
Definition: Dis.f90:790
subroutine write_grb(this, icelltype)
Write a binary grid file.
Definition: Dis.f90:548
subroutine log_griddata(this, found)
Write dimensions to list file.
Definition: Dis.f90:351
subroutine, public dis_cr(dis, name_model, input_mempath, inunit, iout)
Create a new structured discretization object.
Definition: Dis.f90:97
integer(i4b) function get_nodenumber_idx3(this, k, i, j, icheck)
Get reduced node number from layer, row and column indices.
Definition: Dis.f90:746
logical function supports_layers(this)
Indicates whether the grid discretization supports layers.
Definition: Dis.f90:1020
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.
Definition: Dis.f90:1096
integer(i4b) function nodeu_from_string(this, lloc, istart, istop, in, iout, line, flag_string, allow_zero)
Convert a string to a user nodenumber.
Definition: Dis.f90:845
subroutine get_dis_type(this, dis_type)
Get the discretization type.
Definition: Dis.f90:1167
subroutine nodeu_to_string(this, nodeu, str)
Convert a user nodenumber to a string (nodenumber) or (k,i,j)
Definition: Dis.f90:664
subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, icolbnd, aname, inunit, iout)
Read a 2d double array into col icolbnd of darray.
Definition: Dis.f90:1357
subroutine allocate_arrays(this)
Allocate and initialize arrays.
Definition: Dis.f90:814
subroutine source_dimensions(this)
Copy dimensions from IDM into package.
Definition: Dis.f90:244
subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, ipos)
Get normal vector components between the cell and a given neighbor.
Definition: Dis.f90:1041
integer(i4b) function get_nodenumber_idx1(this, nodeu, icheck)
Get reduced node number from user node number.
Definition: Dis.f90:715
subroutine grid_finalize(this)
Finalize grid (check properties, allocate arrays, compute connections)
Definition: Dis.f90:384
subroutine source_options(this)
Copy options from IDM into package.
Definition: Dis.f90:184
subroutine read_int_array(this, line, lloc, istart, istop, iout, in, iarray, aname)
Read an integer array.
Definition: Dis.f90:1236
subroutine source_griddata(this)
Copy grid data from IDM into package.
Definition: Dis.f90:330
subroutine get_polyverts(this, ic, polyverts, closed)
Get a 2D array of polygon vertices, listed in.
Definition: Dis.f90:1188
subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, darray, aname)
Read a double precision array.
Definition: Dis.f90:1295
subroutine record_srcdst_list_header(this, text, textmodel, textpackage, dstmodel, dstpackage, naux, auxtxt, ibdchn, nlist, iout)
Record list header for imeth=6.
Definition: Dis.f90:1492
subroutine log_dimensions(this, found)
Write dimensions to list file.
Definition: Dis.f90:305
subroutine, public line_unit_vector(x0, y0, z0, x1, y1, z1, xcomp, ycomp, zcomp, vmag)
Calculate the vector components (xcomp, ycomp, and zcomp) for a line defined by two points,...
Definition: DisvGeom.f90:508
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:80
subroutine, public get_ijk(nodenumber, nrow, ncol, nlay, irow, icol, ilay)
Get row, column and layer indices from node number and grid dimensions. If nodenumber is invalid,...
Definition: GeomUtil.f90:97
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...
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public ulaprufw(ncol, nrow, kstp, kper, ilay, iout, buf, text, userfmt, nvalues, nwidth, editdesc)
Print 1 layer array with user formatting in wrap format.
subroutine, public ubdsv06(kstp, kper, text, modelnam1, paknam1, modelnam2, paknam2, ibdchn, naux, auxtxt, ncol, nrow, nlay, nlist, iout, delt, pertim, totim)
Write header records for cell-by-cell flow terms for one component of flow.
subroutine, public ulasav(buf, text, kstp, kper, pertim, totim, ncol, nrow, ilay, ichn)
Save 1 layer array on disk.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
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
subroutine, public memorylist_remove(component, subcomponent, context)
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
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
character(len=linelength) idm_context
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
Simplifies tracking parameters sourced from the input context.
Definition: Dis.f90:76
Structured grid discretization.
Definition: Dis.f90:23