MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
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 memorystore_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  use constantsmodule, only: lenbigline
551  ! -- dummy
552  class(distype) :: this
553  integer(I4B), dimension(:), intent(in) :: icelltype
554  ! -- local
555  integer(I4B) :: iunit, ntxt, ncpl, version
556  integer(I4B), parameter :: lentxt = 100
557  character(len=50) :: txthdr
558  character(len=lentxt) :: txt
559  character(len=LINELENGTH) :: fname
560  character(len=LENBIGLINE) :: crs
561  logical(LGP) :: found_crs
562  character(len=*), parameter :: fmtgrdsave = &
563  "(4X,'BINARY GRID INFORMATION WILL BE WRITTEN TO:', &
564  &/,6X,'UNIT NUMBER: ', I0,/,6X, 'FILE NAME: ', A)"
565  !
566  ! -- Initialize
567  version = 1
568  ntxt = 16
569  ncpl = this%nrow * this%ncol
570  !
571  call mem_set_value(crs, 'CRS', this%input_mempath, found_crs)
572  !
573  ! -- set version
574  if (found_crs) then
575  ntxt = ntxt + 1
576  version = 2
577  end if
578  !
579  ! -- Open the file
580  fname = trim(this%output_fname)
581  iunit = getunit()
582  write (this%iout, fmtgrdsave) iunit, trim(adjustl(fname))
583  call openfile(iunit, this%iout, trim(adjustl(fname)), 'DATA(BINARY)', &
584  form, access, 'REPLACE')
585  !
586  ! -- write header information
587  write (txthdr, '(a)') 'GRID DIS'
588  txthdr(50:50) = new_line('a')
589  write (iunit) txthdr
590  write (txthdr, '(a, i0)') 'VERSION ', version
591  txthdr(50:50) = new_line('a')
592  write (iunit) txthdr
593  write (txthdr, '(a, i0)') 'NTXT ', ntxt
594  txthdr(50:50) = new_line('a')
595  write (iunit) txthdr
596  write (txthdr, '(a, i0)') 'LENTXT ', lentxt
597  txthdr(50:50) = new_line('a')
598  write (iunit) txthdr
599  !
600  ! -- write variable definitions
601  write (txt, '(3a, i0)') 'NCELLS ', 'INTEGER ', 'NDIM 0 # ', this%nodesuser
602  txt(lentxt:lentxt) = new_line('a')
603  write (iunit) txt
604  write (txt, '(3a, i0)') 'NLAY ', 'INTEGER ', 'NDIM 0 # ', this%nlay
605  txt(lentxt:lentxt) = new_line('a')
606  write (iunit) txt
607  write (txt, '(3a, i0)') 'NROW ', 'INTEGER ', 'NDIM 0 # ', this%nrow
608  txt(lentxt:lentxt) = new_line('a')
609  write (iunit) txt
610  write (txt, '(3a, i0)') 'NCOL ', 'INTEGER ', 'NDIM 0 # ', this%ncol
611  txt(lentxt:lentxt) = new_line('a')
612  write (iunit) txt
613  write (txt, '(3a, i0)') 'NJA ', 'INTEGER ', 'NDIM 0 # ', this%nja
614  txt(lentxt:lentxt) = new_line('a')
615  write (iunit) txt
616  write (txt, '(3a, 1pg24.15)') 'XORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%xorigin
617  txt(lentxt:lentxt) = new_line('a')
618  write (iunit) txt
619  write (txt, '(3a, 1pg24.15)') 'YORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%yorigin
620  txt(lentxt:lentxt) = new_line('a')
621  write (iunit) txt
622  write (txt, '(3a, 1pg24.15)') 'ANGROT ', 'DOUBLE ', 'NDIM 0 # ', this%angrot
623  txt(lentxt:lentxt) = new_line('a')
624  write (iunit) txt
625  write (txt, '(3a, i0)') 'DELR ', 'DOUBLE ', 'NDIM 1 ', this%ncol
626  txt(lentxt:lentxt) = new_line('a')
627  write (iunit) txt
628  write (txt, '(3a, i0)') 'DELC ', 'DOUBLE ', 'NDIM 1 ', this%nrow
629  txt(lentxt:lentxt) = new_line('a')
630  write (iunit) txt
631  write (txt, '(3a, i0)') 'TOP ', 'DOUBLE ', 'NDIM 1 ', ncpl
632  txt(lentxt:lentxt) = new_line('a')
633  write (iunit) txt
634  write (txt, '(3a, i0)') 'BOTM ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
635  txt(lentxt:lentxt) = new_line('a')
636  write (iunit) txt
637  write (txt, '(3a, i0)') 'IA ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
638  txt(lentxt:lentxt) = new_line('a')
639  write (iunit) txt
640  write (txt, '(3a, i0)') 'JA ', 'INTEGER ', 'NDIM 1 ', size(this%con%jausr)
641  txt(lentxt:lentxt) = new_line('a')
642  write (iunit) txt
643  write (txt, '(3a, i0)') 'IDOMAIN ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
644  txt(lentxt:lentxt) = new_line('a')
645  write (iunit) txt
646  write (txt, '(3a, i0)') 'ICELLTYPE ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
647  txt(lentxt:lentxt) = new_line('a')
648  write (iunit) txt
649  !
650  ! -- if version 2 write character array headers
651  if (version == 2) then
652  if (found_crs) then
653  write (txt, '(3a, i0)') 'CRS ', 'CHARACTER ', 'NDIM 1 ', &
654  len_trim(crs)
655  txt(lentxt:lentxt) = new_line('a')
656  write (iunit) txt
657  end if
658  end if
659  !
660  ! -- write data
661  write (iunit) this%nodesuser ! ncells
662  write (iunit) this%nlay ! nlay
663  write (iunit) this%nrow ! nrow
664  write (iunit) this%ncol ! ncol
665  write (iunit) this%nja ! nja
666  write (iunit) this%xorigin ! xorigin
667  write (iunit) this%yorigin ! yorigin
668  write (iunit) this%angrot ! angrot
669  write (iunit) this%delr ! delr
670  write (iunit) this%delc ! delc
671  write (iunit) this%top2d ! top2d
672  write (iunit) this%bot3d ! bot3d
673  write (iunit) this%con%iausr ! iausr
674  write (iunit) this%con%jausr ! jausr
675  write (iunit) this%idomain ! idomain
676  write (iunit) icelltype ! icelltype
677  !
678  ! -- if version 2 write character array data
679  if (version == 2) then
680  if (found_crs) write (iunit) trim(crs) ! crs user input
681  end if
682  !
683  ! -- Close the file
684  close (iunit)
685  !
686  end subroutine write_grb
687 
688  !> @brief Convert a user nodenumber to a string (nodenumber) or (k,i,j)
689  !<
690  subroutine nodeu_to_string(this, nodeu, str)
691  ! -- dummy
692  class(distype) :: this
693  integer(I4B), intent(in) :: nodeu
694  character(len=*), intent(inout) :: str
695  ! -- local
696  integer(I4B) :: i, j, k
697  character(len=10) :: kstr, istr, jstr
698  !
699  call get_ijk(nodeu, this%nrow, this%ncol, this%nlay, i, j, k)
700  write (kstr, '(i10)') k
701  write (istr, '(i10)') i
702  write (jstr, '(i10)') j
703  str = '('//trim(adjustl(kstr))//','// &
704  trim(adjustl(istr))//','// &
705  trim(adjustl(jstr))//')'
706  !
707  end subroutine nodeu_to_string
708 
709  !> @brief Convert a user nodenumber to an array (nodenumber) or (k,i,j)
710  !<
711  subroutine nodeu_to_array(this, nodeu, arr)
712  ! -- dummy
713  class(distype) :: this
714  integer(I4B), intent(in) :: nodeu
715  integer(I4B), dimension(:), intent(inout) :: arr
716  ! -- local
717  integer(I4B) :: isize
718  integer(I4B) :: i, j, k
719  !
720  ! -- check the size of arr
721  isize = size(arr)
722  if (isize /= this%ndim) then
723  write (errmsg, '(a,i0,a,i0,a)') &
724  'Program error: nodeu_to_array size of array (', isize, &
725  ') is not equal to the discretization dimension (', this%ndim, ')'
726  call store_error(errmsg, terminate=.true.)
727  end if
728  !
729  ! -- get k, i, j
730  call get_ijk(nodeu, this%nrow, this%ncol, this%nlay, i, j, k)
731  !
732  ! -- fill array
733  arr(1) = k
734  arr(2) = i
735  arr(3) = j
736  !
737  end subroutine nodeu_to_array
738 
739  !> @brief Get reduced node number from user node number
740  !<
741  function get_nodenumber_idx1(this, nodeu, icheck) result(nodenumber)
742  ! -- return
743  integer(I4B) :: nodenumber
744  ! -- dummy
745  class(distype), intent(in) :: this
746  integer(I4B), intent(in) :: nodeu
747  integer(I4B), intent(in) :: icheck
748  !
749  ! -- check the node number if requested
750  if (icheck /= 0) then
751  !
752  ! -- If within valid range, convert to reduced nodenumber
753  if (nodeu < 1 .or. nodeu > this%nodesuser) then
754  write (errmsg, '(a,i0,a)') &
755  'Node number (', nodeu, &
756  ') less than 1 or greater than the number of nodes.'
757  call store_error(errmsg)
758  nodenumber = 0
759  else
760  nodenumber = nodeu
761  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
762  end if
763  else
764  nodenumber = nodeu
765  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
766  end if
767  !
768  end function get_nodenumber_idx1
769 
770  !> @brief Get reduced node number from layer, row and column indices
771  !<
772  function get_nodenumber_idx3(this, k, i, j, icheck) result(nodenumber)
773  ! -- return
774  integer(I4B) :: nodenumber
775  ! -- dummy
776  class(distype), intent(in) :: this
777  integer(I4B), intent(in) :: k, i, j
778  integer(I4B), intent(in) :: icheck
779  ! -- local
780  integer(I4B) :: nodeu
781  ! formats
782  character(len=*), parameter :: fmterr = &
783  "('Error in structured-grid cell indices: layer = ',i0,', &
784  &row = ',i0,', column = ',i0)"
785  !
786  nodeu = get_node(k, i, j, this%nlay, this%nrow, this%ncol)
787  if (nodeu < 1) then
788  write (errmsg, fmterr) k, i, j
789  call store_error(errmsg, terminate=.true.)
790  end if
791  nodenumber = nodeu
792  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
793  !
794  ! -- check the node number if requested
795  if (icheck /= 0) then
796  !
797  if (k < 1 .or. k > this%nlay) &
798  call store_error('Layer less than one or greater than nlay')
799  if (i < 1 .or. i > this%nrow) &
800  call store_error('Row less than one or greater than nrow')
801  if (j < 1 .or. j > this%ncol) &
802  call store_error('Column less than one or greater than ncol')
803  !
804  ! -- Error if outside of range
805  if (nodeu < 1 .or. nodeu > this%nodesuser) then
806  write (errmsg, '(a,i0,a)') &
807  'Node number (', nodeu, ')less than 1 or greater than nodes.'
808  call store_error(errmsg)
809  end if
810  end if
811  !
812  end function get_nodenumber_idx3
813 
814  !> @brief Allocate and initialize scalar variables
815  !<
816  subroutine allocate_scalars(this, name_model, input_mempath)
817  ! -- dummy
818  class(distype) :: this
819  character(len=*), intent(in) :: name_model
820  character(len=*), intent(in) :: input_mempath
821  !
822  ! -- Allocate parent scalars
823  call this%DisBaseType%allocate_scalars(name_model, input_mempath)
824  !
825  ! -- Allocate
826  call mem_allocate(this%nlay, 'NLAY', this%memoryPath)
827  call mem_allocate(this%nrow, 'NROW', this%memoryPath)
828  call mem_allocate(this%ncol, 'NCOL', this%memoryPath)
829  !
830  ! -- Initialize
831  this%nlay = 0
832  this%nrow = 0
833  this%ncol = 0
834  this%ndim = 3
835  !
836  end subroutine allocate_scalars
837 
838  !> @brief Allocate and initialize arrays
839  !<
840  subroutine allocate_arrays(this)
841  ! -- dummy
842  class(distype) :: this
843  !
844  ! -- Allocate arrays in DisBaseType (mshape, top, bot, area)
845  call this%DisBaseType%allocate_arrays()
846  !
847  ! -- Allocate arrays for DisType
848  if (this%nodes < this%nodesuser) then
849  call mem_allocate(this%nodeuser, this%nodes, 'NODEUSER', this%memoryPath)
850  call mem_allocate(this%nodereduced, this%nodesuser, 'NODEREDUCED', &
851  this%memoryPath)
852  else
853  call mem_allocate(this%nodeuser, 1, 'NODEUSER', this%memoryPath)
854  call mem_allocate(this%nodereduced, 1, 'NODEREDUCED', this%memoryPath)
855  end if
856  !
857  ! -- Initialize
858  this%mshape(1) = this%nlay
859  this%mshape(2) = this%nrow
860  this%mshape(3) = this%ncol
861  !
862  end subroutine allocate_arrays
863 
864  !> @brief Convert a string to a user nodenumber
865  !!
866  !! Parse layer, row and column and return user nodenumber.
867  !! If flag_string is present and true, the first token may be
868  !! non-numeric (e.g. boundary name). In this case, return -2.
869  !<
870  function nodeu_from_string(this, lloc, istart, istop, in, iout, line, &
871  flag_string, allow_zero) result(nodeu)
872  ! -- dummy
873  class(distype) :: this
874  integer(I4B), intent(inout) :: lloc
875  integer(I4B), intent(inout) :: istart
876  integer(I4B), intent(inout) :: istop
877  integer(I4B), intent(in) :: in
878  integer(I4B), intent(in) :: iout
879  character(len=*), intent(inout) :: line
880  logical, optional, intent(in) :: flag_string
881  logical, optional, intent(in) :: allow_zero
882  integer(I4B) :: nodeu
883  ! -- local
884  integer(I4B) :: k, i, j, nlay, nrow, ncol
885  integer(I4B) :: lloclocal, ndum, istat, n
886  real(dp) :: r
887  !
888  if (present(flag_string)) then
889  if (flag_string) then
890  ! Check to see if first token in line can be read as an integer.
891  lloclocal = lloc
892  call urword(line, lloclocal, istart, istop, 1, ndum, r, iout, in)
893  read (line(istart:istop), *, iostat=istat) n
894  if (istat /= 0) then
895  ! First token in line is not an integer; return flag to this effect.
896  nodeu = -2
897  return
898  end if
899  end if
900  end if
901  !
902  nlay = this%mshape(1)
903  nrow = this%mshape(2)
904  ncol = this%mshape(3)
905  !
906  call urword(line, lloc, istart, istop, 2, k, r, iout, in)
907  call urword(line, lloc, istart, istop, 2, i, r, iout, in)
908  call urword(line, lloc, istart, istop, 2, j, r, iout, in)
909  !
910  if (k == 0 .and. i == 0 .and. j == 0) then
911  if (present(allow_zero)) then
912  if (allow_zero) then
913  nodeu = 0
914  return
915  end if
916  end if
917  end if
918  !
919  errmsg = ""
920  !
921  if (k < 1 .or. k > nlay) then
922  write (errmsg, '(a,i0,a)') &
923  'Layer number in list (', k, ') is outside of the grid.'
924  end if
925  if (i < 1 .or. i > nrow) then
926  write (errmsg, '(a,1x,a,i0,a)') &
927  trim(adjustl(errmsg)), 'Row number in list (', i, &
928  ') is outside of the grid.'
929  end if
930  if (j < 1 .or. j > ncol) then
931  write (errmsg, '(a,1x,a,i0,a)') &
932  trim(adjustl(errmsg)), 'Column number in list (', j, &
933  ') is outside of the grid.'
934  end if
935  !
936  nodeu = get_node(k, i, j, nlay, nrow, ncol)
937  !
938  if (nodeu < 1 .or. nodeu > this%nodesuser) then
939  write (errmsg, '(a,1x,a,i0,a)') &
940  trim(adjustl(errmsg)), &
941  "Node number in list (", nodeu, ") is outside of the grid. "// &
942  "Cell number cannot be determined in line '"// &
943  trim(adjustl(line))//"'."
944  end if
945  !
946  if (len_trim(adjustl(errmsg)) > 0) then
947  call store_error(errmsg)
948  call store_error_unit(in)
949  end if
950  !
951  end function nodeu_from_string
952 
953  !> @brief Convert a cellid string to a user nodenumber
954  !!
955  !! If flag_string is present and true, the first token may be
956  !! non-numeric (e.g. boundary name). In this case, return -2.
957  !!
958  !! If allow_zero is present and true, and all indices are zero, the
959  !! result can be zero. If allow_zero is false, a zero in any index is an error.
960  !<
961  function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, &
962  allow_zero) result(nodeu)
963  ! -- return
964  integer(I4B) :: nodeu
965  ! -- dummy
966  class(distype) :: this
967  character(len=*), intent(inout) :: cellid
968  integer(I4B), intent(in) :: inunit
969  integer(I4B), intent(in) :: iout
970  logical, optional, intent(in) :: flag_string
971  logical, optional, intent(in) :: allow_zero
972  ! -- local
973  integer(I4B) :: lloclocal, istart, istop, ndum, n
974  integer(I4B) :: k, i, j, nlay, nrow, ncol
975  integer(I4B) :: istat
976  real(dp) :: r
977  !
978  if (present(flag_string)) then
979  if (flag_string) then
980  ! Check to see if first token in cellid can be read as an integer.
981  lloclocal = 1
982  call urword(cellid, lloclocal, istart, istop, 1, ndum, r, iout, inunit)
983  read (cellid(istart:istop), *, iostat=istat) n
984  if (istat /= 0) then
985  ! First token in cellid is not an integer; return flag to this effect.
986  nodeu = -2
987  return
988  end if
989  end if
990  end if
991  !
992  nlay = this%mshape(1)
993  nrow = this%mshape(2)
994  ncol = this%mshape(3)
995  !
996  lloclocal = 1
997  call urword(cellid, lloclocal, istart, istop, 2, k, r, iout, inunit)
998  call urword(cellid, lloclocal, istart, istop, 2, i, r, iout, inunit)
999  call urword(cellid, lloclocal, istart, istop, 2, j, r, iout, inunit)
1000  !
1001  if (k == 0 .and. i == 0 .and. j == 0) then
1002  if (present(allow_zero)) then
1003  if (allow_zero) then
1004  nodeu = 0
1005  return
1006  end if
1007  end if
1008  end if
1009  !
1010  errmsg = ""
1011  !
1012  if (k < 1 .or. k > nlay) then
1013  write (errmsg, '(a,i0,a)') &
1014  'Layer number in list (', k, ') is outside of the grid.'
1015  end if
1016  if (i < 1 .or. i > nrow) then
1017  write (errmsg, '(a,1x,a,i0,a)') &
1018  trim(adjustl(errmsg)), 'Row number in list (', i, &
1019  ') is outside of the grid.'
1020  end if
1021  if (j < 1 .or. j > ncol) then
1022  write (errmsg, '(a,1x,a,i0,a)') &
1023  trim(adjustl(errmsg)), 'Column number in list (', j, &
1024  ') is outside of the grid.'
1025  end if
1026  !
1027  nodeu = get_node(k, i, j, nlay, nrow, ncol)
1028  !
1029  if (nodeu < 1 .or. nodeu > this%nodesuser) then
1030  write (errmsg, '(a,1x,a,i0,a)') &
1031  trim(adjustl(errmsg)), &
1032  "Cell number cannot be determined for cellid ("// &
1033  trim(adjustl(cellid))//") and results in a user "// &
1034  "node number (", nodeu, ") that is outside of the grid."
1035  end if
1036  !
1037  if (len_trim(adjustl(errmsg)) > 0) then
1038  call store_error(errmsg)
1039  call store_error_unit(inunit)
1040  end if
1041  !
1042  end function nodeu_from_cellid
1043 
1044  !> @brief Indicates whether the grid discretization supports layers
1045  !<
1046  logical function supports_layers(this)
1047  ! -- dummy
1048  class(distype) :: this
1049  !
1050  supports_layers = .true.
1051  !
1052  end function supports_layers
1053 
1054  !> @brief Return number of cells per layer (nrow * ncol)
1055  !<
1056  function get_ncpl(this)
1057  integer(I4B) :: get_ncpl
1058  class(distype) :: this
1059  get_ncpl = this%nrow * this%ncol
1060  end function get_ncpl
1061 
1062  !> @brief Get normal vector components between the cell and a given neighbor
1063  !!
1064  !! The normal points outward from the shared face between noden and nodem.
1065  !<
1066  subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, &
1067  ipos)
1068  ! -- dummy
1069  class(distype) :: this
1070  integer(I4B), intent(in) :: noden !< cell (reduced nn)
1071  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
1072  integer(I4B), intent(in) :: ihc !< horizontal connection flag
1073  real(DP), intent(inout) :: xcomp
1074  real(DP), intent(inout) :: ycomp
1075  real(DP), intent(inout) :: zcomp
1076  integer(I4B), intent(in) :: ipos
1077  ! -- local
1078  integer(I4B) :: nodeu1, i1, j1, k1
1079  integer(I4B) :: nodeu2, i2, j2, k2
1080  !
1081  ! -- Set vector components based on ihc
1082  if (ihc == 0) then
1083  xcomp = dzero
1084  ycomp = dzero
1085  if (nodem < noden) then
1086  !
1087  ! -- nodem must be above noden, so upward connection
1088  zcomp = done
1089  else
1090  !
1091  ! -- nodem must be below noden, so downward connection
1092  zcomp = -done
1093  end if
1094  else
1095  xcomp = dzero
1096  ycomp = dzero
1097  zcomp = dzero
1098  nodeu1 = this%get_nodeuser(noden)
1099  nodeu2 = this%get_nodeuser(nodem)
1100  call get_ijk(nodeu1, this%nrow, this%ncol, this%nlay, i1, j1, k1)
1101  call get_ijk(nodeu2, this%nrow, this%ncol, this%nlay, i2, j2, k2)
1102  if (i2 < i1) then ! back
1103  ycomp = done
1104  elseif (j2 < j1) then ! left
1105  xcomp = -done
1106  elseif (j2 > j1) then ! right
1107  xcomp = done
1108  else ! front
1109  ycomp = -done
1110  end if
1111  !
1112  end if
1113  !
1114  end subroutine connection_normal
1115 
1116  !> @brief Get unit vector components between the cell and a given neighbor
1117  !!
1118  !! Saturation must be provided to compute cell center vertical coordinates.
1119  !! Also return the straight-line connection length.
1120  !<
1121  subroutine connection_vector(this, noden, nodem, nozee, satn, satm, ihc, &
1122  xcomp, ycomp, zcomp, conlen)
1123  ! -- modules
1124  use disvgeom, only: line_unit_vector
1125  ! -- dummy
1126  class(distype) :: this
1127  integer(I4B), intent(in) :: noden !< cell (reduced nn)
1128  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
1129  logical, intent(in) :: nozee
1130  real(DP), intent(in) :: satn
1131  real(DP), intent(in) :: satm
1132  integer(I4B), intent(in) :: ihc !< horizontal connection flag
1133  real(DP), intent(inout) :: xcomp
1134  real(DP), intent(inout) :: ycomp
1135  real(DP), intent(inout) :: zcomp
1136  real(DP), intent(inout) :: conlen
1137  ! -- local
1138  real(DP) :: z1, z2
1139  real(DP) :: x1, y1, x2, y2
1140  real(DP) :: ds
1141  integer(I4B) :: i1, i2, j1, j2, k1, k2
1142  integer(I4B) :: nodeu1, nodeu2, ipos
1143  !
1144  ! -- Set vector components based on ihc
1145  if (ihc == 0) then
1146  !
1147  ! -- vertical connection; set zcomp positive upward
1148  xcomp = dzero
1149  ycomp = dzero
1150  if (nodem < noden) then
1151  zcomp = done
1152  else
1153  zcomp = -done
1154  end if
1155  z1 = this%bot(noden) + dhalf * (this%top(noden) - this%bot(noden))
1156  z2 = this%bot(nodem) + dhalf * (this%top(nodem) - this%bot(nodem))
1157  conlen = abs(z2 - z1)
1158  else
1159  !
1160  if (nozee) then
1161  z1 = dzero
1162  z2 = dzero
1163  else
1164  z1 = this%bot(noden) + dhalf * satn * (this%top(noden) - this%bot(noden))
1165  z2 = this%bot(nodem) + dhalf * satm * (this%top(nodem) - this%bot(nodem))
1166  end if
1167  ipos = this%con%getjaindex(noden, nodem)
1168  ds = this%con%cl1(this%con%jas(ipos)) + this%con%cl2(this%con%jas(ipos))
1169  nodeu1 = this%get_nodeuser(noden)
1170  nodeu2 = this%get_nodeuser(nodem)
1171  call get_ijk(nodeu1, this%nrow, this%ncol, this%nlay, i1, j1, k1)
1172  call get_ijk(nodeu2, this%nrow, this%ncol, this%nlay, i2, j2, k2)
1173  x1 = dzero
1174  x2 = dzero
1175  y1 = dzero
1176  y2 = dzero
1177  if (i2 < i1) then ! back
1178  y2 = ds
1179  elseif (j2 < j1) then ! left
1180  x2 = -ds
1181  elseif (j2 > j1) then ! right
1182  x2 = ds
1183  else ! front
1184  y2 = -ds
1185  end if
1186  call line_unit_vector(x1, y1, z1, x2, y2, z2, xcomp, ycomp, zcomp, conlen)
1187  end if
1188  !
1189  end subroutine
1190 
1191  !> @brief Get the discretization type
1192  !<
1193  subroutine get_dis_type(this, dis_type)
1194  ! -- dummy
1195  class(distype), intent(in) :: this
1196  character(len=*), intent(out) :: dis_type
1197  !
1198  dis_type = "DIS"
1199  !
1200  end subroutine get_dis_type
1201 
1202  !> @brief Get the discretization type enumeration
1203  function get_dis_enum(this) result(dis_enum)
1204  use constantsmodule, only: dis
1205  class(distype), intent(in) :: this
1206  integer(I4B) :: dis_enum
1207  dis_enum = dis
1208  end function get_dis_enum
1209 
1210  !> @brief Get a 2D array of polygon vertices, listed in
1211  !!
1212  !! clockwise order beginning with the lower left corner
1213  !<
1214  subroutine get_polyverts(this, ic, polyverts, closed)
1215  ! -- dummy
1216  class(distype), intent(inout) :: this
1217  integer(I4B), intent(in) :: ic !< cell number (reduced)
1218  real(DP), allocatable, intent(out) :: polyverts(:, :) !< polygon vertices (column-major indexing)
1219  logical(LGP), intent(in), optional :: closed !< whether to close the polygon, duplicating a vertex
1220  ! -- local
1221  integer(I4B) :: icu, nverts, irow, jcol, klay
1222  real(DP) :: cellx, celly, dxhalf, dyhalf
1223  logical(LGP) :: lclosed
1224  !
1225  nverts = 4
1226  !
1227  ! check closed option
1228  if (.not. (present(closed))) then
1229  lclosed = .false.
1230  else
1231  lclosed = closed
1232  end if
1233  !
1234  ! allocate vertices array
1235  if (lclosed) then
1236  allocate (polyverts(2, nverts + 1))
1237  else
1238  allocate (polyverts(2, nverts))
1239  end if
1240  !
1241  ! set vertices
1242  icu = this%get_nodeuser(ic)
1243  call get_ijk(icu, this%nrow, this%ncol, this%nlay, irow, jcol, klay)
1244  cellx = this%cellx(jcol)
1245  celly = this%celly(irow)
1246  dxhalf = dhalf * this%delr(jcol)
1247  dyhalf = dhalf * this%delc(irow)
1248  polyverts(:, 1) = (/cellx - dxhalf, celly - dyhalf/) ! SW
1249  polyverts(:, 2) = (/cellx - dxhalf, celly + dyhalf/) ! NW
1250  polyverts(:, 3) = (/cellx + dxhalf, celly + dyhalf/) ! NE
1251  polyverts(:, 4) = (/cellx + dxhalf, celly - dyhalf/) ! SE
1252  !
1253  ! close if enabled
1254  if (lclosed) &
1255  polyverts(:, nverts + 1) = polyverts(:, 1)
1256  !
1257  end subroutine
1258 
1259  !> @brief Read an integer array
1260  !<
1261  subroutine read_int_array(this, line, lloc, istart, istop, iout, in, &
1262  iarray, aname)
1263  ! -- dummy
1264  class(distype), intent(inout) :: this
1265  character(len=*), intent(inout) :: line
1266  integer(I4B), intent(inout) :: lloc
1267  integer(I4B), intent(inout) :: istart
1268  integer(I4B), intent(inout) :: istop
1269  integer(I4B), intent(in) :: in
1270  integer(I4B), intent(in) :: iout
1271  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: iarray
1272  character(len=*), intent(in) :: aname
1273  ! -- local
1274  integer(I4B) :: ival
1275  real(DP) :: rval
1276  integer(I4B) :: nlay
1277  integer(I4B) :: nrow
1278  integer(I4B) :: ncol
1279  integer(I4B) :: nval
1280  integer(I4B), dimension(:), pointer, contiguous :: itemp
1281  !
1282  ! -- Point the temporary pointer array, which is passed to the reading
1283  ! subroutine. The temporary array will point to ibuff if it is a
1284  ! reduced structured system, or to iarray if it is an unstructured
1285  ! model.
1286  nlay = this%mshape(1)
1287  nrow = this%mshape(2)
1288  ncol = this%mshape(3)
1289  !
1290  if (this%nodes < this%nodesuser) then
1291  nval = this%nodesuser
1292  itemp => this%ibuff
1293  else
1294  nval = this%nodes
1295  itemp => iarray
1296  end if
1297  !
1298  ! -- Read the array
1299  call urword(line, lloc, istart, istop, 1, ival, rval, iout, in)
1300  if (line(istart:istop) .EQ. 'LAYERED') then
1301  !
1302  ! -- Read layered input
1303  call readarray(in, itemp, aname, this%ndim, ncol, nrow, nlay, nval, &
1304  iout, 1, nlay)
1305  else
1306  !
1307  ! -- Read unstructured input
1308  call readarray(in, itemp, aname, this%ndim, nval, iout, 0)
1309  end if
1310  !
1311  ! -- If reduced model, then need to copy from itemp(=>ibuff) to iarray
1312  if (this%nodes < this%nodesuser) then
1313  call this%fill_grid_array(itemp, iarray)
1314  end if
1315  !
1316  end subroutine read_int_array
1317 
1318  !> @brief Read a double precision array
1319  !<
1320  subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, &
1321  darray, aname)
1322  ! -- dummy
1323  class(distype), intent(inout) :: this
1324  character(len=*), intent(inout) :: line
1325  integer(I4B), intent(inout) :: lloc
1326  integer(I4B), intent(inout) :: istart
1327  integer(I4B), intent(inout) :: istop
1328  integer(I4B), intent(in) :: in
1329  integer(I4B), intent(in) :: iout
1330  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray
1331  character(len=*), intent(in) :: aname
1332  ! -- local
1333  integer(I4B) :: ival
1334  real(DP) :: rval
1335  integer(I4B) :: nlay
1336  integer(I4B) :: nrow
1337  integer(I4B) :: ncol
1338  integer(I4B) :: nval
1339  real(DP), dimension(:), pointer, contiguous :: dtemp
1340  !
1341  ! -- Point the temporary pointer array, which is passed to the reading
1342  ! subroutine. The temporary array will point to dbuff if it is a
1343  ! reduced structured system, or to darray if it is an unstructured
1344  ! model.
1345  nlay = this%mshape(1)
1346  nrow = this%mshape(2)
1347  ncol = this%mshape(3)
1348  !
1349  if (this%nodes < this%nodesuser) then
1350  nval = this%nodesuser
1351  dtemp => this%dbuff
1352  else
1353  nval = this%nodes
1354  dtemp => darray
1355  end if
1356  !
1357  ! -- Read the array
1358  call urword(line, lloc, istart, istop, 1, ival, rval, iout, in)
1359  if (line(istart:istop) .EQ. 'LAYERED') then
1360  !
1361  ! -- Read structured input
1362  call readarray(in, dtemp, aname, this%ndim, ncol, nrow, nlay, nval, &
1363  iout, 1, nlay)
1364  else
1365  !
1366  ! -- Read unstructured input
1367  call readarray(in, dtemp, aname, this%ndim, nval, iout, 0)
1368  end if
1369  !
1370  ! -- If reduced model, then need to copy from dtemp(=>dbuff) to darray
1371  if (this%nodes < this%nodesuser) then
1372  call this%fill_grid_array(dtemp, darray)
1373  end if
1374  !
1375  end subroutine read_dbl_array
1376 
1377  !> @brief Read a 2d double array into col icolbnd of darray
1378  !!
1379  !! For cells that are outside of the active domain,
1380  !! do not copy the array value into darray.
1381  !<
1382  subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, &
1383  icolbnd, aname, inunit, iout)
1384  ! -- dummy
1385  class(distype) :: this
1386  integer(I4B), intent(in) :: maxbnd
1387  integer(I4B), dimension(maxbnd) :: nodelist
1388  integer(I4B), intent(in) :: ncolbnd
1389  real(DP), dimension(ncolbnd, maxbnd), intent(inout) :: darray
1390  integer(I4B), intent(in) :: icolbnd
1391  character(len=*), intent(in) :: aname
1392  integer(I4B), intent(in) :: inunit
1393  integer(I4B), intent(in) :: iout
1394  ! -- local
1395  integer(I4B) :: ir, ic, ncol, nrow, nlay, nval, ipos, nodeu
1396  !
1397  ! -- set variables
1398  nlay = this%mshape(1)
1399  nrow = this%mshape(2)
1400  ncol = this%mshape(3)
1401  !
1402  ! -- Read the array
1403  nval = ncol * nrow
1404  call readarray(inunit, this%dbuff, aname, this%ndim, ncol, nrow, nlay, &
1405  nval, iout, 0, 0)
1406  !
1407  ! -- Copy array into bound. Note that this routine was substantially
1408  ! changed on 9/21/2021 to support changes to READASARRAYS input
1409  ! for recharge and evapotranspiration. nodelist and bound are of
1410  ! size nrow * ncol and correspond directly to dbuff.
1411  ipos = 1
1412  do ir = 1, nrow
1413  do ic = 1, ncol
1414  nodeu = get_node(1, ir, ic, nlay, nrow, ncol)
1415  darray(icolbnd, ipos) = this%dbuff(nodeu)
1416  ipos = ipos + 1
1417  !
1418  end do
1419  end do
1420  !
1421  end subroutine read_layer_array
1422 
1423  !> @brief Record a double precision array.
1424  !!
1425  !! The array is written to a formatted or unformatted external file
1426  !! depending on the arguments.
1427  !<
1428  subroutine record_array(this, darray, iout, iprint, idataun, aname, &
1429  cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
1430  ! -- dummy
1431  class(distype), intent(inout) :: this
1432  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray !< double precision array to record
1433  integer(I4B), intent(in) :: iout !< ascii output unit number
1434  integer(I4B), intent(in) :: iprint !< whether to print the array
1435  integer(I4B), intent(in) :: idataun !< binary output unit number, if negative don't write by layers, write entire array
1436  character(len=*), intent(in) :: aname !< text descriptor
1437  character(len=*), intent(in) :: cdatafmp !< write format
1438  integer(I4B), intent(in) :: nvaluesp !< values per line
1439  integer(I4B), intent(in) :: nwidthp !< number width
1440  character(len=*), intent(in) :: editdesc !< format type (I, G, F, S, E)
1441  real(DP), intent(in) :: dinact !< double precision value for cells excluded from model domain
1442  ! -- local
1443  integer(I4B) :: k, ifirst
1444  integer(I4B) :: nlay
1445  integer(I4B) :: nrow
1446  integer(I4B) :: ncol
1447  integer(I4B) :: nval
1448  integer(I4B) :: nodeu, noder
1449  integer(I4B) :: istart, istop
1450  real(DP), dimension(:), pointer, contiguous :: dtemp
1451  ! -- formats
1452  character(len=*), parameter :: fmthsv = &
1453  "(1X,/1X,a,' WILL BE SAVED ON UNIT ',I4, &
1454  &' AT END OF TIME STEP',I5,', STRESS PERIOD ',I4)"
1455  !
1456  ! -- set variables
1457  nlay = this%mshape(1)
1458  nrow = this%mshape(2)
1459  ncol = this%mshape(3)
1460  !
1461  ! -- If this is a reduced model, then copy the values from darray into
1462  ! dtemp.
1463  if (this%nodes < this%nodesuser) then
1464  nval = this%nodes
1465  dtemp => this%dbuff
1466  do nodeu = 1, this%nodesuser
1467  noder = this%get_nodenumber(nodeu, 0)
1468  if (noder <= 0) then
1469  dtemp(nodeu) = dinact
1470  cycle
1471  end if
1472  dtemp(nodeu) = darray(noder)
1473  end do
1474  else
1475  nval = this%nodes
1476  dtemp => darray
1477  end if
1478  !
1479  ! -- Print to iout if iprint /= 0
1480  if (iprint /= 0) then
1481  istart = 1
1482  do k = 1, nlay
1483  istop = istart + nrow * ncol - 1
1484  call ulaprufw(ncol, nrow, kstp, kper, k, iout, dtemp(istart:istop), &
1485  aname, cdatafmp, nvaluesp, nwidthp, editdesc)
1486  istart = istop + 1
1487  end do
1488  end if
1489  !
1490  ! -- Save array to an external file.
1491  if (idataun > 0) then
1492  ! -- write to binary file by layer
1493  ifirst = 1
1494  istart = 1
1495  do k = 1, nlay
1496  istop = istart + nrow * ncol - 1
1497  if (ifirst == 1) write (iout, fmthsv) &
1498  trim(adjustl(aname)), idataun, &
1499  kstp, kper
1500  ifirst = 0
1501  call ulasav(dtemp(istart:istop), aname, kstp, kper, &
1502  pertim, totim, ncol, nrow, k, idataun)
1503  istart = istop + 1
1504  end do
1505  elseif (idataun < 0) then
1506  !
1507  ! -- write entire array as one record
1508  call ubdsv1(kstp, kper, aname, -idataun, dtemp, ncol, nrow, nlay, &
1509  iout, delt, pertim, totim)
1510  end if
1511  !
1512  end subroutine record_array
1513 
1514  !> @brief Record list header for imeth=6
1515  !<
1516  subroutine record_srcdst_list_header(this, text, textmodel, textpackage, &
1517  dstmodel, dstpackage, naux, auxtxt, &
1518  ibdchn, nlist, iout)
1519  ! -- dummy
1520  class(distype) :: this
1521  character(len=16), intent(in) :: text
1522  character(len=16), intent(in) :: textmodel
1523  character(len=16), intent(in) :: textpackage
1524  character(len=16), intent(in) :: dstmodel
1525  character(len=16), intent(in) :: dstpackage
1526  integer(I4B), intent(in) :: naux
1527  character(len=16), dimension(:), intent(in) :: auxtxt
1528  integer(I4B), intent(in) :: ibdchn
1529  integer(I4B), intent(in) :: nlist
1530  integer(I4B), intent(in) :: iout
1531  ! -- local
1532  integer(I4B) :: nlay, nrow, ncol
1533  !
1534  nlay = this%mshape(1)
1535  nrow = this%mshape(2)
1536  ncol = this%mshape(3)
1537  !
1538  ! -- Use ubdsv06 to write list header
1539  call ubdsv06(kstp, kper, text, textmodel, textpackage, dstmodel, dstpackage, &
1540  ibdchn, naux, auxtxt, ncol, nrow, nlay, &
1541  nlist, iout, delt, pertim, totim)
1542  !
1543  end subroutine record_srcdst_list_header
1544 
1545  !> @brief Convert an integer array (layer numbers) to nodelist
1546  !<
1547  subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
1548  ! -- dummy
1549  class(distype) :: this
1550  integer(I4B), intent(in) :: maxbnd
1551  integer(I4B), dimension(:), pointer, contiguous :: darray
1552  integer(I4B), dimension(maxbnd), intent(inout) :: nodelist
1553  integer(I4B), intent(inout) :: nbound
1554  character(len=*), intent(in) :: aname
1555  ! -- local
1556  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nval, nodeu, noder, ipos, ierr
1557  !
1558  ! -- set variables
1559  nlay = this%mshape(1)
1560  nrow = this%mshape(2)
1561  ncol = this%mshape(3)
1562  !
1563  if (this%ndim > 1) then
1564  !
1565  nval = ncol * nrow
1566  !
1567  ! -- Copy array into nodelist
1568  ipos = 1
1569  ierr = 0
1570  do ir = 1, nrow
1571  do ic = 1, ncol
1572  nodeu = get_node(1, ir, ic, nlay, nrow, ncol)
1573  il = darray(nodeu)
1574  if (il < 1 .or. il > nlay) then
1575  write (errmsg, '(a,1x,i0)') 'Invalid layer number:', il
1576  call store_error(errmsg, terminate=.true.)
1577  end if
1578  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
1579  noder = this%get_nodenumber(nodeu, 0)
1580  if (ipos > maxbnd) then
1581  ierr = ipos
1582  else
1583  nodelist(ipos) = noder
1584  end if
1585  ipos = ipos + 1
1586  end do
1587  end do
1588  !
1589  ! -- Check for errors
1590  nbound = ipos - 1
1591  if (ierr > 0) then
1592  write (errmsg, '(a,1x,i0)') &
1593  'MAXBOUND dimension is too small.'// &
1594  'INCREASE MAXBOUND TO:', ierr
1595  call store_error(errmsg, terminate=.true.)
1596  end if
1597  !
1598  ! -- If nbound < maxbnd, then initialize nodelist to zero in this range
1599  if (nbound < maxbnd) then
1600  do ipos = nbound + 1, maxbnd
1601  nodelist(ipos) = 0
1602  end do
1603  end if
1604  !
1605  else
1606  !
1607  ! -- For unstructured, read nodelist directly, then check node numbers
1608  nodelist = darray
1609  do noder = 1, maxbnd
1610  if (noder < 1 .or. noder > this%nodes) then
1611  write (errmsg, '(a,1x,i0)') 'Invalid node number:', noder
1612  call store_error(errmsg, terminate=.true.)
1613  end if
1614  end do
1615  nbound = maxbnd
1616  !
1617  end if
1618  !
1619  end subroutine nlarray_to_nodelist
1620 
1621 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:45
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15
@ dis
DIS6 discretization.
Definition: Constants.f90:155
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:68
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
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
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:1057
subroutine nodeu_to_array(this, nodeu, arr)
Convert a user nodenumber to an array (nodenumber) or (k,i,j)
Definition: Dis.f90:712
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:1204
subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
Convert an integer array (layer numbers) to nodelist.
Definition: Dis.f90:1548
subroutine record_array(this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
Record a double precision array.
Definition: Dis.f90:1430
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:963
subroutine allocate_scalars(this, name_model, input_mempath)
Allocate and initialize scalar variables.
Definition: Dis.f90:817
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:773
logical function supports_layers(this)
Indicates whether the grid discretization supports layers.
Definition: Dis.f90:1047
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:1123
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:872
subroutine get_dis_type(this, dis_type)
Get the discretization type.
Definition: Dis.f90:1194
subroutine nodeu_to_string(this, nodeu, str)
Convert a user nodenumber to a string (nodenumber) or (k,i,j)
Definition: Dis.f90:691
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:1384
subroutine allocate_arrays(this)
Allocate and initialize arrays.
Definition: Dis.f90:841
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:1068
integer(i4b) function get_nodenumber_idx1(this, nodeu, icheck)
Get reduced node number from user node number.
Definition: Dis.f90:742
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:1263
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:1215
subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, darray, aname)
Read a double precision array.
Definition: Dis.f90:1322
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:1519
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:475
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:83
subroutine, public get_ijk(nodenumber, nrow, ncol, nlay, irow, icol, ilay)
Get row, column and layer indices from node number and grid dimensions. If nodenumber is invalid,...
Definition: GeomUtil.f90:100
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 memorystore_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