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

Data Types

interface  cachedgradienttype
 Decorator that adds caching to any gradient computation implementation. More...
 

Functions/Subroutines

type(cachedgradienttype) function constructor (gradient, dis)
 
subroutine destructor (this)
 
real(dp) function, dimension(3) get (this, n)
 
subroutine set_field (this, phi)
 

Function/Subroutine Documentation

◆ constructor()

type(cachedgradienttype) function cachedgradientmodule::constructor ( class(igradienttype), intent(inout), allocatable  gradient,
class(disbasetype), intent(in), pointer  dis 
)
private

Definition at line 50 of file CachedGradient.f90.

51  ! --dummy
52  class(IGradientType), allocatable, intent(inout) :: gradient
53  class(DisBaseType), pointer, intent(in) :: dis
54  !-- return
55  type(CachedGradientType) :: cached_gradient
56 
57  cached_gradient%dis => dis
58 
59  call move_alloc(gradient, cached_gradient%gradient) ! Take ownership
60  allocate (cached_gradient%cached_gradients(dis%nodes, 3))
61 
Here is the caller graph for this function:

◆ destructor()

subroutine cachedgradientmodule::destructor ( type(cachedgradienttype), intent(inout)  this)
private

Definition at line 64 of file CachedGradient.f90.

65  ! -- dummy
66  type(CachedGradientType), intent(inout) :: this
67 
68  deallocate (this%cached_gradients)
Here is the caller graph for this function:

◆ get()

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

Definition at line 71 of file CachedGradient.f90.

72  ! -- dummy
73  class(CachedGradientType), target :: this
74  integer(I4B), intent(in) :: n
75  !-- return
76  real(DP), dimension(3) :: grad_c
77 
78  grad_c = this%cached_gradients(n, :)

◆ set_field()

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

Definition at line 81 of file CachedGradient.f90.

82  ! -- dummy
83  class(CachedGradientType), target :: this
84  real(DP), dimension(:), pointer, intent(in) :: phi
85  ! -- local
86  integer(I4B) :: n
87 
88  call this%gradient%set_field(phi)
89  do n = 1, this%dis%nodes
90  this%cached_gradients(n, :) = this%gradient%get(n)
91  end do
92