MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
gwt-dsp.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use constantsmodule, only: done, dzero, dhalf, dpi
6  use basedismodule, only: disbasetype
7  use tspfmimodule, only: tspfmitype
8  use xt3dmodule, only: xt3dtype, xt3d_cr
11 
12  implicit none
13  private
14  public :: gwtdsptype
15  public :: dsp_cr
16 
17  type, extends(numericalpackagetype) :: gwtdsptype
18 
19  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() ! pointer to GWT model ibound
20  type(tspfmitype), pointer :: fmi => null() ! pointer to GWT fmi object
21  real(dp), dimension(:), pointer, contiguous :: thetam => null() ! pointer to GWT storage porosity (voids per aquifer volume)
22  real(dp), dimension(:), pointer, contiguous :: diffc => null() ! molecular diffusion coefficient for each cell
23  real(dp), dimension(:), pointer, contiguous :: alh => null() ! longitudinal horizontal dispersivity
24  real(dp), dimension(:), pointer, contiguous :: alv => null() ! longitudinal vertical dispersivity
25  real(dp), dimension(:), pointer, contiguous :: ath1 => null() ! transverse horizontal dispersivity
26  real(dp), dimension(:), pointer, contiguous :: ath2 => null() ! transverse horizontal dispersivity
27  real(dp), dimension(:), pointer, contiguous :: atv => null() ! transverse vertical dispersivity
28  integer(I4B), pointer :: idiffc => null() ! flag indicating diffusion is active
29  integer(I4B), pointer :: idisp => null() ! flag indicating mechanical dispersion is active
30  integer(I4B), pointer :: ialh => null() ! longitudinal horizontal dispersivity data flag
31  integer(I4B), pointer :: ialv => null() ! longitudinal vertical dispersivity data flag
32  integer(I4B), pointer :: iath1 => null() ! transverse horizontal dispersivity data flag
33  integer(I4B), pointer :: iath2 => null() ! transverse horizontal dispersivity data flag
34  integer(I4B), pointer :: iatv => null() ! transverse vertical dispersivity data flag
35  integer(I4B), pointer :: ixt3doff => null() ! xt3d off flag, xt3d is set inactive if 1
36  integer(I4B), pointer :: ixt3drhs => null() ! xt3d rhs flag, xt3d rhs is set active if 1
37  integer(I4B), pointer :: ixt3d => null() ! flag indicating xt3d is active
38  type(xt3dtype), pointer :: xt3d => null() ! xt3d object
39  real(dp), dimension(:), pointer, contiguous :: dispcoef => null() ! disp coefficient (only if xt3d not active)
40  integer(I4B), pointer :: id22 => null() ! flag indicating d22 is available
41  integer(I4B), pointer :: id33 => null() ! flag indicating d33 is available
42  real(dp), dimension(:), pointer, contiguous :: d11 => null() ! dispersion coefficient
43  real(dp), dimension(:), pointer, contiguous :: d22 => null() ! dispersion coefficient
44  real(dp), dimension(:), pointer, contiguous :: d33 => null() ! dispersion coefficient
45  real(dp), dimension(:), pointer, contiguous :: angle1 => null() ! rotation angle 1
46  real(dp), dimension(:), pointer, contiguous :: angle2 => null() ! rotation angle 2
47  real(dp), dimension(:), pointer, contiguous :: angle3 => null() ! rotation angle 3
48  integer(I4B), pointer :: iangle1 => null() ! flag indicating angle1 is available
49  integer(I4B), pointer :: iangle2 => null() ! flag indicating angle2 is available
50  integer(I4B), pointer :: iangle3 => null() ! flag indicating angle3 is available
51 
52  contains
53 
54  procedure :: dsp_df
55  procedure :: dsp_ac
56  procedure :: dsp_mc
57  procedure :: dsp_ar
58  procedure :: dsp_ad
59  procedure :: dsp_fc
60  procedure :: dsp_cq
61  procedure :: dsp_da
62  procedure :: allocate_scalars
63  procedure :: allocate_arrays
64  procedure, private :: source_options
65  procedure, private :: source_griddata
66  procedure, private :: log_options
67  procedure, private :: log_griddata
68  procedure, private :: calcdispellipse
69  procedure, private :: calcdispcoef
70 
71  end type gwtdsptype
72 
73 contains
74 
75  !> @brief Create a DSP object
76  !<
77  subroutine dsp_cr(dspobj, name_model, input_mempath, inunit, iout, fmi)
78  ! -- modules
79  use kindmodule, only: lgp
81  ! -- dummy
82  type(gwtdsptype), pointer :: dspobj
83  character(len=*), intent(in) :: name_model
84  character(len=*), intent(in) :: input_mempath
85  integer(I4B), intent(in) :: inunit
86  integer(I4B), intent(in) :: iout
87  type(tspfmitype), intent(in), target :: fmi
88  ! -- locals
89  ! -- formats
90  character(len=*), parameter :: fmtdsp = &
91  "(1x,/1x,'DSP-- DISPERSION PACKAGE, VERSION 1, 1/24/2018', &
92  &' INPUT READ FROM MEMPATH: ', A, //)"
93  !
94  ! -- Create the object
95  allocate (dspobj)
96  !
97  ! -- create name and memory path
98  call dspobj%set_names(1, name_model, 'DSP', 'DSP', input_mempath)
99  !
100  ! -- Allocate scalars
101  call dspobj%allocate_scalars()
102  !
103  ! -- Set variables
104  dspobj%inunit = inunit
105  dspobj%iout = iout
106  dspobj%fmi => fmi
107  !
108  if (dspobj%inunit > 0) then
109  !
110  ! -- Print a message identifying the dispersion package.
111  if (dspobj%iout > 0) then
112  write (dspobj%iout, fmtdsp) input_mempath
113  end if
114  end if
115  end subroutine dsp_cr
116 
117  !> @brief Define DSP object
118  !<
119  subroutine dsp_df(this, dis, dspOptions)
120  ! -- modules
121  ! -- dummy
122  class(gwtdsptype) :: this
123  class(disbasetype), pointer :: dis
124  type(gwtdspoptionstype), optional, intent(in) :: dspOptions !< the optional DSP options, used when not
125  !! creating DSP from file
126  ! -- local
127  !
128  ! -- Store pointer to dis
129  this%dis => dis
130  !
131  !
132  ! -- set default xt3d representation to on and lhs
133  this%ixt3d = 1
134  !
135  ! -- Read dispersion options
136  if (present(dspoptions)) then
137  this%ixt3d = dspoptions%ixt3d
138  !
139  ! -- Allocate only, grid data will not be read from file
140  call this%allocate_arrays(this%dis%nodes)
141  else
142  !
143  ! -- Source options
144  call this%source_options()
145  call this%allocate_arrays(this%dis%nodes)
146  !
147  ! -- Source dispersion data
148  call this%source_griddata()
149  end if
150  !
151  ! -- xt3d create
152  if (this%ixt3d > 0) then
153  call xt3d_cr(this%xt3d, this%name_model, this%inunit, this%iout, &
154  ldispopt=.true.)
155  this%xt3d%ixt3d = this%ixt3d
156  call this%xt3d%xt3d_df(dis)
157  end if
158  end subroutine dsp_df
159 
160  !> @brief Add connections to DSP
161  !!
162  !! Add connections for extended neighbors to the sparse matrix
163  !<
164  subroutine dsp_ac(this, moffset, sparse)
165  ! -- modules
166  use sparsemodule, only: sparsematrix
168  ! -- dummy
169  class(gwtdsptype) :: this
170  integer(I4B), intent(in) :: moffset
171  type(sparsematrix), intent(inout) :: sparse
172  ! -- local
173  !
174  ! -- Add extended neighbors (neighbors of neighbors)
175  if (this%ixt3d > 0) call this%xt3d%xt3d_ac(moffset, sparse)
176  end subroutine dsp_ac
177 
178  !> @brief Map DSP connections
179  !!
180  !! Map connections and construct iax, jax, and idxglox
181  !<
182  subroutine dsp_mc(this, moffset, matrix_sln)
183  ! -- modules
185  ! -- dummy
186  class(gwtdsptype) :: this
187  integer(I4B), intent(in) :: moffset
188  class(matrixbasetype), pointer :: matrix_sln
189  ! -- local
190  !
191  ! -- Call xt3d map connections
192  if (this%ixt3d > 0) call this%xt3d%xt3d_mc(moffset, matrix_sln)
193  end subroutine dsp_mc
194 
195  !> @brief Allocate and read method for package
196  !!
197  !! Method to allocate and read static data for the package.
198  !<
199  subroutine dsp_ar(this, ibound, thetam)
200  ! -- modules
201  ! -- dummy
202  class(gwtdsptype) :: this
203  integer(I4B), dimension(:), pointer, contiguous :: ibound
204  real(DP), dimension(:), pointer, contiguous :: thetam
205  ! -- local
206  ! -- formats
207  character(len=*), parameter :: fmtdsp = &
208  "(1x,/1x,'DSP-- DISPERSION PACKAGE, VERSION 1, 1/24/2018', &
209  &' INPUT READ FROM UNIT ', i0, //)"
210  !
211  ! -- dsp pointers to arguments that were passed in
212  this%ibound => ibound
213  this%thetam => thetam
214  end subroutine dsp_ar
215 
216  !> @brief Advance method for the package
217  !<
218  subroutine dsp_ad(this)
219  ! -- modules
220  use tdismodule, only: kstp, kper
221  ! -- dummy
222  class(gwtdsptype) :: this
223  ! -- local
224  !
225  ! -- xt3d
226  ! TODO: might consider adding a new mf6 level set pointers method, and
227  ! doing this stuff there instead of in the time step loop.
228  if (kstp * kper == 1) then
229  if (this%ixt3d > 0) then
230  call this%xt3d%xt3d_ar(this%fmi%ibdgwfsat0, this%d11, this%id33, &
231  this%d33, this%fmi%gwfsat, this%id22, this%d22, &
232  this%iangle1, this%iangle2, this%iangle3, &
233  this%angle1, this%angle2, this%angle3)
234  end if
235  end if
236  !
237  ! -- Fill d11, d22, d33, angle1, angle2, angle3 using specific discharge
238  call this%calcdispellipse()
239  !
240  ! -- Recalculate dispersion coefficients if the flows were updated
241  if (this%fmi%iflowsupdated == 1) then
242  if (this%ixt3d == 0) then
243  call this%calcdispcoef()
244  else if (this%ixt3d > 0) then
245  call this%xt3d%xt3d_fcpc(this%dis%nodes, .false.)
246  end if
247  end if
248  end subroutine dsp_ad
249 
250  !> @brief Fill coefficient method for package
251  !!
252  !! Method to calculate and fill coefficients for the package.
253  !<
254  subroutine dsp_fc(this, kiter, nodes, nja, matrix_sln, idxglo, rhs, cnew)
255  ! -- modules
256  ! -- dummy
257  class(gwtdsptype) :: this
258  integer(I4B) :: kiter
259  integer(I4B), intent(in) :: nodes
260  integer(I4B), intent(in) :: nja
261  class(matrixbasetype), pointer :: matrix_sln
262  integer(I4B), intent(in), dimension(nja) :: idxglo
263  real(DP), intent(inout), dimension(nodes) :: rhs
264  real(DP), intent(inout), dimension(nodes) :: cnew
265  ! -- local
266  integer(I4B) :: n, m, idiag, idiagm, ipos, isympos, isymcon
267  real(DP) :: dnm
268  !
269  if (this%ixt3d > 0) then
270  call this%xt3d%xt3d_fc(kiter, matrix_sln, idxglo, rhs, cnew)
271  else
272  do n = 1, nodes
273  if (this%fmi%ibdgwfsat0(n) == 0) cycle
274  idiag = this%dis%con%ia(n)
275  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
276  if (this%dis%con%mask(ipos) == 0) cycle
277  m = this%dis%con%ja(ipos)
278  if (m < n) cycle
279  if (this%fmi%ibdgwfsat0(m) == 0) cycle
280  isympos = this%dis%con%jas(ipos)
281  dnm = this%dispcoef(isympos)
282  !
283  ! -- Contribution to row n
284  call matrix_sln%add_value_pos(idxglo(ipos), dnm)
285  call matrix_sln%add_value_pos(idxglo(idiag), -dnm)
286  !
287  ! -- Contribution to row m
288  idiagm = this%dis%con%ia(m)
289  isymcon = this%dis%con%isym(ipos)
290  call matrix_sln%add_value_pos(idxglo(isymcon), dnm)
291  call matrix_sln%add_value_pos(idxglo(idiagm), -dnm)
292  end do
293  end do
294  end if
295  end subroutine dsp_fc
296 
297  !> @ brief Calculate flows for package
298  !!
299  !! Method to calculate dispersion contribution to flowja
300  !<
301  subroutine dsp_cq(this, cnew, flowja)
302  ! -- modules
303  ! -- dummy
304  class(gwtdsptype) :: this
305  real(DP), intent(inout), dimension(:) :: cnew
306  real(DP), intent(inout), dimension(:) :: flowja
307  ! -- local
308  integer(I4B) :: n, m, ipos, isympos
309  real(DP) :: dnm
310  !
311  ! -- Calculate dispersion and add to flowja
312  if (this%ixt3d > 0) then
313  call this%xt3d%xt3d_flowja(cnew, flowja)
314  else
315  do n = 1, this%dis%nodes
316  if (this%fmi%ibdgwfsat0(n) == 0) cycle
317  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
318  m = this%dis%con%ja(ipos)
319  if (this%fmi%ibdgwfsat0(m) == 0) cycle
320  isympos = this%dis%con%jas(ipos)
321  dnm = this%dispcoef(isympos)
322  flowja(ipos) = flowja(ipos) + dnm * (cnew(m) - cnew(n))
323  end do
324  end do
325  end if
326  end subroutine dsp_cq
327 
328  !> @ brief Allocate scalar variables for package
329  !!
330  !! Method to allocate scalar variables for the package.
331  !<
332  subroutine allocate_scalars(this)
333  ! -- modules
335  use constantsmodule, only: dzero
336  ! -- dummy
337  class(gwtdsptype) :: this
338  ! -- local
339  !
340  ! -- allocate scalars in NumericalPackageType
341  call this%NumericalPackageType%allocate_scalars()
342  !
343  ! -- Allocate
344  call mem_allocate(this%idiffc, 'IDIFFC', this%memoryPath)
345  call mem_allocate(this%idisp, 'IDISP', this%memoryPath)
346  call mem_allocate(this%ialh, 'IALH', this%memoryPath)
347  call mem_allocate(this%ialv, 'IALV', this%memoryPath)
348  call mem_allocate(this%iath1, 'IATH1', this%memoryPath)
349  call mem_allocate(this%iath2, 'IATH2', this%memoryPath)
350  call mem_allocate(this%iatv, 'IATV', this%memoryPath)
351  call mem_allocate(this%ixt3doff, 'IXT3DOFF', this%memoryPath)
352  call mem_allocate(this%ixt3drhs, 'IXT3DRHS', this%memoryPath)
353  call mem_allocate(this%ixt3d, 'IXT3D', this%memoryPath)
354  call mem_allocate(this%id22, 'ID22', this%memoryPath)
355  call mem_allocate(this%id33, 'ID33', this%memoryPath)
356  call mem_allocate(this%iangle1, 'IANGLE1', this%memoryPath)
357  call mem_allocate(this%iangle2, 'IANGLE2', this%memoryPath)
358  call mem_allocate(this%iangle3, 'IANGLE3', this%memoryPath)
359  !
360  ! -- Initialize
361  this%idiffc = 0
362  this%idisp = 0
363  this%ialh = 0
364  this%ialv = 0
365  this%iath1 = 0
366  this%iath2 = 0
367  this%iatv = 0
368  this%ixt3doff = 0
369  this%ixt3drhs = 0
370  this%ixt3d = 0
371  this%id22 = 1
372  this%id33 = 1
373  this%iangle1 = 1
374  this%iangle2 = 1
375  this%iangle3 = 1
376  end subroutine allocate_scalars
377 
378  !> @ brief Allocate arrays for package
379  !!
380  !! Method to allocate arrays for the package.
381  !<
382  subroutine allocate_arrays(this, nodes)
383  ! -- modules
385  use constantsmodule, only: dzero
386  ! -- dummy
387  class(gwtdsptype) :: this
388  integer(I4B), intent(in) :: nodes
389  ! -- local
390  !
391  ! -- Allocate
392  call mem_allocate(this%alh, nodes, 'ALH', trim(this%memoryPath))
393  call mem_allocate(this%alv, nodes, 'ALV', trim(this%memoryPath))
394  call mem_allocate(this%ath1, nodes, 'ATH1', trim(this%memoryPath))
395  call mem_allocate(this%ath2, nodes, 'ATH2', trim(this%memoryPath))
396  call mem_allocate(this%atv, nodes, 'ATV', trim(this%memoryPath))
397  call mem_allocate(this%diffc, nodes, 'DIFFC', trim(this%memoryPath))
398  call mem_allocate(this%d11, nodes, 'D11', trim(this%memoryPath))
399  call mem_allocate(this%d22, nodes, 'D22', trim(this%memoryPath))
400  call mem_allocate(this%d33, nodes, 'D33', trim(this%memoryPath))
401  call mem_allocate(this%angle1, nodes, 'ANGLE1', trim(this%memoryPath))
402  call mem_allocate(this%angle2, nodes, 'ANGLE2', trim(this%memoryPath))
403  call mem_allocate(this%angle3, nodes, 'ANGLE3', trim(this%memoryPath))
404  !
405  ! -- Allocate dispersion coefficient array if xt3d not in use
406  if (this%ixt3d == 0) then
407  call mem_allocate(this%dispcoef, this%dis%njas, 'DISPCOEF', &
408  trim(this%memoryPath))
409  else
410  call mem_allocate(this%dispcoef, 0, 'DISPCOEF', trim(this%memoryPath))
411  end if
412  end subroutine allocate_arrays
413 
414  !> @ brief Deallocate memory
415  !!
416  !! Method to deallocate memory for the package.
417  !<
418  subroutine dsp_da(this)
419  ! -- modules
423  ! -- dummy
424  class(gwtdsptype) :: this
425  ! -- local
426  !
427  ! -- Deallocate input memory
428  call memorystore_remove(this%name_model, 'DSP', idm_context)
429  !
430  ! -- deallocate arrays
431  if (this%inunit /= 0) then
432  call mem_deallocate(this%alh, 'ALH', trim(this%memoryPath))
433  call mem_deallocate(this%alv, 'ALV', trim(this%memoryPath))
434  call mem_deallocate(this%ath1, 'ATH1', trim(this%memoryPath))
435  call mem_deallocate(this%ath2, 'ATH2', trim(this%memoryPath))
436  call mem_deallocate(this%atv, 'ATV', trim(this%memoryPath))
437  call mem_deallocate(this%diffc)
438  call mem_deallocate(this%d11)
439  call mem_deallocate(this%d22)
440  call mem_deallocate(this%d33)
441  call mem_deallocate(this%angle1)
442  call mem_deallocate(this%angle2)
443  call mem_deallocate(this%angle3)
444  call mem_deallocate(this%dispcoef)
445  if (this%ixt3d > 0) call this%xt3d%xt3d_da()
446  end if
447  !
448  ! -- deallocate objects
449  if (this%ixt3d > 0) deallocate (this%xt3d)
450  !
451  ! -- deallocate scalars
452  call mem_deallocate(this%idiffc)
453  call mem_deallocate(this%idisp)
454  call mem_deallocate(this%ialh)
455  call mem_deallocate(this%ialv)
456  call mem_deallocate(this%iath1)
457  call mem_deallocate(this%iath2)
458  call mem_deallocate(this%iatv)
459  call mem_deallocate(this%ixt3doff)
460  call mem_deallocate(this%ixt3drhs)
461  call mem_deallocate(this%ixt3d)
462  call mem_deallocate(this%id22)
463  call mem_deallocate(this%id33)
464  call mem_deallocate(this%iangle1)
465  call mem_deallocate(this%iangle2)
466  call mem_deallocate(this%iangle3)
467  !
468  ! -- deallocate variables in NumericalPackageType
469  call this%NumericalPackageType%da()
470  !
471  ! -- nullify pointers
472  this%thetam => null()
473  end subroutine dsp_da
474 
475  !> @brief Write user options to list file
476  !<
477  subroutine log_options(this, found)
479  class(gwtdsptype) :: this
480  type(gwtdspparamfoundtype), intent(in) :: found
481 
482  write (this%iout, '(1x,a)') 'Setting DSP Options'
483  write (this%iout, '(4x,a,i0)') 'XT3D formulation [0=INACTIVE, 1=ACTIVE, &
484  &3=ACTIVE RHS] set to: ', this%ixt3d
485  write (this%iout, '(1x,a,/)') 'End Setting DSP Options'
486  end subroutine log_options
487 
488  !> @brief Update simulation mempath options
489  !<
490  subroutine source_options(this)
491  ! -- modules
492  !use KindModule, only: LGP
495  ! -- dummy
496  class(gwtdsptype) :: this
497  ! -- locals
498  type(gwtdspparamfoundtype) :: found
499  !
500  ! -- update defaults with idm sourced values
501  call mem_set_value(this%ixt3doff, 'XT3D_OFF', this%input_mempath, &
502  found%xt3d_off)
503  call mem_set_value(this%ixt3drhs, 'XT3D_RHS', this%input_mempath, &
504  found%xt3d_rhs)
505  !
506  ! -- set xt3d state flag
507  if (found%xt3d_off) this%ixt3d = 0
508  if (found%xt3d_rhs) this%ixt3d = 2
509  !
510  ! -- log options
511  if (this%iout > 0) then
512  call this%log_options(found)
513  end if
514  end subroutine source_options
515 
516  !> @brief Write dimensions to list file
517  !<
518  subroutine log_griddata(this, found)
520  class(gwtdsptype) :: this
521  type(gwtdspparamfoundtype), intent(in) :: found
522 
523  write (this%iout, '(1x,a)') 'Setting DSP Griddata'
524 
525  if (found%diffc) then
526  write (this%iout, '(4x,a)') 'DIFFC set from input file'
527  end if
528 
529  if (found%alh) then
530  write (this%iout, '(4x,a)') 'ALH set from input file'
531  end if
532 
533  if (found%alv) then
534  write (this%iout, '(4x,a)') 'ALV set from input file'
535  end if
536 
537  if (found%ath1) then
538  write (this%iout, '(4x,a)') 'ATH1 set from input file'
539  end if
540 
541  if (found%ath2) then
542  write (this%iout, '(4x,a)') 'ATH2 set from input file'
543  end if
544 
545  if (found%atv) then
546  write (this%iout, '(4x,a)') 'ATV set from input file'
547  end if
548 
549  write (this%iout, '(1x,a,/)') 'End Setting DSP Griddata'
550 
551  end subroutine log_griddata
552 
553  !> @brief Update DSP simulation data from input mempath
554  !<
555  subroutine source_griddata(this)
556  ! -- modules
560  use constantsmodule, only: linelength
562  ! -- dummy
563  class(gwtdsptype) :: this
564  ! -- locals
565  character(len=LINELENGTH) :: errmsg
566  type(gwtdspparamfoundtype) :: found
567  integer(I4B), dimension(:), pointer, contiguous :: map
568  ! -- formats
569  !
570  ! -- set map
571  map => null()
572  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
573  !
574  ! -- update defaults with idm sourced values
575  call mem_set_value(this%diffc, 'DIFFC', this%input_mempath, map, found%diffc)
576  call mem_set_value(this%alh, 'ALH', this%input_mempath, map, found%alh)
577  call mem_set_value(this%alv, 'ALV', this%input_mempath, map, found%alv)
578  call mem_set_value(this%ath1, 'ATH1', this%input_mempath, map, found%ath1)
579  call mem_set_value(this%ath2, 'ATH2', this%input_mempath, map, found%ath2)
580  call mem_set_value(this%atv, 'ATV', this%input_mempath, map, found%atv)
581  !
582  ! -- set active flags
583  if (found%diffc) this%idiffc = 1
584  if (found%alh) this%ialh = 1
585  if (found%alv) this%ialv = 1
586  if (found%ath1) this%iath1 = 1
587  if (found%ath2) this%iath2 = 1
588  if (found%atv) this%iatv = 1
589  !
590  ! -- reallocate diffc if not found
591  if (.not. found%diffc) then
592  call mem_reallocate(this%diffc, 0, 'DIFFC', trim(this%memoryPath))
593  end if
594  !
595  ! -- set this%idisp flag
596  if (found%alh) this%idisp = this%idisp + 1
597  if (found%alv) this%idisp = this%idisp + 1
598  if (found%ath1) this%idisp = this%idisp + 1
599  if (found%ath2) this%idisp = this%idisp + 1
600  !
601  ! -- manage dispersion arrays
602  if (this%idisp > 0) then
603  if (.not. (found%alh .and. found%ath1)) then
604  write (errmsg, '(a)') &
605  'if dispersivities are specified then ALH and ATH1 are required.'
606  call store_error(errmsg)
607  end if
608  ! -- If alv not specified then point it to alh
609  if (.not. found%alv) &
610  call mem_reassignptr(this%alv, 'ALV', trim(this%memoryPath), &
611  'ALH', trim(this%memoryPath))
612  ! -- If ath2 not specified then point it to ath1
613  if (.not. found%ath2) &
614  call mem_reassignptr(this%ath2, 'ATH2', trim(this%memoryPath), &
615  'ATH1', trim(this%memoryPath))
616  ! -- If atv not specified then point it to ath2
617  if (.not. found%atv) &
618  call mem_reassignptr(this%atv, 'ATV', trim(this%memoryPath), &
619  'ATH2', trim(this%memoryPath))
620  else
621  call mem_reallocate(this%alh, 0, 'ALH', trim(this%memoryPath))
622  call mem_reallocate(this%alv, 0, 'ALV', trim(this%memoryPath))
623  call mem_reallocate(this%ath1, 0, 'ATH1', trim(this%memoryPath))
624  call mem_reallocate(this%ath2, 0, 'ATH2', trim(this%memoryPath))
625  call mem_reallocate(this%atv, 0, 'ATV', trim(this%memoryPath))
626  end if
627  !
628  ! -- log griddata
629  if (this%iout > 0) then
630  call this%log_griddata(found)
631  end if
632  end subroutine source_griddata
633 
634  !> @brief Calculate dispersion coefficients
635  !<
636  subroutine calcdispellipse(this)
637  ! -- modules
638  ! -- dummy
639  class(gwtdsptype) :: this
640  ! -- local
641  integer(I4B) :: nodes, n
642  real(DP) :: q, qx, qy, qz
643  real(DP) :: alh, alv, ath1, ath2, atv, a
644  real(DP) :: al, at1, at2
645  real(DP) :: qzoqsquared
646  real(DP) :: dstar
647  !
648  ! -- loop through and calculate dispersion coefficients and angles
649  nodes = size(this%d11)
650  do n = 1, nodes
651  !
652  ! -- initialize
653  this%d11(n) = dzero
654  this%d22(n) = dzero
655  this%d33(n) = dzero
656  this%angle1(n) = dzero
657  this%angle2(n) = dzero
658  this%angle3(n) = dzero
659  if (this%fmi%ibdgwfsat0(n) == 0) cycle
660  !
661  ! -- specific discharge
662  qx = dzero
663  qy = dzero
664  qz = dzero
665  q = dzero
666  qx = this%fmi%gwfspdis(1, n)
667  qy = this%fmi%gwfspdis(2, n)
668  qz = this%fmi%gwfspdis(3, n)
669  q = qx**2 + qy**2 + qz**2
670  if (q > dzero) q = sqrt(q)
671  !
672  ! -- dispersion coefficients
673  alh = dzero
674  alv = dzero
675  ath1 = dzero
676  ath2 = dzero
677  atv = dzero
678  if (this%idisp > 0) then
679  alh = this%alh(n)
680  alv = this%alv(n)
681  ath1 = this%ath1(n)
682  ath2 = this%ath2(n)
683  atv = this%atv(n)
684  end if
685  dstar = dzero
686  if (this%idiffc > 0) then
687  ! -- Multiply diffusion coefficient by mobile porosity, defined
688  ! as volume of voids in the mobile domain per aquifer volume
689  dstar = this%diffc(n) * this%thetam(n)
690  end if
691  !
692  ! -- Calculate the longitudal and transverse dispersivities
693  al = dzero
694  at1 = dzero
695  at2 = dzero
696  if (q > dzero) then
697  qzoqsquared = (qz / q)**2
698  al = alh * (done - qzoqsquared) + alv * qzoqsquared
699  at1 = ath1 * (done - qzoqsquared) + atv * qzoqsquared
700  at2 = ath2 * (done - qzoqsquared) + atv * qzoqsquared
701  end if
702  !
703  ! -- Calculate and save the diagonal components of the dispersion tensor
704  this%d11(n) = al * q + dstar
705  this%d22(n) = at1 * q + dstar
706  this%d33(n) = at2 * q + dstar
707  !
708  ! -- Angles of rotation if velocity based dispersion tensor
709  if (this%idisp > 0) then
710  !
711  ! -- angles of rotation from model coordinates to direction of velocity
712  ! qx / q = cos(a1) * cos(a2)
713  ! qy / q = sin(a1) * cos(a2)
714  ! qz / q = sin(a2)
715  !
716  ! -- angle3 is zero
717  this%angle3(n) = dzero
718  !
719  ! -- angle2
720  a = dzero
721  if (q > dzero) a = qz / q
722  this%angle2(n) = asin(a)
723  !
724  ! -- angle1
725  a = q * cos(this%angle2(n))
726  if (a /= dzero) then
727  a = qx / a
728  else
729  a = dzero
730  end if
731  !
732  ! -- acos(1) not defined, so set to zero if necessary
733  if (a <= -done) then
734  this%angle1(n) = dpi
735  elseif (a >= done) then
736  this%angle1(n) = dzero
737  else
738  this%angle1(n) = acos(a)
739  end if
740  !
741  end if
742  end do
743  end subroutine calcdispellipse
744 
745  !> @brief Calculate dispersion coefficients
746  !<
747  subroutine calcdispcoef(this)
748  ! -- modules
749  use hgeoutilmodule, only: hyeff
751  ! -- dummy
752  class(gwtdsptype) :: this
753  ! -- local
754  integer(I4B) :: nodes, n, m, idiag, ipos
755  real(DP) :: clnm, clmn, dn, dm
756  real(DP) :: vg1, vg2, vg3
757  integer(I4B) :: ihc, isympos
758  integer(I4B) :: iavgmeth
759  real(DP) :: satn, satm, topn, topm, botn, botm
760  real(DP) :: hwva, cond, cn, cm, denom
761  real(DP) :: anm, amn, thksatn, thksatm
762  !
763  ! -- set iavgmeth = 1 to use arithmetic averaging for effective dispersion
764  iavgmeth = 1
765  !
766  ! -- Process connections
767  nodes = size(this%d11)
768  do n = 1, nodes
769  if (this%fmi%ibdgwfsat0(n) == 0) cycle
770  idiag = this%dis%con%ia(n)
771  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
772  !
773  ! -- Set m to connected cell
774  m = this%dis%con%ja(ipos)
775  !
776  ! -- skip for lower triangle
777  if (m < n) cycle
778  isympos = this%dis%con%jas(ipos)
779  this%dispcoef(isympos) = dzero
780  if (this%fmi%ibdgwfsat0(m) == 0) cycle
781  !
782  ! -- cell dimensions
783  hwva = this%dis%con%hwva(isympos)
784  clnm = this%dis%con%cl1(isympos)
785  clmn = this%dis%con%cl2(isympos)
786  ihc = this%dis%con%ihc(isympos)
787  topn = this%dis%top(n)
788  topm = this%dis%top(m)
789  botn = this%dis%bot(n)
790  botm = this%dis%bot(m)
791  !
792  ! -- flow model information
793  satn = this%fmi%gwfsat(n)
794  satm = this%fmi%gwfsat(m)
795  !
796  ! -- Calculate dispersion coefficient for cell n in the direction
797  ! normal to the shared n-m face and for cell m in the direction
798  ! normal to the shared n-m face.
799  call this%dis%connection_normal(n, m, ihc, vg1, vg2, vg3, ipos)
800  dn = hyeff(this%d11(n), this%d22(n), this%d33(n), &
801  this%angle1(n), this%angle2(n), this%angle3(n), &
802  vg1, vg2, vg3, iavgmeth)
803  dm = hyeff(this%d11(m), this%d22(m), this%d33(m), &
804  this%angle1(m), this%angle2(m), this%angle3(m), &
805  vg1, vg2, vg3, iavgmeth)
806  !
807  ! -- Calculate dispersion conductance based on NPF subroutines and the
808  ! effective dispersion coefficients dn and dm.
809  if (ihc == 0) then
810  clnm = satn * (topn - botn) * dhalf
811  clmn = satm * (topm - botm) * dhalf
812  anm = hwva
813  !
814  ! -- n is convertible and unsaturated
815  if (satn == dzero) then
816  anm = dzero
817  else if (n > m .and. satn < done) then
818  anm = dzero
819  end if
820  !
821  ! -- m is convertible and unsaturated
822  if (satm == dzero) then
823  anm = dzero
824  else if (m > n .and. satm < done) then
825  anm = dzero
826  end if
827  !
828  ! -- amn is the same as anm for vertical flow
829  amn = anm
830  !
831  else
832  !
833  ! -- horizontal conductance
834  !
835  ! -- handle vertically staggered case
836  if (ihc == 2) then
837  thksatn = staggered_thkfrac(topn, botn, satn, topm, botm)
838  thksatm = staggered_thkfrac(topm, botm, satm, topn, botn)
839  else
840  thksatn = (topn - botn) * satn
841  thksatm = (topm - botm) * satm
842  end if
843  !
844  ! -- calculate the saturated area term
845  anm = thksatn * hwva
846  amn = thksatm * hwva
847  !
848  ! -- n or m is unsaturated, so no dispersion
849  if (satn == dzero .or. satm == dzero) then
850  anm = dzero
851  amn = dzero
852  end if
853  !
854  end if
855  !
856  ! -- calculate conductance using the two half cell conductances
857  cn = dzero
858  if (clnm > dzero) cn = dn * anm / clnm
859  cm = dzero
860  if (clmn > dzero) cm = dm * amn / clmn
861  denom = cn + cm
862  if (denom > dzero) then
863  cond = cn * cm / denom
864  else
865  cond = dzero
866  end if
867  !
868  ! -- Assign the calculated dispersion conductance
869  this%dispcoef(isympos) = cond
870  !
871  end do
872  end do
873  end subroutine calcdispcoef
874 
875 end module gwtdspmodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dpi
real constant
Definition: Constants.f90:128
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module contains stateless conductance functions.
real(dp) function, public staggered_thkfrac(top, bot, sat, topc, botc)
Calculate the thickness fraction for staggered grids.
subroutine log_options(this, found)
Write user options to list file.
Definition: gwt-dsp.f90:478
subroutine allocate_scalars(this)
@ brief Allocate scalar variables for package
Definition: gwt-dsp.f90:333
subroutine dsp_df(this, dis, dspOptions)
Define DSP object.
Definition: gwt-dsp.f90:120
subroutine dsp_da(this)
@ brief Deallocate memory
Definition: gwt-dsp.f90:419
subroutine, public dsp_cr(dspobj, name_model, input_mempath, inunit, iout, fmi)
Create a DSP object.
Definition: gwt-dsp.f90:78
subroutine log_griddata(this, found)
Write dimensions to list file.
Definition: gwt-dsp.f90:519
subroutine dsp_fc(this, kiter, nodes, nja, matrix_sln, idxglo, rhs, cnew)
Fill coefficient method for package.
Definition: gwt-dsp.f90:255
subroutine calcdispcoef(this)
Calculate dispersion coefficients.
Definition: gwt-dsp.f90:748
subroutine dsp_ar(this, ibound, thetam)
Allocate and read method for package.
Definition: gwt-dsp.f90:200
subroutine dsp_ad(this)
Advance method for the package.
Definition: gwt-dsp.f90:219
subroutine dsp_mc(this, moffset, matrix_sln)
Map DSP connections.
Definition: gwt-dsp.f90:183
subroutine source_griddata(this)
Update DSP simulation data from input mempath.
Definition: gwt-dsp.f90:556
subroutine dsp_cq(this, cnew, flowja)
@ brief Calculate flows for package
Definition: gwt-dsp.f90:302
subroutine allocate_arrays(this, nodes)
@ brief Allocate arrays for package
Definition: gwt-dsp.f90:383
subroutine source_options(this)
Update simulation mempath options.
Definition: gwt-dsp.f90:491
subroutine calcdispellipse(this)
Calculate dispersion coefficients.
Definition: gwt-dsp.f90:637
subroutine dsp_ac(this, moffset, sparse)
Add connections to DSP.
Definition: gwt-dsp.f90:165
General-purpose hydrogeologic functions.
Definition: HGeoUtil.f90:2
real(dp) function, public hyeff(k11, k22, k33, ang1, ang2, ang3, vg1, vg2, vg3, iavgmeth)
Calculate the effective horizontal hydraulic conductivity from an ellipse using a specified direction...
Definition: HGeoUtil.f90:31
This module defines variable data types.
Definition: kind.f90:8
subroutine, public memorystore_remove(component, subcomponent, context)
This module contains the base numerical package type.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
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
subroutine, public xt3d_cr(xt3dobj, name_model, inunit, iout, ldispopt)
Create a new xt3d object.
data structure (and helpers) for passing dsp option data