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