MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
LeastSquaresGradient.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
3  use constantsmodule, only: done
4 
5  Use igradient
6  use basedismodule, only: disbasetype
7  use pseudoinversemodule, only: pinv
9 
10  implicit none
11  private
12 
13  public :: leastsquaresgradienttype
14 
15  type array2d
16  real(dp), dimension(:, :), allocatable :: data
17  end type array2d
18 
19  !> @brief Weighted least-squares gradient method for structured and unstructured grids.
20  !!
21  !! This class implements a least-squares gradient reconstruction for use on both structured and unstructured grids.
22  !! For each cell, it precomputes and caches a gradient reconstruction matrix using the Moore-Penrose pseudoinverse,
23  !! based on the geometry and connectivity of the mesh. The operator is created once during initialization
24  !! and can then be efficiently applied to any scalar field to compute the gradient in each cell.
25  !! The gradient can then be computed by multiplying the reconstruction matrix with the difference vector.
26  !! ∇ɸ = R * ∑(ɸ_i - ɸ_up), where i are the neighboring cells.
27  !!
28  !! - The gradient operator is constructed using normalized direction vectors between cell centers,
29  !! scaled by the inverse of the distance.
30  !! - The least-squares approach ensures robust gradients even for irregular or rank-deficient stencils.
31  !! - The operator is cached for each cell, so gradient computation is efficient for repeated queries.
32  !! - The class provides a `get` method to compute the gradient for any cell and scalar field.
33  !!
34  !! @note Boundary cells are not handled in a special manner. This may impact the quality of the gradient
35  !! near boundaries, especially if a cell does not have enough neighbors (fewer than three in 3D).
36  !<
38  class(disbasetype), pointer :: dis
39  type(array2d), allocatable, dimension(:) :: r ! Gradient reconstruction matrix
40  contains
41  procedure :: get
42 
43  procedure, private :: compute_cell_gradient
46 
47  interface leastsquaresgradienttype
48  module procedure constructor
49  end interface leastsquaresgradienttype
50 
51 contains
52  function constructor(dis) Result(gradient)
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
69  end function constructor
70 
71  function create_gradient_reconstruction_matrix(this, n) result(R)
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 
113 
114  function get(this, n, c) result(grad_c)
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)
123  end function get
124 
125  function compute_cell_gradient(this, n, phi_new) result(grad_c)
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 
154  end function compute_cell_gradient
155 
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
integer(i4b) function, public number_connected_faces(dis, n)
Returns the number of connected faces for a given cell.
Definition: DisUtils.f90:22
real(dp) function, dimension(3), public node_distance(dis, n, m)
Returns the vector distance from cell n to cell m.
Definition: DisUtils.f90:39
This module defines variable data types.
Definition: kind.f90:8
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)
type(leastsquaresgradienttype) function constructor(dis)
real(dp) function, dimension(size(a, dim=2), size(a, dim=1)), public pinv(A)
Abstract interface for cell-based gradient computation.
Definition: IGradient.f90:15
Weighted least-squares gradient method for structured and unstructured grids.