MODFLOW 6  version 6.7.0.dev0
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  ! -- dummy
486  class(dis2dtype) :: this
487  integer(I4B), dimension(:), intent(in) :: icelltype
488  ! -- local
489  integer(I4B) :: iunit, ntxt
490  integer(I4B), parameter :: lentxt = 100
491  character(len=50) :: txthdr
492  character(len=lentxt) :: txt
493  character(len=LINELENGTH) :: fname
494  character(len=*), parameter :: fmtgrdsave = &
495  "(4X,'BINARY GRID INFORMATION WILL BE WRITTEN TO:', &
496  &/,6X,'UNIT NUMBER: ', I0,/,6X, 'FILE NAME: ', A)"
497  !
498  ! -- Initialize
499  ntxt = 14
500  !
501  ! -- Open the file
502  fname = trim(this%input_fname)//'.grb'
503  iunit = getunit()
504  write (this%iout, fmtgrdsave) iunit, trim(adjustl(fname))
505  call openfile(iunit, this%iout, trim(adjustl(fname)), 'DATA(BINARY)', &
506  form, access, 'REPLACE')
507  !
508  ! -- write header information
509  write (txthdr, '(a)') 'GRID DIS2D'
510  txthdr(50:50) = new_line('a')
511  write (iunit) txthdr
512  write (txthdr, '(a)') 'VERSION 1'
513  txthdr(50:50) = new_line('a')
514  write (iunit) txthdr
515  write (txthdr, '(a, i0)') 'NTXT ', ntxt
516  txthdr(50:50) = new_line('a')
517  write (iunit) txthdr
518  write (txthdr, '(a, i0)') 'LENTXT ', lentxt
519  txthdr(50:50) = new_line('a')
520  write (iunit) txthdr
521  !
522  ! -- write variable definitions
523  write (txt, '(3a, i0)') 'NCELLS ', 'INTEGER ', 'NDIM 0 # ', this%nodesuser
524  txt(lentxt:lentxt) = new_line('a')
525  write (iunit) txt
526  write (txt, '(3a, i0)') 'NROW ', 'INTEGER ', 'NDIM 0 # ', this%nrow
527  txt(lentxt:lentxt) = new_line('a')
528  write (iunit) txt
529  write (txt, '(3a, i0)') 'NCOL ', 'INTEGER ', 'NDIM 0 # ', this%ncol
530  txt(lentxt:lentxt) = new_line('a')
531  write (iunit) txt
532  write (txt, '(3a, i0)') 'NJA ', 'INTEGER ', 'NDIM 0 # ', this%nja
533  txt(lentxt:lentxt) = new_line('a')
534  write (iunit) txt
535  write (txt, '(3a, 1pg24.15)') 'XORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%xorigin
536  txt(lentxt:lentxt) = new_line('a')
537  write (iunit) txt
538  write (txt, '(3a, 1pg24.15)') 'YORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%yorigin
539  txt(lentxt:lentxt) = new_line('a')
540  write (iunit) txt
541  write (txt, '(3a, 1pg24.15)') 'ANGROT ', 'DOUBLE ', 'NDIM 0 # ', this%angrot
542  txt(lentxt:lentxt) = new_line('a')
543  write (iunit) txt
544  write (txt, '(3a, i0)') 'DELR ', 'DOUBLE ', 'NDIM 1 ', this%ncol
545  txt(lentxt:lentxt) = new_line('a')
546  write (iunit) txt
547  write (txt, '(3a, i0)') 'DELC ', 'DOUBLE ', 'NDIM 1 ', this%nrow
548  txt(lentxt:lentxt) = new_line('a')
549  write (iunit) txt
550  write (txt, '(3a, i0)') 'BOTM ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
551  txt(lentxt:lentxt) = new_line('a')
552  write (iunit) txt
553  write (txt, '(3a, i0)') 'IA ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
554  txt(lentxt:lentxt) = new_line('a')
555  write (iunit) txt
556  write (txt, '(3a, i0)') 'JA ', 'INTEGER ', 'NDIM 1 ', size(this%con%jausr)
557  txt(lentxt:lentxt) = new_line('a')
558  write (iunit) txt
559  write (txt, '(3a, i0)') 'IDOMAIN ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
560  txt(lentxt:lentxt) = new_line('a')
561  write (iunit) txt
562  write (txt, '(3a, i0)') 'ICELLTYPE ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
563  txt(lentxt:lentxt) = new_line('a')
564  write (iunit) txt
565  !
566  ! -- write data
567  write (iunit) this%nodesuser ! ncells
568  write (iunit) this%nrow ! nrow
569  write (iunit) this%ncol ! ncol
570  write (iunit) this%nja ! nja
571  write (iunit) this%xorigin ! xorigin
572  write (iunit) this%yorigin ! yorigin
573  write (iunit) this%angrot ! angrot
574  write (iunit) this%delr ! delr
575  write (iunit) this%delc ! delc
576  write (iunit) this%bottom ! bottom
577  write (iunit) this%con%iausr ! iausr
578  write (iunit) this%con%jausr ! jausr
579  write (iunit) this%idomain ! idomain
580  write (iunit) icelltype ! icelltype
581  !
582  ! -- Close the file
583  close (iunit)
584  !
585  end subroutine write_grb
586 
587  !> @brief Convert a user nodenumber to a string (nodenumber) or (k,i,j)
588  !<
589  subroutine nodeu_to_string(this, nodeu, str)
590  ! -- dummy
591  class(dis2dtype) :: this
592  integer(I4B), intent(in) :: nodeu
593  character(len=*), intent(inout) :: str
594  ! -- local
595  integer(I4B) :: i, j, k
596  character(len=10) :: istr, jstr
597  !
598  call get_ijk(nodeu, this%nrow, this%ncol, 1, i, j, k)
599  write (istr, '(i10)') i
600  write (jstr, '(i10)') j
601  str = '('//trim(adjustl(istr))//','// &
602  trim(adjustl(jstr))//')'
603  !
604  end subroutine nodeu_to_string
605 
606  !> @brief Convert a user nodenumber to an array (nodenumber) or (i,j)
607  !<
608  subroutine nodeu_to_array(this, nodeu, arr)
609  ! -- dummy
610  class(dis2dtype) :: this
611  integer(I4B), intent(in) :: nodeu
612  integer(I4B), dimension(:), intent(inout) :: arr
613  ! -- local
614  integer(I4B) :: isize
615  integer(I4B) :: i, j, k
616  !
617  ! -- check the size of arr
618  isize = size(arr)
619  if (isize /= this%ndim) then
620  write (errmsg, '(a,i0,a,i0,a)') &
621  'Program error: nodeu_to_array size of array (', isize, &
622  ') is not equal to the discretization dimension (', this%ndim, ')'
623  call store_error(errmsg, terminate=.true.)
624  end if
625  !
626  ! -- get k, i, j
627  call get_ijk(nodeu, this%nrow, this%ncol, 1, i, j, k)
628  !
629  ! -- fill array
630  arr(1) = i
631  arr(2) = j
632  !
633  end subroutine nodeu_to_array
634 
635  !> @brief Get reduced node number from user node number
636  !<
637  function get_nodenumber_idx1(this, nodeu, icheck) result(nodenumber)
638  ! -- return
639  integer(I4B) :: nodenumber
640  ! -- dummy
641  class(dis2dtype), intent(in) :: this
642  integer(I4B), intent(in) :: nodeu
643  integer(I4B), intent(in) :: icheck
644  !
645  ! -- check the node number if requested
646  if (icheck /= 0) then
647  !
648  ! -- If within valid range, convert to reduced nodenumber
649  if (nodeu < 1 .or. nodeu > this%nodesuser) then
650  write (errmsg, '(a,i0,a)') &
651  'Node number (', nodeu, &
652  ') less than 1 or greater than the number of nodes.'
653  call store_error(errmsg)
654  nodenumber = 0
655  else
656  nodenumber = nodeu
657  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
658  end if
659  else
660  nodenumber = nodeu
661  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
662  end if
663  !
664  end function get_nodenumber_idx1
665 
666  !> @brief Get reduced node number from layer, row and column indices
667  !<
668  function get_nodenumber_idx2(this, k, j, icheck) result(nodenumber)
669  ! -- return
670  integer(I4B) :: nodenumber
671  ! -- dummy
672  class(dis2dtype), intent(in) :: this
673  integer(I4B), intent(in) :: k, j
674  integer(I4B), intent(in) :: icheck
675  ! -- local
676  integer(I4B) :: nodeu, i
677  ! formats
678  character(len=*), parameter :: fmterr = &
679  "('Error in structured-grid cell indices: row = ',i0,&
680  &', column = ',i0)"
681  !
682  i = k
683  nodeu = get_node(1, i, j, 1, this%nrow, this%ncol)
684  if (nodeu < 1) then
685  write (errmsg, fmterr) i, j
686  call store_error(errmsg, terminate=.true.)
687  end if
688  nodenumber = nodeu
689  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
690  !
691  ! -- check the node number if requested
692  if (icheck /= 0) then
693  !
694  if (i < 1 .or. i > this%nrow) &
695  call store_error('Row less than one or greater than nrow')
696  if (j < 1 .or. j > this%ncol) &
697  call store_error('Column less than one or greater than ncol')
698  !
699  ! -- Error if outside of range
700  if (nodeu < 1 .or. nodeu > this%nodesuser) then
701  write (errmsg, '(a,i0,a)') &
702  'Node number (', nodeu, ')less than 1 or greater than nodes.'
703  call store_error(errmsg)
704  end if
705  end if
706  !
707  end function get_nodenumber_idx2
708 
709  !> @brief Allocate and initialize scalar variables
710  !<
711  subroutine allocate_scalars(this, name_model, input_mempath)
712  ! -- dummy
713  class(dis2dtype) :: this
714  character(len=*), intent(in) :: name_model
715  character(len=*), intent(in) :: input_mempath
716  !
717  ! -- Allocate parent scalars
718  call this%DisBaseType%allocate_scalars(name_model, input_mempath)
719  !
720  ! -- Allocate
721  call mem_allocate(this%nrow, 'NROW', this%memoryPath)
722  call mem_allocate(this%ncol, 'NCOL', this%memoryPath)
723  !
724  ! -- Initialize
725  this%nrow = 0
726  this%ncol = 0
727  this%ndim = 2
728  !
729  end subroutine allocate_scalars
730 
731  !> @brief Allocate and initialize arrays
732  !<
733  subroutine allocate_arrays(this)
734  ! -- dummy
735  class(dis2dtype) :: this
736  !
737  ! -- Allocate arrays in DisBaseType (mshape, top, bot, area)
738  call this%DisBaseType%allocate_arrays()
739  !
740  ! -- Allocate arrays for DisType
741  if (this%nodes < this%nodesuser) then
742  call mem_allocate(this%nodeuser, this%nodes, 'NODEUSER', this%memoryPath)
743  call mem_allocate(this%nodereduced, this%nodesuser, 'NODEREDUCED', &
744  this%memoryPath)
745  else
746  call mem_allocate(this%nodeuser, 1, 'NODEUSER', this%memoryPath)
747  call mem_allocate(this%nodereduced, 1, 'NODEREDUCED', this%memoryPath)
748  end if
749  !
750  ! -- Initialize
751  this%mshape(1) = this%nrow
752  this%mshape(2) = this%ncol
753  !
754  end subroutine allocate_arrays
755 
756  !> @brief Convert a string to a user nodenumber
757  !!
758  !! Parse layer, row and column and return user nodenumber.
759  !! If flag_string is present and true, the first token may be
760  !! non-numeric (e.g. boundary name). In this case, return -2.
761  !<
762  function nodeu_from_string(this, lloc, istart, istop, in, iout, line, &
763  flag_string, allow_zero) result(nodeu)
764  ! -- dummy
765  class(dis2dtype) :: this
766  integer(I4B), intent(inout) :: lloc
767  integer(I4B), intent(inout) :: istart
768  integer(I4B), intent(inout) :: istop
769  integer(I4B), intent(in) :: in
770  integer(I4B), intent(in) :: iout
771  character(len=*), intent(inout) :: line
772  logical, optional, intent(in) :: flag_string
773  logical, optional, intent(in) :: allow_zero
774  integer(I4B) :: nodeu
775  ! -- local
776  integer(I4B) :: i, j, nrow, ncol
777  integer(I4B) :: lloclocal, ndum, istat, n
778  real(dp) :: r
779  !
780  if (present(flag_string)) then
781  if (flag_string) then
782  ! Check to see if first token in line can be read as an integer.
783  lloclocal = lloc
784  call urword(line, lloclocal, istart, istop, 1, ndum, r, iout, in)
785  read (line(istart:istop), *, iostat=istat) n
786  if (istat /= 0) then
787  ! First token in line is not an integer; return flag to this effect.
788  nodeu = -2
789  return
790  end if
791  end if
792  end if
793  !
794  nrow = this%mshape(1)
795  ncol = this%mshape(2)
796  !
797  call urword(line, lloc, istart, istop, 2, i, r, iout, in)
798  call urword(line, lloc, istart, istop, 2, j, r, iout, in)
799  !
800  if (i == 0 .and. j == 0) then
801  if (present(allow_zero)) then
802  if (allow_zero) then
803  nodeu = 0
804  return
805  end if
806  end if
807  end if
808  !
809  errmsg = ""
810  !
811  if (i < 1 .or. i > nrow) then
812  write (errmsg, '(a,1x,a,i0,a)') &
813  trim(adjustl(errmsg)), 'Row number in list (', i, &
814  ') is outside of the grid.'
815  end if
816  if (j < 1 .or. j > ncol) then
817  write (errmsg, '(a,1x,a,i0,a)') &
818  trim(adjustl(errmsg)), 'Column number in list (', j, &
819  ') is outside of the grid.'
820  end if
821  !
822  nodeu = get_node(1, i, j, 1, nrow, ncol)
823  !
824  if (nodeu < 1 .or. nodeu > this%nodesuser) then
825  write (errmsg, '(a,1x,a,i0,a)') &
826  trim(adjustl(errmsg)), &
827  "Node number in list (", nodeu, ") is outside of the grid. "// &
828  "Cell number cannot be determined in line '"// &
829  trim(adjustl(line))//"'."
830  end if
831  !
832  if (len_trim(adjustl(errmsg)) > 0) then
833  call store_error(errmsg)
834  call store_error_unit(in)
835  end if
836  !
837  end function nodeu_from_string
838 
839  !> @brief Convert a cellid string to a user nodenumber
840  !!
841  !! If flag_string is present and true, the first token may be
842  !! non-numeric (e.g. boundary name). In this case, return -2.
843  !!
844  !! If allow_zero is present and true, and all indices are zero, the
845  !! result can be zero. If allow_zero is false, a zero in any index is an error.
846  !<
847  function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, &
848  allow_zero) result(nodeu)
849  ! -- return
850  integer(I4B) :: nodeu
851  ! -- dummy
852  class(dis2dtype) :: this
853  character(len=*), intent(inout) :: cellid
854  integer(I4B), intent(in) :: inunit
855  integer(I4B), intent(in) :: iout
856  logical, optional, intent(in) :: flag_string
857  logical, optional, intent(in) :: allow_zero
858  ! -- local
859  integer(I4B) :: lloclocal, istart, istop, ndum, n
860  integer(I4B) :: i, j, nrow, ncol
861  integer(I4B) :: istat
862  real(dp) :: r
863  !
864  if (present(flag_string)) then
865  if (flag_string) then
866  ! Check to see if first token in cellid can be read as an integer.
867  lloclocal = 1
868  call urword(cellid, lloclocal, istart, istop, 1, ndum, r, iout, inunit)
869  read (cellid(istart:istop), *, iostat=istat) n
870  if (istat /= 0) then
871  ! First token in cellid is not an integer; return flag to this effect.
872  nodeu = -2
873  return
874  end if
875  end if
876  end if
877  !
878  nrow = this%mshape(1)
879  ncol = this%mshape(2)
880  !
881  lloclocal = 1
882  call urword(cellid, lloclocal, istart, istop, 2, i, r, iout, inunit)
883  call urword(cellid, lloclocal, istart, istop, 2, j, r, iout, inunit)
884  !
885  if (i == 0 .and. j == 0) then
886  if (present(allow_zero)) then
887  if (allow_zero) then
888  nodeu = 0
889  return
890  end if
891  end if
892  end if
893  !
894  errmsg = ""
895  !
896  if (i < 1 .or. i > nrow) then
897  write (errmsg, '(a,1x,a,i0,a)') &
898  trim(adjustl(errmsg)), 'Row number in list (', i, &
899  ') is outside of the grid.'
900  end if
901  if (j < 1 .or. j > ncol) then
902  write (errmsg, '(a,1x,a,i0,a)') &
903  trim(adjustl(errmsg)), 'Column number in list (', j, &
904  ') is outside of the grid.'
905  end if
906  !
907  nodeu = get_node(1, i, j, 1, nrow, ncol)
908  !
909  if (nodeu < 1 .or. nodeu > this%nodesuser) then
910  write (errmsg, '(a,1x,a,i0,a)') &
911  trim(adjustl(errmsg)), &
912  "Cell number cannot be determined for cellid ("// &
913  trim(adjustl(cellid))//") and results in a user "// &
914  "node number (", nodeu, ") that is outside of the grid."
915  end if
916  !
917  if (len_trim(adjustl(errmsg)) > 0) then
918  call store_error(errmsg)
919  call store_error_unit(inunit)
920  end if
921  !
922  end function nodeu_from_cellid
923 
924  !> @brief Indicates whether the grid discretization supports layers
925  !<
926  logical function supports_layers(this)
927  ! -- dummy
928  class(dis2dtype) :: this
929  !
930  supports_layers = .true.
931  !
932  end function supports_layers
933 
934  !> @brief Return number of cells per layer (nrow * ncol)
935  !<
936  function get_ncpl(this)
937  integer(I4B) :: get_ncpl
938  class(dis2dtype) :: this
939  get_ncpl = this%nrow * this%ncol
940  end function get_ncpl
941 
942  !> @brief Get normal vector components between the cell and a given neighbor
943  !!
944  !! The normal points outward from the shared face between noden and nodem.
945  !<
946  subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, &
947  ipos)
948  ! -- dummy
949  class(dis2dtype) :: this
950  integer(I4B), intent(in) :: noden !< cell (reduced nn)
951  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
952  integer(I4B), intent(in) :: ihc !< horizontal connection flag
953  real(DP), intent(inout) :: xcomp
954  real(DP), intent(inout) :: ycomp
955  real(DP), intent(inout) :: zcomp
956  integer(I4B), intent(in) :: ipos
957  ! -- local
958  integer(I4B) :: nodeu1, i1, j1, k1
959  integer(I4B) :: nodeu2, i2, j2, k2
960  !
961  ! -- Set vector components based on ihc
962  if (ihc == 0) then
963  xcomp = dzero
964  ycomp = dzero
965  if (nodem < noden) then
966  !
967  ! -- nodem must be above noden, so upward connection
968  zcomp = done
969  else
970  !
971  ! -- nodem must be below noden, so downward connection
972  zcomp = -done
973  end if
974  else
975  xcomp = dzero
976  ycomp = dzero
977  zcomp = dzero
978  nodeu1 = this%get_nodeuser(noden)
979  nodeu2 = this%get_nodeuser(nodem)
980  call get_ijk(nodeu1, this%nrow, this%ncol, 1, i1, j1, k1)
981  call get_ijk(nodeu2, this%nrow, this%ncol, 1, i2, j2, k2)
982  if (i2 < i1) then ! back
983  ycomp = done
984  elseif (j2 < j1) then ! left
985  xcomp = -done
986  elseif (j2 > j1) then ! right
987  xcomp = done
988  else ! front
989  ycomp = -done
990  end if
991  !
992  end if
993  !
994  end subroutine connection_normal
995 
996  !> @brief Get unit vector components between the cell and a given neighbor
997  !!
998  !< For dis2d there is no z component
999  subroutine connection_vector(this, noden, nodem, nozee, satn, satm, ihc, &
1000  xcomp, ycomp, zcomp, conlen)
1001  ! modules
1002  use disvgeom, only: line_unit_vector
1003  ! dummy
1004  class(dis2dtype) :: this
1005  integer(I4B), intent(in) :: noden !< cell (reduced nn)
1006  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
1007  logical, intent(in) :: nozee !< not used for dis2d
1008  real(DP), intent(in) :: satn !< not used for dis2d
1009  real(DP), intent(in) :: satm !< not used for dis2d
1010  integer(I4B), intent(in) :: ihc !< not used for dis2d (always horizontal)
1011  real(DP), intent(inout) :: xcomp !< x component of the connection vector
1012  real(DP), intent(inout) :: ycomp !< y component of the connection vector
1013  real(DP), intent(inout) :: zcomp !< z component, which is always zero
1014  real(DP), intent(inout) :: conlen !< calculated connection length
1015  ! local
1016  real(DP) :: z1, z2
1017  real(DP) :: x1, y1, x2, y2
1018  real(DP) :: ds
1019  integer(I4B) :: i1, i2, j1, j2, k1, k2
1020  integer(I4B) :: nodeu1, nodeu2, ipos
1021 
1022  ! Calculate vector components
1023  z1 = dzero
1024  z2 = dzero
1025  ipos = this%con%getjaindex(noden, nodem)
1026  ds = this%con%cl1(this%con%jas(ipos)) + this%con%cl2(this%con%jas(ipos))
1027  nodeu1 = this%get_nodeuser(noden)
1028  nodeu2 = this%get_nodeuser(nodem)
1029  call get_ijk(nodeu1, this%nrow, this%ncol, 1, i1, j1, k1)
1030  call get_ijk(nodeu2, this%nrow, this%ncol, 1, i2, j2, k2)
1031  x1 = dzero
1032  x2 = dzero
1033  y1 = dzero
1034  y2 = dzero
1035  if (i2 < i1) then ! back
1036  y2 = ds
1037  elseif (j2 < j1) then ! left
1038  x2 = -ds
1039  elseif (j2 > j1) then ! right
1040  x2 = ds
1041  else ! front
1042  y2 = -ds
1043  end if
1044  call line_unit_vector(x1, y1, z1, x2, y2, z2, xcomp, ycomp, zcomp, conlen)
1045  end subroutine connection_vector
1046 
1047  !> @brief Get the discretization type
1048  !<
1049  subroutine get_dis_type(this, dis_type)
1050  ! -- dummy
1051  class(dis2dtype), intent(in) :: this
1052  character(len=*), intent(out) :: dis_type
1053  !
1054  dis_type = "DIS2D"
1055  !
1056  end subroutine get_dis_type
1057 
1058  !> @brief Get the discretization type enumeration
1059  function get_dis_enum(this) result(dis_enum)
1060  use constantsmodule, only: dis2d
1061  class(dis2dtype), intent(in) :: this
1062  integer(I4B) :: dis_enum
1063  dis_enum = dis2d
1064  end function get_dis_enum
1065 
1066  !> @brief Get a 2D array of polygon vertices, listed in
1067  !!
1068  !! clockwise order beginning with the lower left corner
1069  !<
1070  subroutine get_polyverts(this, ic, polyverts, closed)
1071  ! -- dummy
1072  class(dis2dtype), intent(inout) :: this
1073  integer(I4B), intent(in) :: ic !< cell number (reduced)
1074  real(DP), allocatable, intent(out) :: polyverts(:, :) !< polygon vertices (column-major indexing)
1075  logical(LGP), intent(in), optional :: closed !< whether to close the polygon, duplicating a vertex
1076  ! -- local
1077  integer(I4B) :: icu, nverts, irow, jcol, klay
1078  real(DP) :: cellx, celly, dxhalf, dyhalf
1079  logical(LGP) :: lclosed
1080  !
1081  nverts = 4
1082  !
1083  ! check closed option
1084  if (.not. (present(closed))) then
1085  lclosed = .false.
1086  else
1087  lclosed = closed
1088  end if
1089  !
1090  ! allocate vertices array
1091  if (lclosed) then
1092  allocate (polyverts(2, nverts + 1))
1093  else
1094  allocate (polyverts(2, nverts))
1095  end if
1096  !
1097  ! set vertices
1098  icu = this%get_nodeuser(ic)
1099  call get_ijk(icu, this%nrow, this%ncol, 1, irow, jcol, klay)
1100  cellx = this%cellx(jcol)
1101  celly = this%celly(irow)
1102  dxhalf = dhalf * this%delr(jcol)
1103  dyhalf = dhalf * this%delc(irow)
1104  polyverts(:, 1) = (/cellx - dxhalf, celly - dyhalf/) ! SW
1105  polyverts(:, 2) = (/cellx - dxhalf, celly + dyhalf/) ! NW
1106  polyverts(:, 3) = (/cellx + dxhalf, celly + dyhalf/) ! NE
1107  polyverts(:, 4) = (/cellx + dxhalf, celly - dyhalf/) ! SE
1108  !
1109  ! close if enabled
1110  if (lclosed) &
1111  polyverts(:, nverts + 1) = polyverts(:, 1)
1112  !
1113  end subroutine
1114 
1115  !> @brief Read an integer array
1116  !< TODO: REMOVE?
1117  subroutine read_int_array(this, line, lloc, istart, istop, iout, in, &
1118  iarray, aname)
1119  ! -- dummy
1120  class(dis2dtype), intent(inout) :: this
1121  character(len=*), intent(inout) :: line
1122  integer(I4B), intent(inout) :: lloc
1123  integer(I4B), intent(inout) :: istart
1124  integer(I4B), intent(inout) :: istop
1125  integer(I4B), intent(in) :: in
1126  integer(I4B), intent(in) :: iout
1127  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: iarray
1128  character(len=*), intent(in) :: aname
1129  ! -- local
1130 
1131  !
1132  end subroutine read_int_array
1133 
1134  !> @brief Read a double precision array
1135  !< TODO: REMOVE?
1136  subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, &
1137  darray, aname)
1138  ! -- dummy
1139  class(dis2dtype), intent(inout) :: this
1140  character(len=*), intent(inout) :: line
1141  integer(I4B), intent(inout) :: lloc
1142  integer(I4B), intent(inout) :: istart
1143  integer(I4B), intent(inout) :: istop
1144  integer(I4B), intent(in) :: in
1145  integer(I4B), intent(in) :: iout
1146  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray
1147  character(len=*), intent(in) :: aname
1148  !
1149  end subroutine read_dbl_array
1150 
1151  !> @brief Read a 2d double array into col icolbnd of darray
1152  !!
1153  !! For cells that are outside of the active domain,
1154  !! do not copy the array value into darray.
1155  !< TODO: REMOVE?
1156  subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, &
1157  icolbnd, aname, inunit, iout)
1158  ! -- dummy
1159  class(dis2dtype) :: this
1160  integer(I4B), intent(in) :: maxbnd
1161  integer(I4B), dimension(maxbnd) :: nodelist
1162  integer(I4B), intent(in) :: ncolbnd
1163  real(DP), dimension(ncolbnd, maxbnd), intent(inout) :: darray
1164  integer(I4B), intent(in) :: icolbnd
1165  character(len=*), intent(in) :: aname
1166  integer(I4B), intent(in) :: inunit
1167  integer(I4B), intent(in) :: iout
1168  !
1169  end subroutine read_layer_array
1170 
1171  !> @brief Record a double precision array.
1172  !!
1173  !! The array is written to a formatted or unformatted external file
1174  !! depending on the arguments.
1175  !<
1176  subroutine record_array(this, darray, iout, iprint, idataun, aname, &
1177  cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
1178  ! -- dummy
1179  class(dis2dtype), intent(inout) :: this
1180  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray !< double precision array to record
1181  integer(I4B), intent(in) :: iout !< ascii output unit number
1182  integer(I4B), intent(in) :: iprint !< whether to print the array
1183  integer(I4B), intent(in) :: idataun !< binary output unit number, if negative don't write by layers, write entire array
1184  character(len=*), intent(in) :: aname !< text descriptor
1185  character(len=*), intent(in) :: cdatafmp !< write format
1186  integer(I4B), intent(in) :: nvaluesp !< values per line
1187  integer(I4B), intent(in) :: nwidthp !< number width
1188  character(len=*), intent(in) :: editdesc !< format type (I, G, F, S, E)
1189  real(DP), intent(in) :: dinact !< double precision value for cells excluded from model domain
1190  ! -- local
1191  integer(I4B) :: k, ifirst
1192  integer(I4B) :: nlay
1193  integer(I4B) :: nrow
1194  integer(I4B) :: ncol
1195  integer(I4B) :: nval
1196  integer(I4B) :: nodeu, noder
1197  integer(I4B) :: istart, istop
1198  real(DP), dimension(:), pointer, contiguous :: dtemp
1199  ! -- formats
1200  character(len=*), parameter :: fmthsv = &
1201  "(1X,/1X,a,' WILL BE SAVED ON UNIT ',I4, &
1202  &' AT END OF TIME STEP',I5,', STRESS PERIOD ',I4)"
1203  !
1204  ! -- set variables
1205  nlay = 1
1206  nrow = this%mshape(1)
1207  ncol = this%mshape(2)
1208  !
1209  ! -- If this is a reduced model, then copy the values from darray into
1210  ! dtemp.
1211  if (this%nodes < this%nodesuser) then
1212  nval = this%nodes
1213  dtemp => this%dbuff
1214  do nodeu = 1, this%nodesuser
1215  noder = this%get_nodenumber(nodeu, 0)
1216  if (noder <= 0) then
1217  dtemp(nodeu) = dinact
1218  cycle
1219  end if
1220  dtemp(nodeu) = darray(noder)
1221  end do
1222  else
1223  nval = this%nodes
1224  dtemp => darray
1225  end if
1226  !
1227  ! -- Print to iout if iprint /= 0
1228  if (iprint /= 0) then
1229  istart = 1
1230  do k = 1, nlay
1231  istop = istart + nrow * ncol - 1
1232  call ulaprufw(ncol, nrow, kstp, kper, k, iout, dtemp(istart:istop), &
1233  aname, cdatafmp, nvaluesp, nwidthp, editdesc)
1234  istart = istop + 1
1235  end do
1236  end if
1237  !
1238  ! -- Save array to an external file.
1239  if (idataun > 0) then
1240  ! -- write to binary file by layer
1241  ifirst = 1
1242  istart = 1
1243  do k = 1, nlay
1244  istop = istart + nrow * ncol - 1
1245  if (ifirst == 1) write (iout, fmthsv) &
1246  trim(adjustl(aname)), idataun, &
1247  kstp, kper
1248  ifirst = 0
1249  call ulasav(dtemp(istart:istop), aname, kstp, kper, &
1250  pertim, totim, ncol, nrow, k, idataun)
1251  istart = istop + 1
1252  end do
1253  elseif (idataun < 0) then
1254  !
1255  ! -- write entire array as one record
1256  call ubdsv1(kstp, kper, aname, -idataun, dtemp, ncol, nrow, nlay, &
1257  iout, delt, pertim, totim)
1258  end if
1259  !
1260  end subroutine record_array
1261 
1262  !> @brief Record list header for imeth=6
1263  !<
1264  subroutine record_srcdst_list_header(this, text, textmodel, textpackage, &
1265  dstmodel, dstpackage, naux, auxtxt, &
1266  ibdchn, nlist, iout)
1267  ! -- dummy
1268  class(dis2dtype) :: this
1269  character(len=16), intent(in) :: text
1270  character(len=16), intent(in) :: textmodel
1271  character(len=16), intent(in) :: textpackage
1272  character(len=16), intent(in) :: dstmodel
1273  character(len=16), intent(in) :: dstpackage
1274  integer(I4B), intent(in) :: naux
1275  character(len=16), dimension(:), intent(in) :: auxtxt
1276  integer(I4B), intent(in) :: ibdchn
1277  integer(I4B), intent(in) :: nlist
1278  integer(I4B), intent(in) :: iout
1279  ! -- local
1280  integer(I4B) :: nlay, nrow, ncol
1281  !
1282  nlay = 1
1283  nrow = this%mshape(1)
1284  ncol = this%mshape(2)
1285  !
1286  ! -- Use ubdsv06 to write list header
1287  call ubdsv06(kstp, kper, text, textmodel, textpackage, dstmodel, dstpackage, &
1288  ibdchn, naux, auxtxt, ncol, nrow, nlay, &
1289  nlist, iout, delt, pertim, totim)
1290  !
1291  end subroutine record_srcdst_list_header
1292 
1293  !> @brief Convert an integer array (layer numbers) to nodelist
1294  !<
1295  subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
1296  ! -- dummy
1297  class(dis2dtype) :: this
1298  integer(I4B), intent(in) :: maxbnd
1299  integer(I4B), dimension(:), pointer, contiguous :: darray
1300  integer(I4B), dimension(maxbnd), intent(inout) :: nodelist
1301  integer(I4B), intent(inout) :: nbound
1302  character(len=*), intent(in) :: aname
1303  ! -- local
1304  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nval, nodeu, noder, ipos, ierr
1305  !
1306  ! -- set variables
1307  nlay = 1
1308  nrow = this%mshape(1)
1309  ncol = this%mshape(2)
1310  !
1311  if (this%ndim > 1) then
1312  !
1313  nval = ncol * nrow
1314  !
1315  ! -- Copy array into nodelist
1316  ipos = 1
1317  ierr = 0
1318  do ir = 1, nrow
1319  do ic = 1, ncol
1320  nodeu = get_node(1, ir, ic, nlay, nrow, ncol)
1321  il = darray(nodeu)
1322  if (il < 1 .or. il > nlay) then
1323  write (errmsg, '(a,1x,i0)') 'Invalid layer number:', il
1324  call store_error(errmsg, terminate=.true.)
1325  end if
1326  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
1327  noder = this%get_nodenumber(nodeu, 0)
1328  if (ipos > maxbnd) then
1329  ierr = ipos
1330  else
1331  nodelist(ipos) = noder
1332  end if
1333  ipos = ipos + 1
1334  end do
1335  end do
1336  !
1337  ! -- Check for errors
1338  nbound = ipos - 1
1339  if (ierr > 0) then
1340  write (errmsg, '(a,1x,i0)') &
1341  'MAXBOUND dimension is too small.'// &
1342  'INCREASE MAXBOUND TO:', ierr
1343  call store_error(errmsg, terminate=.true.)
1344  end if
1345  !
1346  ! -- If nbound < maxbnd, then initialize nodelist to zero in this range
1347  if (nbound < maxbnd) then
1348  do ipos = nbound + 1, maxbnd
1349  nodelist(ipos) = 0
1350  end do
1351  end if
1352  !
1353  else
1354  !
1355  ! -- For unstructured, read nodelist directly, then check node numbers
1356  nodelist = darray
1357  do noder = 1, maxbnd
1358  if (noder < 1 .or. noder > this%nodes) then
1359  write (errmsg, '(a,1x,i0)') 'Invalid node number:', noder
1360  call store_error(errmsg, terminate=.true.)
1361  end if
1362  end do
1363  nbound = maxbnd
1364  !
1365  end if
1366  !
1367  end subroutine nlarray_to_nodelist
1368 
1369 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
@ 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:1158
subroutine read_int_array(this, line, lloc, istart, istop, iout, in, iarray, aname)
Read an integer array.
Definition: Dis2d.f90:1119
integer(i4b) function get_nodenumber_idx2(this, k, j, icheck)
Get reduced node number from layer, row and column indices.
Definition: Dis2d.f90:669
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:590
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:1267
subroutine get_dis_type(this, dis_type)
Get the discretization type.
Definition: Dis2d.f90:1050
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:849
integer(i4b) function get_dis_enum(this)
Get the discretization type enumeration.
Definition: Dis2d.f90:1060
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:609
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:1001
subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, darray, aname)
Read a double precision array.
Definition: Dis2d.f90:1138
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:937
subroutine record_array(this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
Record a double precision array.
Definition: Dis2d.f90:1178
subroutine allocate_scalars(this, name_model, input_mempath)
Allocate and initialize scalar variables.
Definition: Dis2d.f90:712
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:927
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:734
integer(i4b) function get_nodenumber_idx1(this, nodeu, icheck)
Get reduced node number from user node number.
Definition: Dis2d.f90:638
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:948
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:764
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:1296
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:1071
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