MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 memorylist_remove(this%name_model, 'DISV', idm_context)
175  call memorylist_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  ! -- dummy
713  class(disvtype) :: this
714  integer(I4B), dimension(:), intent(in) :: icelltype
715  ! -- local
716  integer(I4B) :: iunit, i, ntxt
717  integer(I4B), parameter :: lentxt = 100
718  character(len=50) :: txthdr
719  character(len=lentxt) :: txt
720  character(len=LINELENGTH) :: fname
721  ! -- formats
722  character(len=*), parameter :: fmtgrdsave = &
723  "(4X,'BINARY GRID INFORMATION WILL BE WRITTEN TO:', &
724  &/,6X,'UNIT NUMBER: ', I0,/,6X, 'FILE NAME: ', A)"
725  !
726  ! -- Initialize
727  ntxt = 20
728  !
729  ! -- Open the file
730  fname = trim(this%input_fname)//'.grb'
731  iunit = getunit()
732  write (this%iout, fmtgrdsave) iunit, trim(adjustl(fname))
733  call openfile(iunit, this%iout, trim(adjustl(fname)), 'DATA(BINARY)', &
734  form, access, 'REPLACE')
735  !
736  ! -- write header information
737  write (txthdr, '(a)') 'GRID DISV'
738  txthdr(50:50) = new_line('a')
739  write (iunit) txthdr
740  write (txthdr, '(a)') 'VERSION 1'
741  txthdr(50:50) = new_line('a')
742  write (iunit) txthdr
743  write (txthdr, '(a, i0)') 'NTXT ', ntxt
744  txthdr(50:50) = new_line('a')
745  write (iunit) txthdr
746  write (txthdr, '(a, i0)') 'LENTXT ', lentxt
747  txthdr(50:50) = new_line('a')
748  write (iunit) txthdr
749  !
750  ! -- write variable definitions
751  write (txt, '(3a, i0)') 'NCELLS ', 'INTEGER ', 'NDIM 0 # ', this%nodesuser
752  txt(lentxt:lentxt) = new_line('a')
753  write (iunit) txt
754  write (txt, '(3a, i0)') 'NLAY ', 'INTEGER ', 'NDIM 0 # ', this%nlay
755  txt(lentxt:lentxt) = new_line('a')
756  write (iunit) txt
757  write (txt, '(3a, i0)') 'NCPL ', 'INTEGER ', 'NDIM 0 # ', this%ncpl
758  txt(lentxt:lentxt) = new_line('a')
759  write (iunit) txt
760  write (txt, '(3a, i0)') 'NVERT ', 'INTEGER ', 'NDIM 0 # ', this%nvert
761  txt(lentxt:lentxt) = new_line('a')
762  write (iunit) txt
763  write (txt, '(3a, i0)') 'NJAVERT ', 'INTEGER ', 'NDIM 0 # ', size(this%javert)
764  txt(lentxt:lentxt) = new_line('a')
765  write (iunit) txt
766  write (txt, '(3a, i0)') 'NJA ', 'INTEGER ', 'NDIM 0 # ', this%con%nja
767  txt(lentxt:lentxt) = new_line('a')
768  write (iunit) txt
769  write (txt, '(3a, 1pg25.15e3)') &
770  'XORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%xorigin
771  txt(lentxt:lentxt) = new_line('a')
772  write (iunit) txt
773  write (txt, '(3a, 1pg25.15e3)') &
774  'YORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%yorigin
775  txt(lentxt:lentxt) = new_line('a')
776  write (iunit) txt
777  write (txt, '(3a, 1pg25.15e3)') 'ANGROT ', 'DOUBLE ', 'NDIM 0 # ', this%angrot
778  txt(lentxt:lentxt) = new_line('a')
779  write (iunit) txt
780  write (txt, '(3a, i0)') 'TOP ', 'DOUBLE ', 'NDIM 1 ', this%ncpl
781  txt(lentxt:lentxt) = new_line('a')
782  write (iunit) txt
783  write (txt, '(3a, i0)') 'BOTM ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
784  txt(lentxt:lentxt) = new_line('a')
785  write (iunit) txt
786  write (txt, '(3a, i0)') 'VERTICES ', 'DOUBLE ', 'NDIM 2 2 ', this%nvert
787  txt(lentxt:lentxt) = new_line('a')
788  write (iunit) txt
789  write (txt, '(3a, i0)') 'CELLX ', 'DOUBLE ', 'NDIM 1 ', this%ncpl
790  txt(lentxt:lentxt) = new_line('a')
791  write (iunit) txt
792  write (txt, '(3a, i0)') 'CELLY ', 'DOUBLE ', 'NDIM 1 ', this%ncpl
793  txt(lentxt:lentxt) = new_line('a')
794  write (iunit) txt
795  write (txt, '(3a, i0)') 'IAVERT ', 'INTEGER ', 'NDIM 1 ', this%ncpl + 1
796  txt(lentxt:lentxt) = new_line('a')
797  write (iunit) txt
798  write (txt, '(3a, i0)') 'JAVERT ', 'INTEGER ', 'NDIM 1 ', size(this%javert)
799  txt(lentxt:lentxt) = new_line('a')
800  write (iunit) txt
801  write (txt, '(3a, i0)') 'IA ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
802  txt(lentxt:lentxt) = new_line('a')
803  write (iunit) txt
804  write (txt, '(3a, i0)') 'JA ', 'INTEGER ', 'NDIM 1 ', size(this%con%jausr)
805  txt(lentxt:lentxt) = new_line('a')
806  write (iunit) txt
807  write (txt, '(3a, i0)') 'IDOMAIN ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
808  txt(lentxt:lentxt) = new_line('a')
809  write (iunit) txt
810  write (txt, '(3a, i0)') 'ICELLTYPE ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
811  txt(lentxt:lentxt) = new_line('a')
812  write (iunit) txt
813  !
814  ! -- write data
815  write (iunit) this%nodesuser ! ncells
816  write (iunit) this%nlay ! nlay
817  write (iunit) this%ncpl ! ncpl
818  write (iunit) this%nvert ! nvert
819  write (iunit) size(this%javert) ! njavert
820  write (iunit) this%nja ! nja
821  write (iunit) this%xorigin ! xorigin
822  write (iunit) this%yorigin ! yorigin
823  write (iunit) this%angrot ! angrot
824  write (iunit) this%top1d ! top
825  write (iunit) this%bot2d ! botm
826  write (iunit) this%vertices ! vertices
827  write (iunit) (this%cellxy(1, i), i=1, this%ncpl) ! cellx
828  write (iunit) (this%cellxy(2, i), i=1, this%ncpl) ! celly
829  write (iunit) this%iavert ! iavert
830  write (iunit) this%javert ! javert
831  write (iunit) this%con%iausr ! iausr
832  write (iunit) this%con%jausr ! jausr
833  write (iunit) this%idomain ! idomain
834  write (iunit) icelltype ! icelltype
835  !
836  ! -- Close the file
837  close (iunit)
838  !
839  end subroutine write_grb
840 
841  !> @brief Convert a user nodenumber to a string (nodenumber) or (k,j)
842  !<
843  subroutine nodeu_to_string(this, nodeu, str)
844  ! -- dummy
845  class(disvtype) :: this
846  integer(I4B), intent(in) :: nodeu
847  character(len=*), intent(inout) :: str
848  ! -- local
849  integer(I4B) :: i, j, k
850  character(len=10) :: kstr, jstr
851  !
852  call get_ijk(nodeu, 1, this%ncpl, this%nlay, i, j, k)
853  write (kstr, '(i10)') k
854  write (jstr, '(i10)') j
855  str = '('//trim(adjustl(kstr))//','// &
856  trim(adjustl(jstr))//')'
857  !
858  end subroutine nodeu_to_string
859 
860  !> @brief Convert a user nodenumber to an array (nodenumber) or (k,j)
861  !<
862  subroutine nodeu_to_array(this, nodeu, arr)
863  ! -- dummy
864  class(disvtype) :: this
865  integer(I4B), intent(in) :: nodeu
866  integer(I4B), dimension(:), intent(inout) :: arr
867  ! -- local
868  integer(I4B) :: isize
869  integer(I4B) :: i, j, k
870  !
871  ! -- check the size of arr
872  isize = size(arr)
873  if (isize /= this%ndim) then
874  write (errmsg, '(a,i0,a,i0,a)') &
875  'Program error: nodeu_to_array size of array (', isize, &
876  ') is not equal to the discretization dimension (', this%ndim, ').'
877  call store_error(errmsg, terminate=.true.)
878  end if
879  !
880  ! -- get k, i, j
881  call get_ijk(nodeu, 1, this%ncpl, this%nlay, i, j, k)
882  !
883  ! -- fill array
884  arr(1) = k
885  arr(2) = j
886  !
887  end subroutine nodeu_to_array
888 
889  !> @brief Get reduced node number from user node number
890  !<
891  function get_nodenumber_idx1(this, nodeu, icheck) result(nodenumber)
892  ! -- return
893  integer(I4B) :: nodenumber
894  ! -- dummy
895  class(disvtype), intent(in) :: this
896  integer(I4B), intent(in) :: nodeu
897  integer(I4B), intent(in) :: icheck
898  ! -- local
899  !
900  ! -- check the node number if requested
901  if (icheck /= 0) then
902  !
903  ! -- If within valid range, convert to reduced nodenumber
904  if (nodeu < 1 .or. nodeu > this%nodesuser) then
905  nodenumber = 0
906  write (errmsg, '(a,i0,a,i0,a)') &
907  'Node number (', nodeu, ') is less than 1 or greater than nodes (', &
908  this%nodesuser, ').'
909  call store_error(errmsg)
910  else
911  nodenumber = nodeu
912  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
913  end if
914  else
915  nodenumber = nodeu
916  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
917  end if
918  !
919  end function get_nodenumber_idx1
920 
921  !> @brief Get reduced node number from layer and within-layer node indices
922  !<
923  function get_nodenumber_idx2(this, k, j, icheck) result(nodenumber)
924  ! -- return
925  integer(I4B) :: nodenumber
926  ! -- dummy
927  class(disvtype), intent(in) :: this
928  integer(I4B), intent(in) :: k, j
929  integer(I4B), intent(in) :: icheck
930  ! -- local
931  integer(I4B) :: nodeu
932  ! -- formats
933  character(len=*), parameter :: fmterr = &
934  &"('Error in disv grid cell indices: layer = ',i0,', node = ',i0)"
935  !
936  nodeu = get_node(k, 1, j, this%nlay, 1, this%ncpl)
937  if (nodeu < 1) then
938  write (errmsg, fmterr) k, j
939  call store_error(errmsg, terminate=.true.)
940  end if
941  nodenumber = nodeu
942  if (this%nodes < this%nodesuser) nodenumber = this%nodereduced(nodeu)
943  !
944  ! -- check the node number if requested
945  if (icheck /= 0) then
946  !
947  errmsg = ""
948  !
949  if (k < 1 .or. k > this%nlay) then
950  write (errmsg, '(a,i0,a)') &
951  'Layer number in list (', k, ') is outside of the grid.'
952  end if
953  if (j < 1 .or. j > this%ncpl) then
954  write (errmsg, '(a,1x,a,i0,a)') &
955  trim(adjustl(errmsg)), 'Node number in list (', j, &
956  ') is outside of the grid.'
957  end if
958  !
959  ! -- Error if outside of range
960  if (nodeu < 1 .or. nodeu > this%nodesuser) then
961  write (errmsg, '(a,1x,a,i0,a,i0,a)') &
962  trim(adjustl(errmsg)), &
963  'Node number (', nodeu, ') is less than 1 or greater '// &
964  'than nodes (', this%nodesuser, ').'
965  end if
966  !
967  if (len_trim(adjustl(errmsg)) > 0) then
968  call store_error(errmsg)
969  end if
970  !
971  end if
972  !
973  end function get_nodenumber_idx2
974 
975  !> @brief Get normal vector components between the cell and a given neighbor
976  !!
977  !! The normal points outward from the shared face between noden and nodem.
978  !<
979  subroutine connection_normal(this, noden, nodem, ihc, xcomp, ycomp, zcomp, &
980  ipos)
981  ! -- dummy
982  class(disvtype) :: this
983  integer(I4B), intent(in) :: noden !< cell (reduced nn)
984  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
985  integer(I4B), intent(in) :: ihc !< horizontal connection flag
986  real(DP), intent(inout) :: xcomp
987  real(DP), intent(inout) :: ycomp
988  real(DP), intent(inout) :: zcomp
989  integer(I4B), intent(in) :: ipos
990  ! -- local
991  real(DP) :: angle, dmult
992  !
993  ! -- Set vector components based on ihc
994  if (ihc == 0) then
995  xcomp = dzero
996  ycomp = dzero
997  if (nodem < noden) then
998  !
999  ! -- nodem must be above noden, so upward connection
1000  zcomp = done
1001  else
1002  !
1003  ! -- nodem must be below noden, so downward connection
1004  zcomp = -done
1005  end if
1006  else
1007  ! -- find from anglex, since anglex is symmetric, need to flip vector
1008  ! for lower triangle (nodem < noden)
1009  !ipos = this%con%getjaindex(noden, nodem)
1010  angle = this%con%anglex(this%con%jas(ipos))
1011  dmult = done
1012  if (nodem < noden) dmult = -done
1013  xcomp = cos(angle) * dmult
1014  ycomp = sin(angle) * dmult
1015  zcomp = dzero
1016  end if
1017  !
1018  end subroutine connection_normal
1019 
1020  !> @brief Get unit vector components between the cell and a given neighbor
1021  !!
1022  !! Saturation must be provided to compute cell center vertical coordinates.
1023  !! Also return the straight-line connection length.
1024  !<
1025  subroutine connection_vector(this, noden, nodem, nozee, satn, satm, ihc, &
1026  xcomp, ycomp, zcomp, conlen)
1027  ! -- dummy
1028  class(disvtype) :: this
1029  integer(I4B), intent(in) :: noden !< cell (reduced nn)
1030  integer(I4B), intent(in) :: nodem !< neighbor (reduced nn)
1031  logical, intent(in) :: nozee
1032  real(DP), intent(in) :: satn
1033  real(DP), intent(in) :: satm
1034  integer(I4B), intent(in) :: ihc !< horizontal connection flag
1035  real(DP), intent(inout) :: xcomp
1036  real(DP), intent(inout) :: ycomp
1037  real(DP), intent(inout) :: zcomp
1038  real(DP), intent(inout) :: conlen
1039  ! -- local
1040  integer(I4B) :: nodeu, ncell2d, mcell2d, k
1041  real(DP) :: xn, xm, yn, ym, zn, zm
1042  !
1043  ! -- Set vector components based on ihc
1044  if (ihc == 0) then
1045  !
1046  ! -- vertical connection; set zcomp positive upward
1047  xcomp = dzero
1048  ycomp = dzero
1049  if (nodem < noden) then
1050  zcomp = done
1051  else
1052  zcomp = -done
1053  end if
1054  zn = this%bot(noden) + dhalf * (this%top(noden) - this%bot(noden))
1055  zm = this%bot(nodem) + dhalf * (this%top(nodem) - this%bot(nodem))
1056  conlen = abs(zm - zn)
1057  else
1058  !
1059  ! -- horizontal connection, with possible z component due to cell offsets
1060  ! and/or water table conditions
1061  if (nozee) then
1062  zn = dzero
1063  zm = dzero
1064  else
1065  zn = this%bot(noden) + dhalf * satn * (this%top(noden) - this%bot(noden))
1066  zm = this%bot(nodem) + dhalf * satm * (this%top(nodem) - this%bot(nodem))
1067  end if
1068  nodeu = this%get_nodeuser(noden)
1069  call get_jk(nodeu, this%ncpl, this%nlay, ncell2d, k)
1070  nodeu = this%get_nodeuser(nodem)
1071  call get_jk(nodeu, this%ncpl, this%nlay, mcell2d, k)
1072  xn = this%cellxy(1, ncell2d)
1073  yn = this%cellxy(2, ncell2d)
1074  xm = this%cellxy(1, mcell2d)
1075  ym = this%cellxy(2, mcell2d)
1076  call line_unit_vector(xn, yn, zn, xm, ym, zm, xcomp, ycomp, zcomp, &
1077  conlen)
1078  end if
1079  !
1080  end subroutine connection_vector
1081 
1082  !> @brief Get the discretization type
1083  !<
1084  subroutine get_dis_type(this, dis_type)
1085  ! -- dummy
1086  class(disvtype), intent(in) :: this
1087  character(len=*), intent(out) :: dis_type
1088  !
1089  dis_type = "DISV"
1090  !
1091  end subroutine get_dis_type
1092 
1093  !> @brief Get the discretization type enumeration
1094  function get_dis_enum(this) result(dis_enum)
1095  use constantsmodule, only: disv
1096  class(disvtype), intent(in) :: this
1097  integer(I4B) :: dis_enum
1098  dis_enum = disv
1099  end function get_dis_enum
1100 
1101  !> @brief Allocate and initialize scalars
1102  !<
1103  subroutine allocate_scalars(this, name_model, input_mempath)
1104  ! -- dummy
1105  class(disvtype) :: this
1106  character(len=*), intent(in) :: name_model
1107  character(len=*), intent(in) :: input_mempath
1108  !
1109  ! -- Allocate parent scalars
1110  call this%DisBaseType%allocate_scalars(name_model, input_mempath)
1111  !
1112  ! -- Allocate
1113  call mem_allocate(this%nlay, 'NLAY', this%memoryPath)
1114  call mem_allocate(this%ncpl, 'NCPL', this%memoryPath)
1115  call mem_allocate(this%nvert, 'NVERT', this%memoryPath)
1116  !
1117  ! -- Initialize
1118  this%nlay = 0
1119  this%ncpl = 0
1120  this%nvert = 0
1121  this%ndim = 2
1122  !
1123  end subroutine allocate_scalars
1124 
1125  !> @brief Allocate and initialize arrays
1126  !<
1127  subroutine allocate_arrays(this)
1128  ! -- dummy
1129  class(disvtype) :: this
1130  !
1131  ! -- Allocate arrays in DisBaseType (mshape, top, bot, area)
1132  call this%DisBaseType%allocate_arrays()
1133  !
1134  ! -- Allocate arrays for DisvType
1135  if (this%nodes < this%nodesuser) then
1136  call mem_allocate(this%nodeuser, this%nodes, 'NODEUSER', this%memoryPath)
1137  call mem_allocate(this%nodereduced, this%nodesuser, 'NODEREDUCED', &
1138  this%memoryPath)
1139  else
1140  call mem_allocate(this%nodeuser, 1, 'NODEUSER', this%memoryPath)
1141  call mem_allocate(this%nodereduced, 1, 'NODEREDUCED', this%memoryPath)
1142  end if
1143  ! -- Initialize
1144  this%mshape(1) = this%nlay
1145  this%mshape(2) = this%ncpl
1146  !
1147  end subroutine allocate_arrays
1148 
1149  !> @brief Get the signed area of the cell
1150  !!
1151  !! A negative result means points are in counter-clockwise orientation.
1152  !! Area is computed from the formula:
1153  !! a = 1/2 *[(x1*y2 + x2*y3 + x3*y4 + ... + xn*y1) -
1154  !! (x2*y1 + x3*y2 + x4*y3 + ... + x1*yn)]
1155  !<
1156  function get_cell2d_area(this, icell2d) result(area)
1157  ! -- dummy
1158  class(disvtype) :: this
1159  integer(I4B), intent(in) :: icell2d
1160  ! -- return
1161  real(dp) :: area
1162  ! -- local
1163  integer(I4B) :: ivert
1164  integer(I4B) :: nvert
1165  integer(I4B) :: icount
1166  integer(I4B) :: iv1
1167  real(dp) :: x
1168  real(dp) :: y
1169  real(dp) :: x1
1170  real(dp) :: y1
1171  !
1172  area = dzero
1173  nvert = this%iavert(icell2d + 1) - this%iavert(icell2d)
1174  icount = 1
1175  iv1 = this%javert(this%iavert(icell2d))
1176  x1 = this%vertices(1, iv1)
1177  y1 = this%vertices(2, iv1)
1178  do ivert = this%iavert(icell2d), this%iavert(icell2d + 1) - 1
1179  x = this%vertices(1, this%javert(ivert))
1180  if (icount < nvert) then
1181  y = this%vertices(2, this%javert(ivert + 1))
1182  else
1183  y = this%vertices(2, this%javert(this%iavert(icell2d)))
1184  end if
1185  area = area + (x - x1) * (y - y1)
1186  icount = icount + 1
1187  end do
1188  !
1189  icount = 1
1190  do ivert = this%iavert(icell2d), this%iavert(icell2d + 1) - 1
1191  y = this%vertices(2, this%javert(ivert))
1192  if (icount < nvert) then
1193  x = this%vertices(1, this%javert(ivert + 1))
1194  else
1195  x = this%vertices(1, this%javert(this%iavert(icell2d)))
1196  end if
1197  area = area - (x - x1) * (y - y1)
1198  icount = icount + 1
1199  end do
1200  !
1201  area = -done * area * dhalf
1202  !
1203  end function get_cell2d_area
1204 
1205  !> @brief Convert a string to a user nodenumber
1206  !!
1207  !! Parse layer and within-layer cell number and return user nodenumber.
1208  !! If flag_string is present and true, the first token may be
1209  !! non-numeric (e.g. boundary name). In this case, return -2.
1210  !<
1211  function nodeu_from_string(this, lloc, istart, istop, in, iout, line, &
1212  flag_string, allow_zero) result(nodeu)
1213  ! -- dummy
1214  class(disvtype) :: this
1215  integer(I4B), intent(inout) :: lloc
1216  integer(I4B), intent(inout) :: istart
1217  integer(I4B), intent(inout) :: istop
1218  integer(I4B), intent(in) :: in
1219  integer(I4B), intent(in) :: iout
1220  character(len=*), intent(inout) :: line
1221  logical, optional, intent(in) :: flag_string
1222  logical, optional, intent(in) :: allow_zero
1223  integer(I4B) :: nodeu
1224  ! -- local
1225  integer(I4B) :: j, k, nlay, nrow, ncpl
1226  integer(I4B) :: lloclocal, ndum, istat, n
1227  real(dp) :: r
1228  !
1229  if (present(flag_string)) then
1230  if (flag_string) then
1231  ! Check to see if first token in line can be read as an integer.
1232  lloclocal = lloc
1233  call urword(line, lloclocal, istart, istop, 1, ndum, r, iout, in)
1234  read (line(istart:istop), *, iostat=istat) n
1235  if (istat /= 0) then
1236  ! First token in line is not an integer; return flag to this effect.
1237  nodeu = -2
1238  return
1239  end if
1240  end if
1241  end if
1242  !
1243  nlay = this%mshape(1)
1244  nrow = 1
1245  ncpl = this%mshape(2)
1246  !
1247  call urword(line, lloc, istart, istop, 2, k, r, iout, in)
1248  call urword(line, lloc, istart, istop, 2, j, r, iout, in)
1249  !
1250  if (k == 0 .and. j == 0) then
1251  if (present(allow_zero)) then
1252  if (allow_zero) then
1253  nodeu = 0
1254  return
1255  end if
1256  end if
1257  end if
1258  !
1259  errmsg = ''
1260  !
1261  if (k < 1 .or. k > nlay) then
1262  write (errmsg, '(a,i0,a)') &
1263  'Layer number in list (', k, ') is outside of the grid.'
1264  end if
1265  if (j < 1 .or. j > ncpl) then
1266  write (errmsg, '(a,1x,a,i0,a)') &
1267  trim(adjustl(errmsg)), 'Cell2d number in list (', j, &
1268  ') is outside of the grid.'
1269  end if
1270  !
1271  nodeu = get_node(k, 1, j, nlay, nrow, ncpl)
1272  !
1273  if (nodeu < 1 .or. nodeu > this%nodesuser) then
1274  write (errmsg, '(a,1x,a,i0,a)') &
1275  trim(adjustl(errmsg)), &
1276  "Node number in list (", nodeu, ") is outside of the grid. "// &
1277  "Cell number cannot be determined in line '"// &
1278  trim(adjustl(line))//"'."
1279  end if
1280  !
1281  if (len_trim(adjustl(errmsg)) > 0) then
1282  call store_error(errmsg)
1283  call store_error_unit(in)
1284  end if
1285  !
1286  end function nodeu_from_string
1287 
1288  !> @brief Convert a cellid string to a user nodenumber
1289  !!
1290  !! If flag_string is present and true, the first token may be
1291  !! non-numeric (e.g. boundary name). In this case, return -2.
1292  !!
1293  !! If allow_zero is present and true, and all indices are zero, the
1294  !! result can be zero. If allow_zero is false, a zero in any index is an error.
1295  !<
1296  function nodeu_from_cellid(this, cellid, inunit, iout, flag_string, &
1297  allow_zero) result(nodeu)
1298  ! -- return
1299  integer(I4B) :: nodeu
1300  ! -- dummy
1301  class(disvtype) :: this
1302  character(len=*), intent(inout) :: cellid
1303  integer(I4B), intent(in) :: inunit
1304  integer(I4B), intent(in) :: iout
1305  logical, optional, intent(in) :: flag_string
1306  logical, optional, intent(in) :: allow_zero
1307  ! -- local
1308  integer(I4B) :: j, k, nlay, nrow, ncpl
1309  integer(I4B) :: lloclocal, ndum, istat, n
1310  integer(I4B) :: istart, istop
1311  real(dp) :: r
1312  !
1313  if (present(flag_string)) then
1314  if (flag_string) then
1315  ! Check to see if first token in cellid can be read as an integer.
1316  lloclocal = 1
1317  call urword(cellid, lloclocal, istart, istop, 1, ndum, r, iout, inunit)
1318  read (cellid(istart:istop), *, iostat=istat) n
1319  if (istat /= 0) then
1320  ! First token in cellid is not an integer; return flag to this effect.
1321  nodeu = -2
1322  return
1323  end if
1324  end if
1325  end if
1326  !
1327  nlay = this%mshape(1)
1328  nrow = 1
1329  ncpl = this%mshape(2)
1330  !
1331  lloclocal = 1
1332  call urword(cellid, lloclocal, istart, istop, 2, k, r, iout, inunit)
1333  call urword(cellid, lloclocal, istart, istop, 2, j, r, iout, inunit)
1334  !
1335  if (k == 0 .and. j == 0) then
1336  if (present(allow_zero)) then
1337  if (allow_zero) then
1338  nodeu = 0
1339  return
1340  end if
1341  end if
1342  end if
1343  !
1344  errmsg = ''
1345  !
1346  if (k < 1 .or. k > nlay) then
1347  write (errmsg, '(a,i0,a)') &
1348  'Layer number in list (', k, ') is outside of the grid.'
1349  end if
1350  if (j < 1 .or. j > ncpl) then
1351  write (errmsg, '(a,1x,a,i0,a)') &
1352  trim(adjustl(errmsg)), 'Cell2d number in list (', j, &
1353  ') is outside of the grid.'
1354  end if
1355  !
1356  nodeu = get_node(k, 1, j, nlay, nrow, ncpl)
1357  !
1358  if (nodeu < 1 .or. nodeu > this%nodesuser) then
1359  write (errmsg, '(a,1x,a,i0,a)') &
1360  trim(adjustl(errmsg)), &
1361  "Cell number cannot be determined for cellid ("// &
1362  trim(adjustl(cellid))//") and results in a user "// &
1363  "node number (", nodeu, ") that is outside of the grid."
1364  end if
1365  !
1366  if (len_trim(adjustl(errmsg)) > 0) then
1367  call store_error(errmsg)
1368  call store_error_unit(inunit)
1369  end if
1370  !
1371  end function nodeu_from_cellid
1372 
1373  !> @brief Indicates whether the grid discretization supports layers
1374  !<
1375  logical function supports_layers(this)
1376  ! -- dummy
1377  class(disvtype) :: this
1378  !
1379  supports_layers = .true.
1380  !
1381  end function supports_layers
1382 
1383  !> @brief Get number of cells per layer (ncpl)
1384  !<
1385  function get_ncpl(this)
1386  ! -- return
1387  integer(I4B) :: get_ncpl
1388  ! -- dummy
1389  class(disvtype) :: this
1390  !
1391  get_ncpl = this%ncpl
1392  !
1393  end function get_ncpl
1394 
1395  !> @brief Get a 2D array of polygon vertices, listed in clockwise order
1396  !! beginning with the lower left corner
1397  !<
1398  subroutine get_polyverts(this, ic, polyverts, closed)
1399  ! -- dummy
1400  class(disvtype), intent(inout) :: this
1401  integer(I4B), intent(in) :: ic !< cell number (reduced)
1402  real(DP), allocatable, intent(out) :: polyverts(:, :) !< polygon vertices (column-major indexing)
1403  logical(LGP), intent(in), optional :: closed !< whether to close the polygon, duplicating a vertex (default false)
1404  ! -- local
1405  integer(I4B) :: icu, icu2d, iavert, ncpl, nverts, m, j
1406  logical(LGP) :: lclosed
1407  !
1408  ! count vertices
1409  ncpl = this%get_ncpl()
1410  icu = this%get_nodeuser(ic)
1411  icu2d = icu - ((icu - 1) / ncpl) * ncpl
1412  nverts = this%iavert(icu2d + 1) - this%iavert(icu2d) - 1
1413  if (nverts .le. 0) nverts = nverts + size(this%javert)
1414  !
1415  ! check closed option
1416  if (.not. (present(closed))) then
1417  lclosed = .false.
1418  else
1419  lclosed = closed
1420  end if
1421  !
1422  ! allocate vertices array
1423  if (lclosed) then
1424  allocate (polyverts(2, nverts + 1))
1425  else
1426  allocate (polyverts(2, nverts))
1427  end if
1428  !
1429  ! set vertices
1430  iavert = this%iavert(icu2d)
1431  do m = 1, nverts
1432  j = this%javert(iavert - 1 + m)
1433  polyverts(:, m) = (/this%vertices(1, j), this%vertices(2, j)/)
1434  end do
1435  !
1436  ! close if enabled
1437  if (lclosed) &
1438  polyverts(:, nverts + 1) = polyverts(:, 1)
1439  !
1440  end subroutine
1441 
1442  !> @brief Read an integer array
1443  !<
1444  subroutine read_int_array(this, line, lloc, istart, istop, iout, in, &
1445  iarray, aname)
1446  ! -- dummy
1447  class(disvtype), intent(inout) :: this
1448  character(len=*), intent(inout) :: line
1449  integer(I4B), intent(inout) :: lloc
1450  integer(I4B), intent(inout) :: istart
1451  integer(I4B), intent(inout) :: istop
1452  integer(I4B), intent(in) :: in
1453  integer(I4B), intent(in) :: iout
1454  integer(I4B), dimension(:), pointer, contiguous, intent(inout) :: iarray
1455  character(len=*), intent(in) :: aname
1456  ! -- local
1457  integer(I4B) :: ival
1458  real(DP) :: rval
1459  integer(I4B) :: nlay
1460  integer(I4B) :: nrow
1461  integer(I4B) :: ncol
1462  integer(I4B) :: nval
1463  integer(I4B), dimension(:), pointer, contiguous :: itemp
1464  !
1465  ! -- Point the temporary pointer array, which is passed to the reading
1466  ! subroutine. The temporary array will point to ibuff if it is a
1467  ! reduced structured system, or to iarray if it is an unstructured
1468  ! model.
1469  nlay = this%mshape(1)
1470  nrow = 1
1471  ncol = this%mshape(2)
1472  !
1473  if (this%nodes < this%nodesuser) then
1474  nval = this%nodesuser
1475  itemp => this%ibuff
1476  else
1477  nval = this%nodes
1478  itemp => iarray
1479  end if
1480  !
1481  ! -- Read the array
1482  call urword(line, lloc, istart, istop, 1, ival, rval, iout, in)
1483  if (line(istart:istop) .EQ. 'LAYERED') then
1484  !
1485  ! -- Read layered input
1486  call readarray(in, itemp, aname, this%ndim, ncol, nrow, nlay, nval, &
1487  iout, 1, nlay)
1488  else
1489  !
1490  ! -- Read unstructured input
1491  call readarray(in, itemp, aname, this%ndim, nval, iout, 0)
1492  end if
1493  !
1494  ! -- If reduced model, then need to copy from itemp(=>ibuff) to iarray
1495  if (this%nodes < this%nodesuser) then
1496  call this%fill_grid_array(itemp, iarray)
1497  end if
1498  !
1499  end subroutine read_int_array
1500 
1501  !> @brief Read a double precision array
1502  !<
1503  subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, &
1504  darray, aname)
1505  ! -- dummy
1506  class(disvtype), intent(inout) :: this
1507  character(len=*), intent(inout) :: line
1508  integer(I4B), intent(inout) :: lloc
1509  integer(I4B), intent(inout) :: istart
1510  integer(I4B), intent(inout) :: istop
1511  integer(I4B), intent(in) :: in
1512  integer(I4B), intent(in) :: iout
1513  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray
1514  character(len=*), intent(in) :: aname
1515  ! -- local
1516  integer(I4B) :: ival
1517  real(DP) :: rval
1518  integer(I4B) :: nlay
1519  integer(I4B) :: nrow
1520  integer(I4B) :: ncol
1521  integer(I4B) :: nval
1522  real(DP), dimension(:), pointer, contiguous :: dtemp
1523  !
1524  ! -- Point the temporary pointer array, which is passed to the reading
1525  ! subroutine. The temporary array will point to dbuff if it is a
1526  ! reduced structured system, or to darray if it is an unstructured
1527  ! model.
1528  nlay = this%mshape(1)
1529  nrow = 1
1530  ncol = this%mshape(2)
1531  !
1532  if (this%nodes < this%nodesuser) then
1533  nval = this%nodesuser
1534  dtemp => this%dbuff
1535  else
1536  nval = this%nodes
1537  dtemp => darray
1538  end if
1539  !
1540  ! -- Read the array
1541  call urword(line, lloc, istart, istop, 1, ival, rval, iout, in)
1542  if (line(istart:istop) .EQ. 'LAYERED') then
1543  !
1544  ! -- Read structured input
1545  call readarray(in, dtemp, aname, this%ndim, ncol, nrow, nlay, nval, &
1546  iout, 1, nlay)
1547  else
1548  !
1549  ! -- Read unstructured input
1550  call readarray(in, dtemp, aname, this%ndim, nval, iout, 0)
1551  end if
1552  !
1553  ! -- If reduced model, then need to copy from dtemp(=>dbuff) to darray
1554  if (this%nodes < this%nodesuser) then
1555  call this%fill_grid_array(dtemp, darray)
1556  end if
1557  !
1558  end subroutine read_dbl_array
1559 
1560  !> @brief Read a 2d double array into col icolbnd of darray
1561  !!
1562  !! For cells that are outside of the active domain, do not copy the array
1563  !! value into darray.
1564  !<
1565  subroutine read_layer_array(this, nodelist, darray, ncolbnd, maxbnd, &
1566  icolbnd, aname, inunit, iout)
1567  ! -- dummy
1568  class(disvtype) :: this
1569  integer(I4B), intent(in) :: ncolbnd
1570  integer(I4B), intent(in) :: maxbnd
1571  integer(I4B), dimension(maxbnd) :: nodelist
1572  real(DP), dimension(ncolbnd, maxbnd), intent(inout) :: darray
1573  integer(I4B), intent(in) :: icolbnd
1574  character(len=*), intent(in) :: aname
1575  integer(I4B), intent(in) :: inunit
1576  integer(I4B), intent(in) :: iout
1577  ! -- local
1578  integer(I4B) :: ir, ic, ncol, nrow, nlay, nval, ipos, nodeu
1579  !
1580  ! -- set variables
1581  nlay = this%mshape(1)
1582  nrow = 1
1583  ncol = this%mshape(2)
1584  !
1585  ! -- Read the array
1586  nval = ncol * nrow
1587  call readarray(inunit, this%dbuff, aname, this%ndim, nval, iout, 0)
1588  !
1589  ! -- Copy array into bound. Note that this routine was substantially
1590  ! changed on 9/21/2021 to support changes to READASARRAYS input
1591  ! for recharge and evapotranspiration. nodelist and bound are of
1592  ! size nrow * ncol and correspond directly to dbuff.
1593  ipos = 1
1594  do ir = 1, nrow
1595  do ic = 1, ncol
1596  nodeu = get_node(1, ir, ic, nlay, nrow, ncol)
1597  darray(icolbnd, ipos) = this%dbuff(nodeu)
1598  ipos = ipos + 1
1599  end do
1600  end do
1601  !
1602  end subroutine read_layer_array
1603 
1604  !> @brief Record a double precision array
1605  !!
1606  !! The array is written to a formatted or unformatted external file depending
1607  !! on the arguments.
1608  !<
1609  subroutine record_array(this, darray, iout, iprint, idataun, aname, &
1610  cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
1611  ! -- dummy
1612  class(disvtype), intent(inout) :: this
1613  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray !< double precision array to record
1614  integer(I4B), intent(in) :: iout !< ascii output unit number
1615  integer(I4B), intent(in) :: iprint !< whether to print the array
1616  integer(I4B), intent(in) :: idataun !< binary output unit number, if negative don't write by layers, write entire array
1617  character(len=*), intent(in) :: aname !< text descriptor
1618  character(len=*), intent(in) :: cdatafmp !< write format
1619  integer(I4B), intent(in) :: nvaluesp !< values per line
1620  integer(I4B), intent(in) :: nwidthp !< number width
1621  character(len=*), intent(in) :: editdesc !< format type (I, G, F, S, E)
1622  real(DP), intent(in) :: dinact !< double precision value for cells excluded from model domain
1623  ! -- local
1624  integer(I4B) :: k, ifirst
1625  integer(I4B) :: nlay
1626  integer(I4B) :: nrow
1627  integer(I4B) :: ncol
1628  integer(I4B) :: nval
1629  integer(I4B) :: nodeu, noder
1630  integer(I4B) :: istart, istop
1631  real(DP), dimension(:), pointer, contiguous :: dtemp
1632  ! -- formats
1633  character(len=*), parameter :: fmthsv = &
1634  "(1X,/1X,a,' WILL BE SAVED ON UNIT ',I4, &
1635  &' AT END OF TIME STEP',I5,', STRESS PERIOD ',I4)"
1636  !
1637  ! -- set variables
1638  nlay = this%mshape(1)
1639  nrow = 1
1640  ncol = this%mshape(2)
1641  !
1642  ! -- If this is a reduced model, then copy the values from darray into
1643  ! dtemp.
1644  if (this%nodes < this%nodesuser) then
1645  nval = this%nodes
1646  dtemp => this%dbuff
1647  do nodeu = 1, this%nodesuser
1648  noder = this%get_nodenumber(nodeu, 0)
1649  if (noder <= 0) then
1650  dtemp(nodeu) = dinact
1651  cycle
1652  end if
1653  dtemp(nodeu) = darray(noder)
1654  end do
1655  else
1656  nval = this%nodes
1657  dtemp => darray
1658  end if
1659  !
1660  ! -- Print to iout if iprint /= 0
1661  if (iprint /= 0) then
1662  istart = 1
1663  do k = 1, nlay
1664  istop = istart + nrow * ncol - 1
1665  call ulaprufw(ncol, nrow, kstp, kper, k, iout, dtemp(istart:istop), &
1666  aname, cdatafmp, nvaluesp, nwidthp, editdesc)
1667  istart = istop + 1
1668  end do
1669  end if
1670  !
1671  ! -- Save array to an external file.
1672  if (idataun > 0) then
1673  ! -- write to binary file by layer
1674  ifirst = 1
1675  istart = 1
1676  do k = 1, nlay
1677  istop = istart + nrow * ncol - 1
1678  if (ifirst == 1) write (iout, fmthsv) &
1679  trim(adjustl(aname)), idataun, &
1680  kstp, kper
1681  ifirst = 0
1682  call ulasav(dtemp(istart:istop), aname, kstp, kper, &
1683  pertim, totim, ncol, nrow, k, idataun)
1684  istart = istop + 1
1685  end do
1686  elseif (idataun < 0) then
1687  !
1688  ! -- write entire array as one record
1689  call ubdsv1(kstp, kper, aname, -idataun, dtemp, ncol, nrow, nlay, &
1690  iout, delt, pertim, totim)
1691  end if
1692  !
1693  end subroutine record_array
1694 
1695  !> @brief Record list header for imeth=6
1696  !<
1697  subroutine record_srcdst_list_header(this, text, textmodel, textpackage, &
1698  dstmodel, dstpackage, naux, auxtxt, &
1699  ibdchn, nlist, iout)
1700  ! -- dummy
1701  class(disvtype) :: this
1702  character(len=16), intent(in) :: text
1703  character(len=16), intent(in) :: textmodel
1704  character(len=16), intent(in) :: textpackage
1705  character(len=16), intent(in) :: dstmodel
1706  character(len=16), intent(in) :: dstpackage
1707  integer(I4B), intent(in) :: naux
1708  character(len=16), dimension(:), intent(in) :: auxtxt
1709  integer(I4B), intent(in) :: ibdchn
1710  integer(I4B), intent(in) :: nlist
1711  integer(I4B), intent(in) :: iout
1712  ! -- local
1713  integer(I4B) :: nlay, nrow, ncol
1714  !
1715  nlay = this%mshape(1)
1716  nrow = 1
1717  ncol = this%mshape(2)
1718  !
1719  ! -- Use ubdsv06 to write list header
1720  call ubdsv06(kstp, kper, text, textmodel, textpackage, dstmodel, dstpackage, &
1721  ibdchn, naux, auxtxt, ncol, nrow, nlay, &
1722  nlist, iout, delt, pertim, totim)
1723  !
1724  end subroutine record_srcdst_list_header
1725 
1726  !> @brief Convert an integer array (layer numbers) to nodelist
1727  !<
1728  subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
1729  ! -- dummy
1730  class(disvtype) :: this
1731  integer(I4B), intent(in) :: maxbnd
1732  integer(I4B), dimension(:), pointer, contiguous :: darray
1733  integer(I4B), dimension(maxbnd), intent(inout) :: nodelist
1734  integer(I4B), intent(inout) :: nbound
1735  character(len=*), intent(in) :: aname
1736  ! -- local
1737  integer(I4B) :: il, ir, ic, ncol, nrow, nlay, nval, nodeu, noder, ipos, ierr
1738  !
1739  ! -- set variables
1740  nlay = this%mshape(1)
1741  nrow = 1
1742  ncol = this%mshape(2)
1743  !
1744  nval = ncol * nrow
1745  !
1746  ! -- Copy array into nodelist
1747  ipos = 1
1748  ierr = 0
1749  do ir = 1, nrow
1750  do ic = 1, ncol
1751  nodeu = get_node(1, ir, ic, nlay, nrow, ncol)
1752  il = darray(nodeu)
1753  if (il < 1 .or. il > nlay) then
1754  write (errmsg, '(a,i0,a)') &
1755  'Invalid layer number (', il, ').'
1756  call store_error(errmsg, terminate=.true.)
1757  end if
1758  nodeu = get_node(il, ir, ic, nlay, nrow, ncol)
1759  noder = this%get_nodenumber(nodeu, 0)
1760  if (ipos > maxbnd) then
1761  ierr = ipos
1762  else
1763  nodelist(ipos) = noder
1764  end if
1765  ipos = ipos + 1
1766  end do
1767  end do
1768  !
1769  ! -- Check for errors
1770  nbound = ipos - 1
1771  if (ierr > 0) then
1772  write (errmsg, '(a,i0,a)') &
1773  'MAXBOUND dimension is too small. Increase MAXBOUND to ', ierr, '.'
1774  call store_error(errmsg, terminate=.true.)
1775  end if
1776  !
1777  ! -- If nbound < maxbnd, then initialize nodelist to zero in this range
1778  if (nbound < maxbnd) then
1779  do ipos = nbound + 1, maxbnd
1780  nodelist(ipos) = 0
1781  end do
1782  end if
1783  !
1784  end subroutine nlarray_to_nodelist
1785 
1786 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:44
@ disv
DISU6 discretization.
Definition: Constants.f90:155
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:67
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
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:508
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:1157
integer(i4b) function get_ncpl(this)
Get number of cells per layer (ncpl)
Definition: Disv.f90:1386
subroutine record_array(this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
Record a double precision array.
Definition: Disv.f90:1611
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:1027
integer(i4b) function get_nodenumber_idx2(this, k, j, icheck)
Get reduced node number from layer and within-layer node indices.
Definition: Disv.f90:924
subroutine nodeu_to_string(this, nodeu, str)
Convert a user nodenumber to a string (nodenumber) or (k,j)
Definition: Disv.f90:844
subroutine nlarray_to_nodelist(this, darray, nodelist, maxbnd, nbound, aname)
Convert an integer array (layer numbers) to nodelist.
Definition: Disv.f90:1729
integer(i4b) function get_nodenumber_idx1(this, nodeu, icheck)
Get reduced node number from user node number.
Definition: Disv.f90:892
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:1399
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:1700
subroutine read_dbl_array(this, line, lloc, istart, istop, iout, in, darray, aname)
Read a double precision array.
Definition: Disv.f90:1505
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:863
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:1298
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:1446
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:1567
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:981
subroutine allocate_arrays(this)
Allocate and initialize arrays.
Definition: Disv.f90:1128
subroutine get_dis_type(this, dis_type)
Get the discretization type.
Definition: Disv.f90:1085
subroutine allocate_scalars(this, name_model, input_mempath)
Allocate and initialize scalars.
Definition: Disv.f90:1104
logical function supports_layers(this)
Indicates whether the grid discretization supports layers.
Definition: Disv.f90:1376
integer(i4b) function get_dis_enum(this)
Get the discretization type enumeration.
Definition: Disv.f90:1095
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:1213
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:80
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:97
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:125
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 memorylist_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