MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
CentralDifferenceScheme.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
3  use constantsmodule, only: dhalf, done
6  use basedismodule, only: disbasetype
7  use tspfmimodule, only: tspfmitype
8 
9  implicit none
10  private
11 
13 
15  private
16  class(disbasetype), pointer :: dis
17  type(tspfmitype), pointer :: fmi
18  real(dp), dimension(:), pointer :: phi
19  contains
20  procedure :: compute
21  procedure :: set_field
23 
25  module procedure constructor
26  end interface centraldifferenceschemetype
27 
28 contains
29  function constructor(dis, fmi) result(interpolation_scheme)
30  ! -- return
31  type(centraldifferenceschemetype) :: interpolation_scheme
32  ! --dummy
33  class(disbasetype), pointer, intent(in) :: dis
34  type(tspfmitype), pointer, intent(in) :: fmi
35 
36  interpolation_scheme%dis => dis
37  interpolation_scheme%fmi => fmi
38 
39  end function constructor
40 
41  !> @brief Set the scalar field for which interpolation will be computed
42  !!
43  !! This method establishes a pointer to the scalar field data for
44  !! subsequent central difference interpolation computations.
45  !<
46  subroutine set_field(this, phi)
47  ! -- dummy
48  class(centraldifferenceschemetype), target :: this
49  real(DP), intent(in), dimension(:), pointer :: phi
50 
51  this%phi => phi
52  end subroutine set_field
53 
54  function compute(this, n, m, iposnm) result(phi_face)
55  !-- return
56  type(coefficientstype) :: phi_face
57  ! -- dummy
58  class(centraldifferenceschemetype), target :: this
59  integer(I4B), intent(in) :: n
60  integer(I4B), intent(in) :: m
61  integer(I4B), intent(in) :: iposnm
62  ! -- local
63  real(dp) :: lnm, lmn, omega
64 
65  ! -- calculate weight based on distances between nodes and the shared
66  ! face of the connection
67  if (this%dis%con%ihc(this%dis%con%jas(iposnm)) == 0) then
68  ! -- vertical connection; assume cell is fully saturated
69  lnm = dhalf * (this%dis%top(n) - this%dis%bot(n))
70  lmn = dhalf * (this%dis%top(m) - this%dis%bot(m))
71  else
72  ! -- horizontal connection
73  lnm = this%dis%con%cl1(this%dis%con%jas(iposnm))
74  lmn = this%dis%con%cl2(this%dis%con%jas(iposnm))
75  end if
76 
77  omega = lmn / (lnm + lmn)
78 
79  phi_face%c_n = omega
80  phi_face%c_m = done - omega
81 
82  end function compute
83 
type(coefficientstype) function compute(this, n, m, iposnm)
subroutine set_field(this, phi)
Set the scalar field for which interpolation will be computed.
type(centraldifferenceschemetype) function constructor(dis, fmi)
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module defines variable data types.
Definition: kind.f90:8