MODFLOW 6  version 6.7.0.dev2
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  end subroutine source_vertices
529 
530  !> @brief Copy cell1d information from input data context
531  !! to model context
532  !<
533  subroutine source_cell1d(this)
534  ! -- modules
537  ! -- dummy
538  class(disv1dtype) :: this
539  ! -- locals
540  character(len=LENMEMPATH) :: idmMemoryPath
541  integer(I4B), dimension(:), contiguous, pointer :: icell1d => null()
542  integer(I4B), dimension(:), contiguous, pointer :: ncvert => null()
543  integer(I4B), dimension(:), contiguous, pointer :: icvert => null()
544  real(DP), dimension(:), contiguous, pointer :: fdc => null()
545  integer(I4B) :: i
546  ! -- formats
547  !
548  ! -- set memory path
549  idmmemorypath = create_mem_path(this%name_model, 'DISV1D', idm_context)
550  !
551  ! -- set pointers to input path ncvert and icvert
552  call mem_setptr(icell1d, 'ICELL1D', idmmemorypath)
553  call mem_setptr(ncvert, 'NCVERT', idmmemorypath)
554  call mem_setptr(icvert, 'ICVERT', idmmemorypath)
555  !
556  ! --
557  if (associated(icell1d) .and. associated(ncvert) &
558  .and. associated(icvert)) then
559  call this%define_cellverts(icell1d, ncvert, icvert)
560  else
561  call store_error('Required cell vertex arrays not found.')
562  end if
563  !
564  ! -- set pointers to cell center arrays
565  call mem_setptr(fdc, 'FDC', idmmemorypath)
566  !
567  ! -- set fractional distance to cell center
568  if (associated(fdc)) then
569  do i = 1, this%nodesuser
570  this%fdc(i) = fdc(i)
571  end do
572  else
573  call store_error('Required fdc array not found.')
574  end if
575 
576  ! calculate length from vertices
577  call calculate_cell_length(this%vertices, this%iavert, this%javert, &
578  this%length)
579 
580  ! calculate cellxy from vertices and fdc
581  call calculate_cellxy(this%vertices, this%fdc, this%iavert, &
582  this%javert, this%length, this%cellxy)
583 
584  ! log
585  if (this%iout > 0) then
586  write (this%iout, '(1x,a)') 'Setting Discretization CELL1D'
587  write (this%iout, '(1x,a,/)') 'End Setting Discretization CELL1D'
588  end if
589  end subroutine source_cell1d
590 
591  !> @brief Construct the iavert and javert integer vectors which
592  !! are compressed sparse row index arrays that relate the vertices
593  !! to reaches
594  !<
595  subroutine define_cellverts(this, icell1d, ncvert, icvert)
596  ! -- modules
597  use sparsemodule, only: sparsematrix
598  ! -- dummy
599  class(disv1dtype) :: this
600  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: icell1d
601  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: ncvert
602  integer(I4B), dimension(:), contiguous, pointer, intent(in) :: icvert
603  ! -- locals
604  type(sparsematrix) :: vert_spm
605  integer(I4B) :: i, j, ierr
606  integer(I4B) :: icv_idx, startvert, maxnnz = 2
607  !
608  ! -- initialize sparse matrix
609  call vert_spm%init(this%nodesuser, this%nvert, maxnnz)
610  !
611  ! -- add sparse matrix connections from input memory paths
612  icv_idx = 1
613  do i = 1, this%nodesuser
614  if (icell1d(i) /= i) call store_error('ICELL1D input sequence violation.')
615  do j = 1, ncvert(i)
616  call vert_spm%addconnection(i, icvert(icv_idx), 0)
617  if (j == 1) then
618  startvert = icvert(icv_idx)
619  end if
620  icv_idx = icv_idx + 1
621  end do
622  end do
623  !
624  ! -- allocate and fill iavert and javert
625  call mem_allocate(this%iavert, this%nodesuser + 1, 'IAVERT', this%memoryPath)
626  call mem_allocate(this%javert, vert_spm%nnz, 'JAVERT', this%memoryPath)
627  call vert_spm%filliaja(this%iavert, this%javert, ierr)
628  call vert_spm%destroy()
629  end subroutine define_cellverts
630 
631  !> @brief Calculate x, y, coordinates of reach midpoint
632  !<
633  subroutine calculate_cellxy(vertices, fdc, iavert, javert, length, cellxy)
634  ! -- dummy
635  real(DP), dimension(:, :), intent(in) :: vertices !< 2d array of vertices with x, y as columns
636  real(DP), dimension(:), intent(in) :: fdc !< fractional distance to reach midpoint (normally 0.5)
637  integer(I4B), dimension(:), intent(in) :: iavert !< csr mapping of vertices to cell reaches
638  integer(I4B), dimension(:), intent(in) :: javert !< csr mapping of vertices to cell reaches
639  real(DP), dimension(:), intent(in) :: length !< vector of cell lengths
640  real(DP), dimension(:, :), intent(inout) :: cellxy !< 2d array of reach midpoint with x, y as columns
641  ! -- local
642  integer(I4B) :: nodes !< number of nodes
643  integer(I4B) :: n !< node index
644  integer(I4B) :: j !< vertex index
645  integer(I4B) :: iv0 !< index for line reach start
646  integer(I4B) :: iv1 !< index for linen reach end
647  integer(I4B) :: ixy !< x, y column index
648  real(DP) :: fd0 !< fractional distance to start of this line reach
649  real(DP) :: fd1 !< fractional distance to end of this line reach
650  real(DP) :: fd !< fractional distance where midpoint (defined by fdc) is located
651  real(DP) :: d !< distance
652 
653  nodes = size(iavert) - 1
654  do n = 1, nodes
655 
656  ! find vertices that span midpoint
657  iv0 = 0
658  iv1 = 0
659  fd0 = dzero
660  do j = iavert(n), iavert(n + 1) - 2
661  d = calcdist(vertices, javert(j), javert(j + 1))
662  fd1 = fd0 + d / length(n)
663 
664  ! if true, then we know the midpoint is some fractional distance (fd)
665  ! from vertex j to vertex j + 1
666  if (fd1 >= fdc(n)) then
667  iv0 = javert(j)
668  iv1 = javert(j + 1)
669  fd = (fdc(n) - fd0) / (fd1 - fd0)
670  exit
671  end if
672  fd0 = fd1
673  end do
674 
675  ! find x, y position of point on line
676  do ixy = 1, 2
677  cellxy(ixy, n) = (done - fd) * vertices(ixy, iv0) + &
678  fd * vertices(ixy, iv1)
679  end do
680 
681  end do
682  end subroutine calculate_cellxy
683 
684  !> @brief Calculate x, y, coordinates of reach midpoint
685  !<
686  subroutine calculate_cell_length(vertices, iavert, javert, length)
687  ! -- dummy
688  real(DP), dimension(:, :), intent(in) :: vertices !< 2d array of vertices with x, y as columns
689  integer(I4B), dimension(:), intent(in) :: iavert !< csr mapping of vertices to cell reaches
690  integer(I4B), dimension(:), intent(in) :: javert !< csr mapping of vertices to cell reaches
691  real(DP), dimension(:), intent(inout) :: length !< 2d array of reach midpoint with x, y as columns
692  ! -- local
693  integer(I4B) :: nodes !< number of nodes
694  integer(I4B) :: n !< node index
695  integer(I4B) :: j !< vertex index
696  real(DP) :: dlen !< length
697 
698  nodes = size(iavert) - 1
699  do n = 1, nodes
700 
701  ! calculate length of this reach
702  dlen = dzero
703  do j = iavert(n), iavert(n + 1) - 2
704  dlen = dlen + calcdist(vertices, javert(j), javert(j + 1))
705  end do
706  length(n) = dlen
707 
708  end do
709  end subroutine calculate_cell_length
710 
711  !> @brief Finalize grid construction
712  !<
713  subroutine grid_finalize(this)
714  ! modules
716  use constantsmodule, only: linelength, dzero, done
717  ! dummy
718  class(disv1dtype) :: this
719  ! local
720  integer(I4B) :: node, noder, k
721  ! format
722  character(len=*), parameter :: fmtnr = &
723  "(/1x, 'The specified IDOMAIN results in a reduced number of cells.',&
724  &/1x, 'Number of user nodes: ',I0,&
725  &/1X, 'Number of nodes in solution: ', I0, //)"
726 
727  ! count active cells
728  this%nodes = 0
729  do k = 1, this%nodesuser
730  if (this%idomain(k) > 0) this%nodes = this%nodes + 1
731  end do
732  !
733  ! Check to make sure nodes is a valid number
734  if (this%nodes == 0) then
735  call store_error('Model does not have any active nodes. Make sure &
736  &IDOMAIN has some values greater than zero.')
737  call store_error_filename(this%input_fname)
738  end if
739 
740  ! Write message if reduced grid
741  if (this%nodes < this%nodesuser) then
742  write (this%iout, fmtnr) this%nodesuser, this%nodes
743  end if
744 
745  ! Array size is now known, so allocate
746  call this%allocate_arrays()
747 
748  ! Fill the nodereduced array with the reduced nodenumber, or
749  ! a negative number to indicate it is a pass-through cell, or
750  ! a zero to indicate that the cell is excluded from the
751  ! solution.
752  if (this%nodes < this%nodesuser) then
753  node = 1
754  noder = 1
755  do k = 1, this%nodesuser
756  if (this%idomain(k) > 0) then
757  this%nodereduced(node) = noder
758  noder = noder + 1
759  else
760  this%nodereduced(node) = 0
761  end if
762  node = node + 1
763  end do
764  end if
765 
766  ! allocate and fill nodeuser if a reduced grid
767  if (this%nodes < this%nodesuser) then
768  node = 1
769  noder = 1
770  do k = 1, this%nodesuser
771  if (this%idomain(k) > 0) then
772  this%nodeuser(noder) = node
773  noder = noder + 1
774  end if
775  node = node + 1
776  end do
777  end if
778 
779  ! Move bottom into bot and put length into disbase%area
780  ! and set x and y center coordinates
781  do node = 1, this%nodesuser
782  noder = node
783  if (this%nodes < this%nodesuser) noder = this%nodereduced(node)
784  if (noder <= 0) cycle
785  this%bot(noder) = this%bottom(node)
786  this%area(noder) = this%length(node)
787  end do
788 
789  ! create connectivity using vertices and cell1d
790  call this%create_connections()
791  end subroutine grid_finalize
792 
793  subroutine allocate_arrays(this)
794  ! -- modules
796  ! -- dummy
797  class(disv1dtype) :: this
798  !
799  ! -- Allocate arrays in DisBaseType (mshape, top, bot, area)
800  ! todo: disbasetype will have memory allocated for unneeded arrays
801  call this%DisBaseType%allocate_arrays()
802  !
803  ! -- Allocate arrays
804  if (this%nodes < this%nodesuser) then
805  call mem_allocate(this%nodeuser, this%nodes, 'NODEUSER', this%memoryPath)
806  call mem_allocate(this%nodereduced, this%nodesuser, 'NODEREDUCED', &
807  this%memoryPath)
808  else
809  call mem_allocate(this%nodeuser, 1, 'NODEUSER', this%memoryPath)
810  call mem_allocate(this%nodereduced, 1, 'NODEREDUCED', this%memoryPath)
811  end if
812  !
813  ! -- Initialize
814  this%mshape(1) = this%nodesuser
815  end subroutine allocate_arrays
816 
817  subroutine create_connections(this)
818  ! modules
819  ! dummy
820  class(disv1dtype) :: this
821  ! local
822  integer(I4B) :: nrsize
823 
824  ! create and fill the connections object
825  nrsize = 0
826  if (this%nodes < this%nodesuser) nrsize = this%nodes
827 
828  ! Allocate connections object
829  allocate (this%con)
830 
831  ! Build connectivity based on vertices
832  call this%con%disv1dconnections_verts(this%name_model, this%nodes, &
833  this%nodesuser, nrsize, this%nvert, &
834  this%vertices, this%iavert, &
835  this%javert, this%cellxy, this%fdc, &
836  this%nodereduced, this%nodeuser, &
837  this%length)
838 
839  this%nja = this%con%nja
840  this%njas = this%con%njas
841 
842  end subroutine create_connections
843 
844  !> @brief Write binary grid file
845  !<
846  subroutine write_grb(this, icelltype)
847  ! -- modules
849  use openspecmodule, only: access, form
850  use constantsmodule, only: lenbigline
851  ! -- dummy
852  class(disv1dtype) :: this
853  integer(I4B), dimension(:), intent(in) :: icelltype
854  ! -- local
855  integer(I4B) :: i, iunit, ntxt, version
856  integer(I4B), parameter :: lentxt = 100
857  character(len=50) :: txthdr
858  character(len=lentxt) :: txt
859  character(len=LINELENGTH) :: fname
860  character(len=LENBIGLINE) :: crs
861  logical(LGP) :: found_crs
862  character(len=*), parameter :: fmtgrdsave = &
863  "(4X,'BINARY GRID INFORMATION WILL BE WRITTEN TO:', &
864  &/,6X,'UNIT NUMBER: ', I0,/,6X, 'FILE NAME: ', A)"
865  !
866  ! -- Initialize
867  version = 1
868  ntxt = 10
869  if (this%nvert > 0) ntxt = ntxt + 5
870  !
871  call mem_set_value(crs, 'CRS', this%input_mempath, found_crs)
872  !
873  ! -- set version
874  if (found_crs) then
875  ntxt = ntxt + 1
876  version = 2
877  end if
878  !
879  ! -- Open the file
880  fname = trim(this%output_fname)
881  iunit = getunit()
882  write (this%iout, fmtgrdsave) iunit, trim(adjustl(fname))
883  call openfile(iunit, this%iout, trim(adjustl(fname)), 'DATA(BINARY)', &
884  form, access, 'REPLACE')
885  !
886  ! -- write header information
887  write (txthdr, '(a)') 'GRID DISV1D'
888  txthdr(50:50) = new_line('a')
889  write (iunit) txthdr
890  write (txthdr, '(a)') 'VERSION 1'
891  txthdr(50:50) = new_line('a')
892  write (iunit) txthdr
893  write (txthdr, '(a, i0)') 'NTXT ', ntxt
894  txthdr(50:50) = new_line('a')
895  write (iunit) txthdr
896  write (txthdr, '(a, i0)') 'LENTXT ', lentxt
897  txthdr(50:50) = new_line('a')
898  write (iunit) txthdr
899  !
900  ! -- write variable definitions
901  write (txt, '(3a, i0)') 'NCELLS ', 'INTEGER ', 'NDIM 0 # ', this%nodesuser
902  txt(lentxt:lentxt) = new_line('a')
903  write (iunit) txt
904  write (txt, '(3a, i0)') 'NJA ', 'INTEGER ', 'NDIM 0 # ', this%con%nja
905  txt(lentxt:lentxt) = new_line('a')
906  write (iunit) txt
907  write (txt, '(3a, 1pg24.15)') 'XORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%xorigin
908  txt(lentxt:lentxt) = new_line('a')
909  write (iunit) txt
910  write (txt, '(3a, 1pg24.15)') 'YORIGIN ', 'DOUBLE ', 'NDIM 0 # ', this%yorigin
911  txt(lentxt:lentxt) = new_line('a')
912  write (iunit) txt
913  write (txt, '(3a, 1pg24.15)') 'ANGROT ', 'DOUBLE ', 'NDIM 0 # ', this%angrot
914  txt(lentxt:lentxt) = new_line('a')
915  write (iunit) txt
916  write (txt, '(3a, i0)') 'BOTM ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
917  txt(lentxt:lentxt) = new_line('a')
918  write (iunit) txt
919  write (txt, '(3a, i0)') 'IA ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
920  txt(lentxt:lentxt) = new_line('a')
921  write (iunit) txt
922  write (txt, '(3a, i0)') 'JA ', 'INTEGER ', 'NDIM 1 ', this%con%nja
923  txt(lentxt:lentxt) = new_line('a')
924  write (iunit) txt
925  write (txt, '(3a, i0)') 'ICELLTYPE ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
926  txt(lentxt:lentxt) = new_line('a')
927  write (iunit) txt
928  write (txt, '(3a, i0)') 'IDOMAIN ', 'INTEGER ', 'NDIM 1 ', this%nodesuser
929  txt(lentxt:lentxt) = new_line('a')
930  write (iunit) txt
931  !
932  ! -- if vertices have been read then write additional header information
933  if (this%nvert > 0) then
934  write (txt, '(3a, i0)') 'VERTICES ', 'DOUBLE ', 'NDIM 2 2 ', this%nvert
935  txt(lentxt:lentxt) = new_line('a')
936  write (iunit) txt
937  write (txt, '(3a, i0)') 'CELLX ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
938  txt(lentxt:lentxt) = new_line('a')
939  write (iunit) txt
940  write (txt, '(3a, i0)') 'CELLY ', 'DOUBLE ', 'NDIM 1 ', this%nodesuser
941  txt(lentxt:lentxt) = new_line('a')
942  write (iunit) txt
943  write (txt, '(3a, i0)') 'IAVERT ', 'INTEGER ', 'NDIM 1 ', this%nodesuser + 1
944  txt(lentxt:lentxt) = new_line('a')
945  write (iunit) txt
946  write (txt, '(3a, i0)') 'JAVERT ', 'INTEGER ', 'NDIM 1 ', size(this%javert)
947  txt(lentxt:lentxt) = new_line('a')
948  write (iunit) txt
949  end if
950  !
951  ! -- if version 2 write character array headers
952  if (version == 2) then
953  if (found_crs) then
954  write (txt, '(3a, i0)') 'CRS ', 'CHARACTER ', 'NDIM 1 ', &
955  len_trim(crs)
956  txt(lentxt:lentxt) = new_line('a')
957  write (iunit) txt
958  end if
959  end if
960  !
961  ! -- write data
962  write (iunit) this%nodesuser ! nodes
963  write (iunit) this%nja ! nja
964  write (iunit) this%xorigin ! xorigin
965  write (iunit) this%yorigin ! yorigin
966  write (iunit) this%angrot ! angrot
967  write (iunit) this%bottom ! botm
968  write (iunit) this%con%iausr ! ia
969  write (iunit) this%con%jausr ! ja
970  write (iunit) icelltype ! icelltype
971  write (iunit) this%idomain ! idomain
972  !
973  ! -- if vertices have been read then write additional data
974  if (this%nvert > 0) then
975  write (iunit) this%vertices ! vertices
976  write (iunit) (this%cellxy(1, i), i=1, this%nodesuser) ! cellx
977  write (iunit) (this%cellxy(2, i), i=1, this%nodesuser) ! celly
978  write (iunit) this%iavert ! iavert
979  write (iunit) this%javert ! javert
980  end if
981  !
982  ! -- if version 2 write character array data
983  if (version == 2) then
984  if (found_crs) write (iunit) trim(crs) ! crs user input
985  end if
986  !
987  ! -- Close the file
988  close (iunit)
989  end subroutine write_grb
990 
991  !>
992  !! Return a nodenumber from the user specified node number with an
993  !! option to perform a check. This subroutine can be overridden by
994  !! child classes to perform mapping to a model node number
995  !<
996  function get_nodenumber_idx1(this, nodeu, icheck) result(nodenumber)
997  class(disv1dtype), intent(in) :: this
998  integer(I4B), intent(in) :: nodeu
999  integer(I4B), intent(in) :: icheck
1000  integer(I4B) :: nodenumber
1001  !
1002  if (icheck /= 0) then
1003  if (nodeu < 1 .or. nodeu > this%nodes) then
1004  write (errmsg, '(a,i10)') &
1005  'Nodenumber less than 1 or greater than nodes:', nodeu
1006  call store_error(errmsg)
1007  end if
1008  end if
1009  !
1010  ! -- set node number based on whether it is reduced or not
1011  if (this%nodes == this%nodesuser) then
1012  nodenumber = nodeu
1013  else
1014  nodenumber = this%nodereduced(nodeu)
1015  end if
1016  end function get_nodenumber_idx1
1017 
1018  subroutine nodeu_to_string(this, nodeu, str)
1019  ! -- dummy
1020  class(disv1dtype) :: this
1021  integer(I4B), intent(in) :: nodeu
1022  character(len=*), intent(inout) :: str
1023  ! -- local
1024  character(len=10) :: nstr
1025  !
1026  write (nstr, '(i0)') nodeu
1027  str = '('//trim(adjustl(nstr))//')'
1028  end subroutine nodeu_to_string
1029 
1030  !>
1031  !! nodeu_from_string -- Receive a string and convert the string to a user
1032  !! nodenumber. The model is unstructured; just read user nodenumber.
1033  !! If flag_string argument is present and true, the first token in string
1034  !! is allowed to be a string (e.g. boundary name). In this case, if a string
1035  !! is encountered, return value as -2.
1036  !<
1037  function nodeu_from_string(this, lloc, istart, istop, in, iout, line, &
1038  flag_string, allow_zero) result(nodeu)
1039  ! -- dummy
1040  class(disv1dtype) :: this
1041  integer(I4B), intent(inout) :: lloc
1042  integer(I4B), intent(inout) :: istart
1043  integer(I4B), intent(inout) :: istop
1044  integer(I4B), intent(in) :: in
1045  integer(I4B), intent(in) :: iout
1046  character(len=*), intent(inout) :: line
1047  logical, optional, intent(in) :: flag_string
1048  logical, optional, intent(in) :: allow_zero
1049  integer(I4B) :: nodeu
1050  ! -- local
1051  integer(I4B) :: lloclocal, ndum, istat, n
1052  real(dp) :: r
1053  !
1054  if (present(flag_string)) then
1055  if (flag_string) then
1056  ! Check to see if first token in line can be read as an integer.
1057  lloclocal = lloc
1058  call urword(line, lloclocal, istart, istop, 1, ndum, r, iout, in)
1059  read (line(istart:istop), *, iostat=istat) n
1060  if (istat /= 0) then
1061  ! First token in line is not an integer; return flag to this effect.
1062  nodeu = -2
1063  return
1064  end if
1065  end if
1066  end if
1067  !
1068  call urword(line, lloc, istart, istop, 2, nodeu, r, iout, in)
1069  !
1070  if (nodeu == 0) then
1071  if (present(allow_zero)) then
1072  if (allow_zero) then
1073  return
1074  end if
1075  end if
1076  end if
1077  !
1078  if (nodeu < 1 .or. nodeu > this%nodesuser) then
1079  write (errmsg, '(a,i0,a)') &
1080  "Node number in list (", nodeu, ") is outside of the grid. "// &
1081  "Cell number cannot be determined in line '"// &
1082  trim(adjustl(line))//"'."
1083  call store_error(errmsg)
1084  call store_error_filename(this%input_fname)
1085  end if
1086  end function nodeu_from_string
1087 
1088  subroutine disv1d_da(this)
1089  ! -- modules
1092  use simvariablesmodule, only: idm_context
1093  ! -- dummy
1094  class(disv1dtype) :: this
1095  ! -- local
1096  logical(LGP) :: deallocate_vertices
1097  !
1098  ! -- Deallocate idm memory
1099  call memorystore_remove(this%name_model, 'DISV1D', idm_context)
1100  !
1101  ! -- scalars
1102  deallocate_vertices = (this%nvert > 0)
1103  call mem_deallocate(this%nvert)
1104  !
1105  ! -- arrays
1106  call mem_deallocate(this%nodeuser)
1107  call mem_deallocate(this%nodereduced)
1108  call mem_deallocate(this%length)
1109  call mem_deallocate(this%width)
1110  call mem_deallocate(this%bottom)
1111  call mem_deallocate(this%idomain)
1112  !
1113  ! -- cdl hack for arrays for vertices and cell1d blocks
1114  if (deallocate_vertices) then
1115  call mem_deallocate(this%vertices)
1116  call mem_deallocate(this%fdc)
1117  call mem_deallocate(this%cellxy)
1118  call mem_deallocate(this%iavert)
1119  call mem_deallocate(this%javert)
1120  end if
1121  !
1122  ! -- DisBaseType deallocate
1123  call this%DisBaseType%dis_da()
1124  end subroutine disv1d_da
1125 
1126  !> @brief Record a double precision array
1127  !!
1128  !! Record a double precision array. The array will be
1129  !! printed to an external file and/or written to an unformatted external file
1130  !! depending on the argument specifications.
1131  !<
1132  subroutine record_array(this, darray, iout, iprint, idataun, aname, &
1133  cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
1134  ! -- modules
1135  use tdismodule, only: kstp, kper, pertim, totim, delt
1137  ! -- dummy
1138  class(disv1dtype), intent(inout) :: this
1139  real(DP), dimension(:), pointer, contiguous, intent(inout) :: darray !< double precision array to record
1140  integer(I4B), intent(in) :: iout !< unit number for ascii output
1141  integer(I4B), intent(in) :: iprint !< flag indicating whether or not to print the array
1142  integer(I4B), intent(in) :: idataun !< unit number to which the array will be written in binary
1143  character(len=*), intent(in) :: aname !< text descriptor of the array
1144  character(len=*), intent(in) :: cdatafmp ! fortran format for writing the array
1145  integer(I4B), intent(in) :: nvaluesp !< number of values per line for printing
1146  integer(I4B), intent(in) :: nwidthp !< width of the number for printing
1147  character(len=*), intent(in) :: editdesc !< format type (I, G, F, S, E)
1148  real(DP), intent(in) :: dinact !< double precision value to use for cells that are excluded from model domain
1149  ! -- local
1150  integer(I4B) :: k, ifirst
1151  integer(I4B) :: nlay
1152  integer(I4B) :: nrow
1153  integer(I4B) :: ncol
1154  integer(I4B) :: nval
1155  integer(I4B) :: nodeu, noder
1156  integer(I4B) :: istart, istop
1157  real(DP), dimension(:), pointer, contiguous :: dtemp
1158  ! -- formats
1159  character(len=*), parameter :: fmthsv = &
1160  "(1X,/1X,a,' WILL BE SAVED ON UNIT ',I4, &
1161  &' AT END OF TIME STEP',I5,', STRESS PERIOD ',I4)"
1162  !
1163  ! -- set variables
1164  nlay = 1
1165  nrow = 1
1166  ncol = this%mshape(1)
1167  !
1168  ! -- If this is a reduced model, then copy the values from darray into
1169  ! dtemp.
1170  if (this%nodes < this%nodesuser) then
1171  nval = this%nodes
1172  dtemp => this%dbuff
1173  do nodeu = 1, this%nodesuser
1174  noder = this%get_nodenumber(nodeu, 0)
1175  if (noder <= 0) then
1176  dtemp(nodeu) = dinact
1177  cycle
1178  end if
1179  dtemp(nodeu) = darray(noder)
1180  end do
1181  else
1182  nval = this%nodes
1183  dtemp => darray
1184  end if
1185  !
1186  ! -- Print to iout if iprint /= 0
1187  if (iprint /= 0) then
1188  istart = 1
1189  do k = 1, nlay
1190  istop = istart + nrow * ncol - 1
1191  call ulaprufw(ncol, nrow, kstp, kper, k, iout, dtemp(istart:istop), &
1192  aname, cdatafmp, nvaluesp, nwidthp, editdesc)
1193  istart = istop + 1
1194  end do
1195  end if
1196  !
1197  ! -- Save array to an external file.
1198  if (idataun > 0) then
1199  ! -- write to binary file by layer
1200  ifirst = 1
1201  istart = 1
1202  do k = 1, nlay
1203  istop = istart + nrow * ncol - 1
1204  if (ifirst == 1) write (iout, fmthsv) &
1205  trim(adjustl(aname)), idataun, &
1206  kstp, kper
1207  ifirst = 0
1208  call ulasav(dtemp(istart:istop), aname, kstp, kper, &
1209  pertim, totim, ncol, nrow, k, idataun)
1210  istart = istop + 1
1211  end do
1212  elseif (idataun < 0) then
1213  !
1214  ! -- write entire array as one record
1215  call ubdsv1(kstp, kper, aname, -idataun, dtemp, ncol, nrow, nlay, &
1216  iout, delt, pertim, totim)
1217  end if
1218  end subroutine record_array
1219 
1220  !> @brief Record list header using ubdsv06
1221  !<
1222  subroutine record_srcdst_list_header(this, text, textmodel, textpackage, &
1223  dstmodel, dstpackage, naux, auxtxt, &
1224  ibdchn, nlist, iout)
1225  ! -- module
1226  use tdismodule, only: kstp, kper, pertim, totim, delt
1227  use inputoutputmodule, only: ubdsv06
1228  ! -- dummy
1229  class(disv1dtype) :: this
1230  character(len=16), intent(in) :: text
1231  character(len=16), intent(in) :: textmodel
1232  character(len=16), intent(in) :: textpackage
1233  character(len=16), intent(in) :: dstmodel
1234  character(len=16), intent(in) :: dstpackage
1235  integer(I4B), intent(in) :: naux
1236  character(len=16), dimension(:), intent(in) :: auxtxt
1237  integer(I4B), intent(in) :: ibdchn
1238  integer(I4B), intent(in) :: nlist
1239  integer(I4B), intent(in) :: iout
1240  ! -- local
1241  integer(I4B) :: nlay, nrow, ncol
1242  !
1243  nlay = 1
1244  nrow = 1
1245  ncol = this%mshape(1)
1246  !
1247  ! -- Use ubdsv06 to write list header
1248  call ubdsv06(kstp, kper, text, textmodel, textpackage, dstmodel, dstpackage, &
1249  ibdchn, naux, auxtxt, ncol, nrow, nlay, &
1250  nlist, iout, delt, pertim, totim)
1251  end subroutine record_srcdst_list_header
1252 
1253  !> @ brief Calculate the flow width between two cells
1254  !!
1255  !! This should only be called for connections with IHC > 0.
1256  !! Routine is needed, so it can be overridden by the linear
1257  !! network discretization, which allows for a separate flow
1258  !< width for each cell.
1259  subroutine get_flow_width(this, n, m, idx_conn, width_n, width_m)
1260  ! dummy
1261  class(disv1dtype) :: this
1262  integer(I4B), intent(in) :: n !< cell node number
1263  integer(I4B), intent(in) :: m !< cell node number
1264  integer(I4B), intent(in) :: idx_conn !< connection index
1265  real(DP), intent(out) :: width_n !< flow width for cell n
1266  real(DP), intent(out) :: width_m !< flow width for cell m
1267 
1268  ! For disv1d case, width_n and width_m can be different
1269  width_n = this%width(n)
1270  width_m = this%width(m)
1271 
1272  end subroutine get_flow_width
1273 
1274 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:1019
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:1225
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:596
subroutine calculate_cellxy(vertices, fdc, iavert, javert, length, cellxy)
Calculate x, y, coordinates of reach midpoint.
Definition: Disv1d.f90:634
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:534
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:997
subroutine grid_finalize(this)
Finalize grid construction.
Definition: Disv1d.f90:714
subroutine get_flow_width(this, n, m, idx_conn, width_n, width_m)
@ brief Calculate the flow width between two cells
Definition: Disv1d.f90:1260
subroutine record_array(this, darray, iout, iprint, idataun, aname, cdatafmp, nvaluesp, nwidthp, editdesc, dinact)
Record a double precision array.
Definition: Disv1d.f90:1134
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:687
subroutine disv1d_da(this)
Definition: Disv1d.f90:1089
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:818
subroutine allocate_arrays(this)
Definition: Disv1d.f90:794
subroutine write_grb(this, icelltype)
Write binary grid file.
Definition: Disv1d.f90:847
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:1039
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)
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: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
Simplifies tracking parameters sourced from the input context.
Definition: Disv1d.f90:67