MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
UpstreamScheme.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
3  use constantsmodule, only: dzero
6  use basedismodule, only: disbasetype
7  use tspfmimodule, only: tspfmitype
8 
9  implicit none
10  private
11 
12  public :: upstreamschemetype
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
22  end type upstreamschemetype
23 
24  interface upstreamschemetype
25  module procedure constructor
26  end interface upstreamschemetype
27 
28 contains
29  function constructor(dis, fmi) result(interpolation_scheme)
30  ! -- return
31  type(upstreamschemetype) :: 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 upstream interpolation computations.
45  !<
46  subroutine set_field(this, phi)
47  ! -- dummy
48  class(upstreamschemetype), 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(upstreamschemetype), 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) :: qnm
64 
65  ! -- Compute the coefficients for the upwind scheme
66  qnm = this%fmi%gwfflowja(iposnm)
67 
68  if (qnm < dzero) then
69  phi_face%c_n = 1.0_dp
70  else
71  phi_face%c_m = 1.0_dp
72  end if
73 
74  end function compute
75 
76 end module upstreamschememodule
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
This module defines variable data types.
Definition: kind.f90:8
type(coefficientstype) function compute(this, n, m, iposnm)
type(upstreamschemetype) function constructor(dis, fmi)
subroutine set_field(this, phi)
Set the scalar field for which interpolation will be computed.