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

Data Types

interface  localcellextrematype
 Computes and caches local extrema for each cell and its connected neighbors. More...
 

Functions/Subroutines

type(localcellextrematype) function constructor (dis)
 
subroutine destructor (this)
 
subroutine set_field (this, phi)
 
real(dp) function, pointer get_min (this, n)
 
real(dp) function, pointer get_max (this, n)
 
subroutine compute_local_extrema (this, phi)
 
subroutine find_local_extrema (this, n, phi, min_phi, max_phi)
 

Function/Subroutine Documentation

◆ compute_local_extrema()

subroutine localcellextremamodule::compute_local_extrema ( class(localcellextrematype), target  this,
real(dp), dimension(:), intent(in)  phi 
)
private

Definition at line 94 of file LocalCellExtrema.f90.

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

◆ constructor()

type(localcellextrematype) function localcellextremamodule::constructor ( class(disbasetype), intent(in), pointer  dis)
private

Definition at line 46 of file LocalCellExtrema.f90.

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))
Here is the caller graph for this function:

◆ destructor()

subroutine localcellextremamodule::destructor ( type(localcellextrematype), intent(inout)  this)
private

Definition at line 58 of file LocalCellExtrema.f90.

59  ! -- dummy
60  type(LocalCellExtremaType), intent(inout) :: this
61 
62  deallocate (this%min)
63  deallocate (this%max)
Here is the caller graph for this function:

◆ find_local_extrema()

subroutine localcellextremamodule::find_local_extrema ( class(localcellextrematype), target  this,
integer(i4b), intent(in)  n,
real(dp), dimension(:), intent(in)  phi,
real(dp), intent(out)  min_phi,
real(dp), intent(out)  max_phi 
)
private

Definition at line 109 of file LocalCellExtrema.f90.

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

◆ get_max()

real(dp) function, pointer localcellextremamodule::get_max ( class(localcellextrematype), target  this,
integer(i4b), intent(in)  n 
)
private

Definition at line 84 of file LocalCellExtrema.f90.

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)

◆ get_min()

real(dp) function, pointer localcellextremamodule::get_min ( class(localcellextrematype), target  this,
integer(i4b), intent(in)  n 
)
private

Definition at line 74 of file LocalCellExtrema.f90.

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)

◆ set_field()

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

Definition at line 66 of file LocalCellExtrema.f90.

67  ! -- dummy
68  class(LocalCellExtremaType), target :: this
69  real(DP), intent(in), dimension(:), pointer :: phi
70 
71  call this%compute_local_extrema(phi)