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