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