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