MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
leastsquaresgradientmodule Module Reference

Data Types

type  array2d
 
interface  leastsquaresgradienttype
 Weighted least-squares gradient method for structured and unstructured grids. More...
 

Functions/Subroutines

type(leastsquaresgradienttype) function constructor (dis)
 
real(dp) function, dimension(:, :), allocatable create_gradient_reconstruction_matrix (this, n)
 
real(dp) function, dimension(3) get (this, n)
 
subroutine set_field (this, phi)
 
real(dp) function, dimension(3) compute_cell_gradient (this, n)
 

Function/Subroutine Documentation

◆ compute_cell_gradient()

real(dp) function, dimension(3) leastsquaresgradientmodule::compute_cell_gradient ( class(leastsquaresgradienttype), target  this,
integer(i4b), intent(in)  n 
)
private

Definition at line 140 of file LeastSquaresGradient.f90.

141  ! -- return
142  real(DP), dimension(3) :: grad_c
143  ! -- dummy
144  class(LeastSquaresGradientType), target :: this
145  integer(I4B), intent(in) :: n
146  ! -- local
147  real(DP), dimension(:, :), pointer :: R
148  integer(I4B) :: ipos, local_pos
149  integer(I4B) :: number_connections
150 
151  integer(I4B) :: m
152  real(DP), dimension(:), allocatable :: dc
153 
154  ! Assemble the concentration difference vector
155  number_connections = number_connected_faces(this%dis, n)
156  allocate (dc(number_connections))
157  local_pos = 1
158  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
159  m = this%dis%con%ja(ipos)
160  dc(local_pos) = this%phi(m) - this%phi(n)
161  local_pos = local_pos + 1
162  end do
163 
164  ! Compute the cells gradient
165  r => this%R(n)%data
166  grad_c = matmul(r, dc)
167 
Here is the call graph for this function:

◆ constructor()

type(leastsquaresgradienttype) function leastsquaresgradientmodule::constructor ( class(disbasetype), intent(in), pointer  dis)
private

Definition at line 60 of file LeastSquaresGradient.f90.

61  ! --dummy
62  class(DisBaseType), pointer, intent(in) :: dis
63  !-- return
64  type(LeastSquaresGradientType) :: gradient
65  ! -- local
66  integer(I4B) :: n, nodes
67 
68  gradient%dis => dis
69  nodes = dis%nodes
70 
71  ! -- Compute the gradient rec
72  nodes = dis%nodes
73  allocate (gradient%R(dis%nodes))
74  do n = 1, nodes
75  gradient%R(n)%data = gradient%create_gradient_reconstruction_matrix(n)
76  end do
Here is the caller graph for this function:

◆ create_gradient_reconstruction_matrix()

real(dp) function, dimension(:, :), allocatable leastsquaresgradientmodule::create_gradient_reconstruction_matrix ( class(leastsquaresgradienttype this,
integer(i4b), intent(in)  n 
)
private

Definition at line 79 of file LeastSquaresGradient.f90.

80  ! -- dummy
81  class(LeastSquaresGradientType) :: this
82  integer(I4B), intent(in) :: n ! Cell index for which to create the operator
83  real(DP), dimension(:, :), allocatable :: R ! The resulting gradient reconstruction matrix (3 x number_connections)
84  ! -- local
85  integer(I4B) :: number_connections ! Number of connected neighboring cells
86  integer(I4B) :: ipos, local_pos, m ! Loop indices and neighbor cell index
87  real(DP) :: length ! Distance between cell centers
88  real(DP), dimension(3) :: dnm ! Vector from cell n to neighbor m
89  real(DP), dimension(:, :), allocatable :: d ! Matrix of normalized direction vectors (number_connections x 3)
90  real(DP), dimension(:, :), allocatable :: inverse_distance ! Diagonal scaling matrix (number_connections x number_connections),
91  ! where each diagonal entry is the inverse of the distance between
92 
93  number_connections = number_connected_faces(this%dis, n)
94 
95  allocate (d(number_connections, 3))
96  allocate (r(3, number_connections))
97  allocate (inverse_distance(number_connections, number_connections))
98 
99  inverse_distance = 0
100  d = 0
101 
102  ! Assemble the distance matrix
103  ! Handle the internal connections
104  local_pos = 1
105  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
106  m = this%dis%con%ja(ipos)
107 
108  dnm = node_distance(this%dis, n, m)
109  length = norm2(dnm)
110 
111  d(local_pos, :) = dnm / length
112  inverse_distance(local_pos, local_pos) = 1.0_dp / length
113 
114  local_pos = local_pos + 1
115  end do
116 
117  ! Compute the gradient reconstructions matrix
118  r = matmul(pinv(d), inverse_distance)
119 
Here is the call graph for this function:

◆ get()

real(dp) function, dimension(3) leastsquaresgradientmodule::get ( class(leastsquaresgradienttype), target  this,
integer(i4b), intent(in)  n 
)
private

Definition at line 122 of file LeastSquaresGradient.f90.

123  ! -- dummy
124  class(LeastSquaresGradientType), target :: this
125  integer(I4B), intent(in) :: n
126  !-- return
127  real(DP), dimension(3) :: grad_c
128 
129  grad_c = this%compute_cell_gradient(n)

◆ set_field()

subroutine leastsquaresgradientmodule::set_field ( class(leastsquaresgradienttype), target  this,
real(dp), dimension(:), intent(in), pointer  phi 
)
private

Definition at line 132 of file LeastSquaresGradient.f90.

133  ! -- dummy
134  class(LeastSquaresGradientType), target :: this
135  real(DP), dimension(:), pointer, intent(in) :: phi
136 
137  this%phi => phi