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