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