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