Calculates the Jacobian of the local residual for finite volume spatial discretizations using a finite difference method. More...

#include <fvbasefdlocallinearizer.hh>

Public Member Functions

 FvBaseFdLocalLinearizer ()
 
 ~FvBaseFdLocalLinearizer ()
 
void init (Simulator &simulator)
 Initialize the local Jacobian object. More...
 
void linearize (const Element &element)
 Compute an element's local Jacobian matrix and evaluate its residual. More...
 
void linearize (ElementContext &elemCtx, const Element &elem)
 Compute an element's local Jacobian matrix and evaluate its residual. More...
 
Scalar numericEpsilon (const ElementContext &elemCtx, unsigned dofIdx, unsigned pvIdx) const
 Returns the epsilon value which is added and removed from the current solution. More...
 
LocalResidual & localResidual ()
 Return reference to the local residual. More...
 
const LocalResidual & localResidual () const
 Return reference to the local residual. More...
 
const ScalarMatrixBlock & jacobian (unsigned domainScvIdx, unsigned rangeScvIdx) const
 Returns the local Jacobian matrix of the residual of a sub-control volume. More...
 
const ScalarVectorBlock & residual (unsigned dofIdx) const
 Returns the local residual of a sub-control volume. More...
 

Static Public Member Functions

static void registerParameters ()
 Register all run-time parameters for the local jacobian. More...
 
static Scalar baseEpsilon ()
 Returns the unweighted epsilon value used to calculate the local derivatives. More...
 

Protected Member Functions

Implementation & asImp_ ()
 
const Implementation & asImp_ () const
 
const Simulator & simulator_ () const
 
const Problem & problem_ () const
 
const Model & model_ () const
 
void resize_ (const ElementContext &elemCtx)
 Resize all internal attributes to the size of the element. More...
 
void reset_ (const ElementContext &elemCtx)
 Reset the all relevant internal attributes to 0. More...
 
void evalPartialDerivative_ (ElementContext &elemCtx, unsigned dofIdx, unsigned pvIdx)
 Compute the partial derivatives of a context's residual functions. More...
 
void updateLocalJacobian_ (const ElementContext &elemCtx, unsigned focusDofIdx, unsigned pvIdx)
 Updates the current local Jacobian matrix with the partial derivatives of all equations for primary variable 'pvIdx' at the degree of freedom associated with 'focusDofIdx'. More...
 

Static Protected Member Functions

static int numericDifferenceMethod_ ()
 Returns the numeric difference method which is applied. More...
 

Protected Attributes

Simulator * simulatorPtr_
 
Model * modelPtr_
 
ElementContext * internalElemContext_
 
LocalEvalBlockVector residual_
 
LocalEvalBlockVector derivResidual_
 
ScalarLocalBlockMatrix jacobian_
 
LocalResidual localResidual_
 

Detailed Description

template<class TypeTag>
class Opm::FvBaseFdLocalLinearizer< TypeTag >

Calculates the Jacobian of the local residual for finite volume spatial discretizations using a finite difference method.

The local Jacobian for a given context is defined as the derivatives of the residuals of all degrees of freedom featured by the stencil with regard to the primary variables of the stencil's "primary" degrees of freedom.

This class implements numeric differentiation using finite difference methods, i.e. forward or backward differences (2nd order), or central differences (3rd order). The method used is determined by the "NumericDifferenceMethod" property:

  • If the value of this property is smaller than 0, backward differences are used, i.e.:

    \[
    \frac{\partial f(x)}{\partial x} \approx \frac{f(x) - f(x - \epsilon)}{\epsilon}
  \]

  • If the value of this property is 0, central differences are used, i.e.:

    \[
    \frac{\partial f(x)}{\partial x} \approx
         \frac{f(x + \epsilon) - f(x - \epsilon)}{2 \epsilon}
  \]

  • if the value of this property is larger than 0, forward differences are used, i.e.:

    \[
    \frac{\partial f(x)}{\partial x} \approx
         \frac{f(x + \epsilon) - f(x)}{\epsilon}
  \]

Here, $ f $ is the residual function for all equations, $x$ is the value of a sub-control volume's primary variable at the evaluation point and $\epsilon$ is a small scalar value larger than 0.

Constructor & Destructor Documentation

◆ FvBaseFdLocalLinearizer()

template<class TypeTag >
Opm::FvBaseFdLocalLinearizer< TypeTag >::FvBaseFdLocalLinearizer ( )
inline

◆ ~FvBaseFdLocalLinearizer()

Member Function Documentation

◆ asImp_() [1/2]

template<class TypeTag >
Implementation & Opm::FvBaseFdLocalLinearizer< TypeTag >::asImp_ ( )
inlineprotected

◆ asImp_() [2/2]

template<class TypeTag >
const Implementation & Opm::FvBaseFdLocalLinearizer< TypeTag >::asImp_ ( ) const
inlineprotected

◆ baseEpsilon()

template<class TypeTag >
static Scalar Opm::FvBaseFdLocalLinearizer< TypeTag >::baseEpsilon ( )
inlinestatic

Returns the unweighted epsilon value used to calculate the local derivatives.

Referenced by Opm::FvBaseFdLocalLinearizer< TypeTag >::numericEpsilon().

◆ evalPartialDerivative_()

template<class TypeTag >
void Opm::FvBaseFdLocalLinearizer< TypeTag >::evalPartialDerivative_ ( ElementContext &  elemCtx,
unsigned  dofIdx,
unsigned  pvIdx 
)
inlineprotected

Compute the partial derivatives of a context's residual functions.

This method can be overwritten by the implementation if a better scheme than numerical differentiation is available.

The default implementation of this method uses numeric differentiation, i.e. forward or backward differences (2nd order), or central differences (3rd order). The method used is determined by the "NumericDifferenceMethod" property:

  • If the value of this property is smaller than 0, backward differences are used, i.e.:

    \[
    \frac{\partial f(x)}{\partial x} \approx \frac{f(x) - f(x - \epsilon)}{\epsilon}
  \]

  • If the value of this property is 0, central differences are used, i.e.:

    \[
    \frac{\partial f(x)}{\partial x} \approx
         \frac{f(x + \epsilon) - f(x - \epsilon)}{2 \epsilon}
  \]

  • if the value of this property is larger than 0, forward differences are used, i.e.:

    \[
    \frac{\partial f(x)}{\partial x} \approx \frac{f(x + \epsilon) - f(x)}{\epsilon}
  \]

Here, $ f $ is the residual function for all equations, $x$ is the value of a sub-control volume's primary variable at the evaluation point and $\epsilon$ is a small value larger than 0.

Parameters
elemCtxThe element context for which the local partial derivative ought to be calculated
dofIdxThe sub-control volume index of the current finite element for which the partial derivative ought to be calculated
pvIdxThe index of the primary variable at the dofIdx' sub-control volume of the current finite element for which the partial derivative ought to be calculated

References Opm::FvBaseFdLocalLinearizer< TypeTag >::asImp_(), Opm::FvBaseFdLocalLinearizer< TypeTag >::derivResidual_, Opm::FvBaseFdLocalLinearizer< TypeTag >::localResidual_, Opm::FvBaseFdLocalLinearizer< TypeTag >::numericDifferenceMethod_(), and Opm::FvBaseFdLocalLinearizer< TypeTag >::residual_.

◆ init()

template<class TypeTag >
void Opm::FvBaseFdLocalLinearizer< TypeTag >::init ( Simulator &  simulator)
inline

Initialize the local Jacobian object.

At this point we can assume that everything has been allocated, although some objects may not yet be completely initialized.

Parameters
simulatorThe simulator object of the simulation.

References Opm::FvBaseFdLocalLinearizer< TypeTag >::internalElemContext_, and Opm::FvBaseFdLocalLinearizer< TypeTag >::simulatorPtr_.

◆ jacobian()

template<class TypeTag >
const ScalarMatrixBlock & Opm::FvBaseFdLocalLinearizer< TypeTag >::jacobian ( unsigned  domainScvIdx,
unsigned  rangeScvIdx 
) const
inline

Returns the local Jacobian matrix of the residual of a sub-control volume.

Parameters
domainScvIdxThe local index of the sub control volume which contains the independents
rangeScvIdxThe local index of the sub control volume which contains the local residual

References Opm::FvBaseFdLocalLinearizer< TypeTag >::jacobian_.

◆ linearize() [1/2]

template<class TypeTag >
void Opm::FvBaseFdLocalLinearizer< TypeTag >::linearize ( const Element &  element)
inline

Compute an element's local Jacobian matrix and evaluate its residual.

The local Jacobian for a given context is defined as the derivatives of the residuals of all degrees of freedom featured by the stencil with regard to the primary variables of the stencil's "primary" degrees of freedom. Adding the local Jacobians for all elements in the grid will give the global Jacobian 'grad f(x)'.

Parameters
elementThe grid element for which the local residual and its local Jacobian should be calculated.

References Opm::FvBaseFdLocalLinearizer< TypeTag >::internalElemContext_, and Opm::FvBaseFdLocalLinearizer< TypeTag >::linearize().

Referenced by Opm::FvBaseFdLocalLinearizer< TypeTag >::linearize().

◆ linearize() [2/2]

template<class TypeTag >
void Opm::FvBaseFdLocalLinearizer< TypeTag >::linearize ( ElementContext &  elemCtx,
const Element &  elem 
)
inline

Compute an element's local Jacobian matrix and evaluate its residual.

The local Jacobian for a given context is defined as the derivatives of the residuals of all degrees of freedom featured by the stencil with regard to the primary variables of the stencil's "primary" degrees of freedom. Adding the local Jacobians for all elements in the grid will give the global Jacobian 'grad f(x)'.

After calling this method the ElementContext is in an undefined state, so do not use it anymore!

Parameters
elemCtxThe element execution context for which the local residual and its local Jacobian should be calculated.

References Opm::FvBaseFdLocalLinearizer< TypeTag >::asImp_(), Opm::FvBaseFdLocalLinearizer< TypeTag >::localResidual_, Opm::FvBaseFdLocalLinearizer< TypeTag >::model_(), Opm::FvBaseFdLocalLinearizer< TypeTag >::reset_(), Opm::FvBaseFdLocalLinearizer< TypeTag >::residual_, Opm::FvBaseFdLocalLinearizer< TypeTag >::resize_(), and Opm::FvBaseFdLocalLinearizer< TypeTag >::updateLocalJacobian_().

◆ localResidual() [1/2]

template<class TypeTag >
LocalResidual & Opm::FvBaseFdLocalLinearizer< TypeTag >::localResidual ( )
inline

Return reference to the local residual.

References Opm::FvBaseFdLocalLinearizer< TypeTag >::localResidual_.

◆ localResidual() [2/2]

template<class TypeTag >
const LocalResidual & Opm::FvBaseFdLocalLinearizer< TypeTag >::localResidual ( ) const
inline

Return reference to the local residual.

References Opm::FvBaseFdLocalLinearizer< TypeTag >::localResidual_.

◆ model_()

template<class TypeTag >
const Model & Opm::FvBaseFdLocalLinearizer< TypeTag >::model_ ( ) const
inlineprotected

◆ numericDifferenceMethod_()

template<class TypeTag >
static int Opm::FvBaseFdLocalLinearizer< TypeTag >::numericDifferenceMethod_ ( )
inlinestaticprotected

Returns the numeric difference method which is applied.

Referenced by Opm::FvBaseFdLocalLinearizer< TypeTag >::evalPartialDerivative_().

◆ numericEpsilon()

template<class TypeTag >
Scalar Opm::FvBaseFdLocalLinearizer< TypeTag >::numericEpsilon ( const ElementContext &  elemCtx,
unsigned  dofIdx,
unsigned  pvIdx 
) const
inline

Returns the epsilon value which is added and removed from the current solution.

Parameters
elemCtxThe element execution context for which the local residual and its gradient should be calculated.
dofIdxThe local index of the element's vertex for which the local derivative ought to be calculated.
pvIdxThe index of the primary variable which gets varied

References Opm::FvBaseFdLocalLinearizer< TypeTag >::baseEpsilon().

◆ problem_()

template<class TypeTag >
const Problem & Opm::FvBaseFdLocalLinearizer< TypeTag >::problem_ ( ) const
inlineprotected

◆ registerParameters()

template<class TypeTag >
static void Opm::FvBaseFdLocalLinearizer< TypeTag >::registerParameters ( )
inlinestatic

Register all run-time parameters for the local jacobian.

◆ reset_()

template<class TypeTag >
void Opm::FvBaseFdLocalLinearizer< TypeTag >::reset_ ( const ElementContext &  elemCtx)
inlineprotected

◆ residual()

template<class TypeTag >
const ScalarVectorBlock & Opm::FvBaseFdLocalLinearizer< TypeTag >::residual ( unsigned  dofIdx) const
inline

Returns the local residual of a sub-control volume.

Parameters
dofIdxThe local index of the sub control volume

References Opm::FvBaseFdLocalLinearizer< TypeTag >::residual_.

◆ resize_()

template<class TypeTag >
void Opm::FvBaseFdLocalLinearizer< TypeTag >::resize_ ( const ElementContext &  elemCtx)
inlineprotected

◆ simulator_()

template<class TypeTag >
const Simulator & Opm::FvBaseFdLocalLinearizer< TypeTag >::simulator_ ( ) const
inlineprotected

◆ updateLocalJacobian_()

template<class TypeTag >
void Opm::FvBaseFdLocalLinearizer< TypeTag >::updateLocalJacobian_ ( const ElementContext &  elemCtx,
unsigned  focusDofIdx,
unsigned  pvIdx 
)
inlineprotected

Updates the current local Jacobian matrix with the partial derivatives of all equations for primary variable 'pvIdx' at the degree of freedom associated with 'focusDofIdx'.

References Opm::FvBaseFdLocalLinearizer< TypeTag >::derivResidual_, and Opm::FvBaseFdLocalLinearizer< TypeTag >::jacobian_.

Referenced by Opm::FvBaseFdLocalLinearizer< TypeTag >::linearize().

Member Data Documentation

◆ derivResidual_

◆ internalElemContext_

◆ jacobian_

◆ localResidual_

◆ modelPtr_

template<class TypeTag >
Model* Opm::FvBaseFdLocalLinearizer< TypeTag >::modelPtr_
protected

◆ residual_

◆ simulatorPtr_


The documentation for this class was generated from the following file: