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