MODFLOW 6  version 6.7.0.dev2
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, c)
 
real(dp) function, dimension(3) compute_cell_gradient (this, n, phi_new)
 

Function/Subroutine Documentation

◆ compute_cell_gradient()

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

Definition at line 125 of file LeastSquaresGradient.f90.

126  ! -- return
127  real(DP), dimension(3) :: grad_c
128  ! -- dummy
129  class(LeastSquaresGradientType), target :: this
130  integer(I4B), intent(in) :: n
131  real(DP), dimension(:), intent(in) :: phi_new
132  ! -- local
133  real(DP), dimension(:, :), pointer :: R
134  integer(I4B) :: ipos, local_pos
135  integer(I4B) :: number_connections
136 
137  integer(I4B) :: m
138  real(DP), dimension(:), allocatable :: dc
139 
140  ! Assemble the concentration difference vector
141  number_connections = number_connected_faces(this%dis, n)
142  allocate (dc(number_connections))
143  local_pos = 1
144  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
145  m = this%dis%con%ja(ipos)
146  dc(local_pos) = phi_new(m) - phi_new(n)
147  local_pos = local_pos + 1
148  end do
149 
150  ! Compute the cells gradient
151  r => this%R(n)%data
152  grad_c = matmul(r, dc)
153 
Here is the call graph for this function:

◆ constructor()

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

Definition at line 52 of file LeastSquaresGradient.f90.

53  ! --dummy
54  class(DisBaseType), pointer, intent(in) :: dis
55  !-- return
56  type(LeastSquaresGradientType) :: gradient
57  ! -- local
58  integer(I4B) :: n, nodes
59 
60  gradient%dis => dis
61  nodes = dis%nodes
62 
63  ! -- Compute the gradient rec
64  nodes = dis%nodes
65  allocate (gradient%R(dis%nodes))
66  do n = 1, nodes
67  gradient%R(n)%data = gradient%create_gradient_reconstruction_matrix(n)
68  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 71 of file LeastSquaresGradient.f90.

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

Definition at line 114 of file LeastSquaresGradient.f90.

115  ! -- dummy
116  class(LeastSquaresGradientType), target :: this
117  integer(I4B), intent(in) :: n
118  real(DP), dimension(:), intent(in) :: c
119  !-- return
120  real(DP), dimension(3) :: grad_c
121 
122  grad_c = this%compute_cell_gradient(n, c)