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