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