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