MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
LocalCellExtrema.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b
3  use basedismodule, only: disbasetype
4 
5  implicit none
6  private
7 
8  public :: localcellextrematype
9 
10  !> @brief Computes and caches local extrema for each cell and its connected neighbors
11  !!
12  !! This class computes the minimum and maximum values within the local
13  !! stencil of each cell (the cell itself plus all directly connected neighboring cells).
14  !! The extrema are computed once when the scalar field is set and then cached for
15  !! fast retrieval. This is particularly useful for TVD limiters and slope
16  !! limiting algorithms that need to enforce monotonicity constraints.
17  !!
18  !! The local extrema computation follows the connectivity pattern defined by the
19  !! discretization object, examining all cells that share a face with the target cell.
20  !! This creates a computational stencil that includes the central cell and its
21  !! immediate neighbors.
22  !!
23  !! @note The get_min() and get_max() methods return pointers for zero-copy performance
24  !! in performance-critical loops. Callers should treat returned values as read-only.
25  !<
27  private
28  class(disbasetype), pointer :: dis
29  real(dp), dimension(:), allocatable :: min
30  real(dp), dimension(:), allocatable :: max
31  contains
32  procedure :: set_field
33  procedure :: get_min
34  procedure :: get_max
35  final :: destructor
36 
37  procedure, private :: compute_local_extrema
38  procedure, private :: find_local_extrema
39  end type localcellextrematype
40 
41  interface localcellextrematype
42  module procedure constructor
43  end interface localcellextrematype
44 
45 contains
46  function constructor(dis) result(local_extrema)
47  ! -- return
48  type(localcellextrematype) :: local_extrema
49  ! -- dummy
50  class(disbasetype), pointer, intent(in) :: dis
51 
52  local_extrema%dis => dis
53 
54  allocate (local_extrema%min(dis%nodes))
55  allocate (local_extrema%max(dis%nodes))
56  end function constructor
57 
58  subroutine destructor(this)
59  ! -- dummy
60  type(localcellextrematype), intent(inout) :: this
61 
62  deallocate (this%min)
63  deallocate (this%max)
64  end subroutine destructor
65 
66  subroutine set_field(this, phi)
67  ! -- dummy
68  class(localcellextrematype), target :: this
69  real(DP), intent(in), dimension(:), pointer :: phi
70 
71  call this%compute_local_extrema(phi)
72  end subroutine set_field
73 
74  function get_min(this, n) result(min_val)
75  ! -- dummy
76  class(localcellextrematype), target :: this
77  integer(I4B), intent(in) :: n ! Node index
78  !-- return
79  real(dp), pointer :: min_val
80 
81  min_val => this%min(n)
82  end function get_min
83 
84  function get_max(this, n) result(max_val)
85  ! -- dummy
86  class(localcellextrematype), target :: this
87  integer(I4B), intent(in) :: n ! Node index
88  !-- return
89  real(dp), pointer :: max_val
90 
91  max_val => this%max(n)
92  end function get_max
93 
94  subroutine compute_local_extrema(this, phi)
95  ! -- dummy
96  class(localcellextrematype), target :: this
97  real(DP), intent(in), dimension(:) :: phi
98  ! -- local
99  integer(I4B) :: n
100  real(DP) :: min_phi, max_phi
101 
102  do n = 1, this%dis%nodes
103  call this%find_local_extrema(n, phi, min_phi, max_phi)
104  this%min(n) = min_phi
105  this%max(n) = max_phi
106  end do
107  end subroutine compute_local_extrema
108 
109  subroutine find_local_extrema(this, n, phi, min_phi, max_phi)
110  ! -- dummy
111  class(localcellextrematype), target :: this
112  integer(I4B), intent(in) :: n
113  real(DP), intent(in), dimension(:) :: phi
114  real(DP), intent(out) :: min_phi, max_phi
115  ! -- local
116  integer(I4B) :: ipos, m
117 
118  min_phi = phi(n)
119  max_phi = phi(n)
120  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
121  m = this%dis%con%ja(ipos)
122  min_phi = min(min_phi, phi(m))
123  max_phi = max(max_phi, phi(m))
124  end do
125  end subroutine
126 
127 end module localcellextremamodule
This module defines variable data types.
Definition: kind.f90:8
type(localcellextrematype) function constructor(dis)
subroutine set_field(this, phi)
subroutine compute_local_extrema(this, phi)
real(dp) function, pointer get_min(this, n)
real(dp) function, pointer get_max(this, n)
subroutine destructor(this)
subroutine find_local_extrema(this, n, phi, min_phi, max_phi)
Computes and caches local extrema for each cell and its connected neighbors.