DOpE
Public Member Functions | Protected Member Functions
DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim > Class Template Reference

#include <integrator.h>

Public Member Functions

 Integrator (INTEGRATORDATACONT &idc1)
 
 Integrator (INTEGRATORDATACONT &idc1, INTEGRATORDATACONT &idc2)
 
 ~Integrator ()
 
void ReInit ()
 
template<typename PROBLEM >
void ComputeNonlinearResidual (PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
 
template<typename PROBLEM >
void ComputeNonlinearLhs (PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
 
template<typename PROBLEM >
void ComputeNonlinearRhs (PROBLEM &pde, VECTOR &residual, bool apply_boundary_values=true)
 
template<typename PROBLEM , typename MATRIX >
void ComputeMatrix (PROBLEM &pde, MATRIX &matrix)
 
template<typename PROBLEM >
void ComputeNonlinearAlgebraicResidual (PROBLEM &pde, VECTOR &residual)
 
template<typename PROBLEM >
void ComputeLocalControlConstraints (PROBLEM &pde, VECTOR &constraints)
 
template<typename PROBLEM >
SCALAR ComputeDomainScalar (PROBLEM &pde)
 
template<typename PROBLEM >
SCALAR ComputePointScalar (PROBLEM &pde)
 
template<typename PROBLEM >
SCALAR ComputeBoundaryScalar (PROBLEM &pde)
 
template<typename PROBLEM >
SCALAR ComputeFaceScalar (PROBLEM &pde)
 
template<typename PROBLEM >
SCALAR ComputeAlgebraicScalar (PROBLEM &pde)
 
template<typename PROBLEM >
void ApplyInitialBoundaryValues (PROBLEM &pde, VECTOR &u)
 
template<typename PROBLEM >
void ApplyTransposedInitialBoundaryValues (PROBLEM &pde, VECTOR &u)
 
template<typename PROBLEM >
void ApplyNewtonBoundaryValues (PROBLEM &pde, VECTOR &u)
 
template<typename PROBLEM , typename MATRIX >
void ApplyNewtonBoundaryValues (PROBLEM &pde, MATRIX &matrix, VECTOR &rhs, VECTOR &sol)
 
void AddDomainData (std::string name, const VECTOR *new_data)
 
void DeleteDomainData (std::string name)
 
void AddParamData (std::string name, const dealii::Vector< SCALAR > *new_data)
 
void DeleteParamData (std::string name)
 
template<typename PROBLEM , class STH , class EDC , class FDC >
void ComputeRefinementIndicators (PROBLEM &pde, DWRDataContainer< STH, INTEGRATORDATACONT, EDC, FDC, VECTOR > &dwrc)
 
template<typename PROBLEM >
void ComputeRefinementIndicators (PROBLEM &pde, ResidualErrorContainer< VECTOR > &dwrc)
 
INTEGRATORDATACONT & GetIntegratorDataContainer () const
 
INTEGRATORDATACONT & GetIntegratorDataContainerFunc () const
 

Protected Member Functions

const std::map< std::string,
const VECTOR * > & 
GetDomainData () const
 
const std::map< std::string,
const dealii::Vector< SCALAR > * > & 
GetParamData () const
 
void AddPresetRightHandSide (double s, VECTOR &residual) const
 

Detailed Description

template<typename INTEGRATORDATACONT, typename VECTOR, typename SCALAR, int dim>
class DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >

This class is used to integrate the righthand side, matrix and so on. It assumes that one uses the same triangulation for the control and state variable.

INTEGRATORDATACONT The type of the integratordatacontainer, which has manages the basic data for integration (quadrature, elementdatacontainer, facedatacontainer etc.) VECTOR Class of the vectors which we use in the integrator. SCALAR Type of the scalars we use in the integrator. dim dimesion of the domain

Constructor & Destructor Documentation

template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::Integrator ( INTEGRATORDATACONT &  idc1)
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::Integrator ( INTEGRATORDATACONT &  idc1,
INTEGRATORDATACONT &  idc2 
)

This constructor gets two INTEGRATORDATACONT. One for the evaluation of the integrals during the assembling and residual computation, the second one for the evaluation of functionals.

template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::~Integrator ( )

Member Function Documentation

template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::AddDomainData ( std::string  name,
const VECTOR *  new_data 
)
inline

This function can be used to pass domain data, i.e., finite element functions, to the problem. The added data will automatically be initialized on the elements where the integration takes place.

Adding multiple data with the same name is prohibited.

Parameters
nameAn identifier by which the PROBLEM in the corresponding Compute* method can access the data.
new_dataA pointer to the data to be added.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::AddParamData ( std::string  name,
const dealii::Vector< SCALAR > *  new_data 
)
inline

This function can be used to pass parameter data, i.e., data independent of the spatial position, to the problem.

Adding multiple data with the same name is prohibited.

Parameters
nameAn identifier by which the PROBLEM in the corresponding Compute* method can access the data.
new_dataA pointer to the data to be added.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::AddPresetRightHandSide ( double  s,
VECTOR &  residual 
) const
inlineprotected

This function is used to add usergiven righthandsides to the residual. This is usefull if the rhs is expensive to compute.

The given rhs needs to be given to the integrator as domain data with the name "fixed_rhs"

Parameters
sA scaling parameter
residualThe given residual. Upon exit from this method This vector is residual = residual + s* "fixed_rhs"
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ApplyInitialBoundaryValues ( PROBLEM &  pde,
VECTOR &  u 
)

This method applies inhomogeneous dirichlet boundary values.

It is assumed that PROBLEM provides functions GetDirichletCompMask, GetDirichletColors, and GetDirichletValues to find the appropriate boundary values.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
uThe vector to which the boundary values are applied.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ApplyNewtonBoundaryValues ( PROBLEM &  pde,
VECTOR &  u 
)

This method applies homogeneous dirichlet boundary values on all Dirichlet boundaries

It is assumed that PROBLEM provides functions GetDirichletCompMask and GetDirichletColors to find the appropriate boundary values.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
uThe vector to which the boundary values are applied.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM , typename MATRIX >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ApplyNewtonBoundaryValues ( PROBLEM &  pde,
MATRIX &  matrix,
VECTOR &  rhs,
VECTOR &  sol 
)

This method applies homogeneous dirichlet boundary values on all Dirichlet boundaries. This method applies the required values To allparts of the linear equation Ax=b, i.e., to A, x, and b.

It is assumed that PROBLEM provides functions GetDirichletCompMask and GetDirichletColors to find the appropriate boundary values.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
matrixThe matrix A
solThe vector x
rhsThe vector b
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ApplyTransposedInitialBoundaryValues ( PROBLEM &  pde,
VECTOR &  u 
)

This method is used in optimization problems with control in the dirichlet values, it then needs to calculate the transposed to the control-to-dirichletvalues mapping. Note that this is currently only usable of the control is 0-dimensional, i.e., a fixed number of parameters. Thus it is only implemented in IntegratorMixedDimensions and provided here for compatibility reasons only.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
uThe vector to which the boundary values are applied.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
SCALAR DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeAlgebraicScalar ( PROBLEM &  pde)

This methods evaluates functionals that are given algebraic manipulation of the unknowns.

It is assumed that PROBLEM provides a method ComputeAlgebraicScalar.

This routine assumes that a corresponding calculation method is provided by PROBLEM

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
Returns
The value of the functional
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
SCALAR DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeBoundaryScalar ( PROBLEM &  pde)

This methods evaluates functionals that are given by evaluation of integrals over parts of the boundary.

It is assumed that PROBLEM provides a method BoundaryFunctional.

This routine assumes that a corresponding calculation method is provided by PROBLEM

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
Returns
The value of the functional
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
SCALAR DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeDomainScalar ( PROBLEM &  pde)

This methods evaluates functionals that are given by an integration over the spatial domain.

It is assumed that PROBLEM provides a method ElementFunctional.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
Returns
The value of the functional
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
SCALAR DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeFaceScalar ( PROBLEM &  pde)

This methods evaluates functionals that are given by evaluation of integrals over certain faces in the domain.

It is assumed that PROBLEM provides a method FaceFunctional.

This routine assumes that a corresponding calculation method is provided by PROBLEM

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
Returns
The value of the functional
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeLocalControlConstraints ( PROBLEM &  pde,
VECTOR &  constraints 
)

This method is used to calculate local constraints, i.e., constraints that do not need any integration but can be calculated directly from the values in the unknowns. Typical examples are box constraints in Lagrange elements.

It is assumed that a method ComputeLocalControlConstraints in PROBLEM provides a method for the calculation of the residual.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the nonlinear pde.
constraintsA vector which contains the constraints.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM , typename MATRIX >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeMatrix ( PROBLEM &  pde,
MATRIX &  matrix 
)

This method is used to calculate the matrix corresponding to the linearized equation.

It assumes that the PROBLEM provides methods named ElementMatrix, BoundaryMatrix, FaceMatrix, and InterfaceRhs For the calculation of the respective quantities, see, e.g., OptProblemContainer for details on these methods.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the nonlinear pde.
matrixA matrix which contains the matrix after completing the method.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeNonlinearAlgebraicResidual ( PROBLEM &  pde,
VECTOR &  residual 
)

This routine is used as a dummy to allow for the solutions of problems that do not need any integration, i.e., that don't involve integration.

It is assumed that a method AlgebraicResidual in PROBLEM provides a method for the calculation of the residual.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the nonlinear pde.
resiualA vector which contains the residual after completing the method.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeNonlinearLhs ( PROBLEM &  pde,
VECTOR &  residual,
bool  apply_boundary_values = true 
)

This method is used to evaluate the left hand side of the residual of the nonlinear equation This are all terms depending on the nonlinear variable. The use of this function can be advantageous to avoid recalculation of terms that do not change with update of the evaluation point.

It assumes that the PROBLEM provides methods named ElementEquation,BoundaryEquation, FaceEquation, and InterfaceEquation For the calculation of the respective quantities, see, e.g., OptProblemContainer for details on these methods.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the nonlinear pde.
resiualA vector which contains the residual after completing the method.
apply_boundary_valuesA flag to decide whether or not to set components belonging to constrained DoFs to zero (e.g., at Dirichlet boundaries)
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeNonlinearResidual ( PROBLEM &  pde,
VECTOR &  residual,
bool  apply_boundary_values = true 
)

This method is used to evaluate the residual of the nonlinear equation.

It assumes that the PROBLEM provides methods named ElementEquation, ElementRhs, BoundaryEquation, BoundaryRhs, FaceEquation, FaceRhs, InterfaceEquation, and PointRhs For the calculation of the respective quantities, see, e.g., OptProblemContainer for details on these methods.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the nonlinear pde.
resiualA vector which contains the residual after completing the method.
apply_boundary_valuesA flag to decide whether or not to set components belonging to constrained DoFs to zero (e.g., at Dirichlet boundaries)
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeNonlinearRhs ( PROBLEM &  pde,
VECTOR &  residual,
bool  apply_boundary_values = true 
)

This method is used to evaluate the righ hand side of the residual of the nonlinear equation This are all terms independent of the nonlinear variable. The use of this function can be advantageous to avoid recalculation of terms that do not change with update of the evaluation point.

It assumes that the PROBLEM provides methods named ElementRhs, BoundaryRhs, FaceRhs, and PointRhs For the calculation of the respective quantities, see, e.g., OptProblemContainer for details on these methods.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the nonlinear pde.
resiualA vector which contains the residual after completing the method.
apply_boundary_valuesA flag to decide whether or not to set components belonging to constrained DoFs to zero (e.g., at Dirichlet boundaries)
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
SCALAR DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputePointScalar ( PROBLEM &  pde)

This methods evaluates functionals that are given by evaluation in certain fixed points in the domain.

It is assumed that PROBLEM provides a method PointFunctional.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe object containing the description of the functional.
Returns
The value of the functional
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM , class STH , class EDC , class FDC >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeRefinementIndicators ( PROBLEM &  pde,
DWRDataContainer< STH, INTEGRATORDATACONT, EDC, FDC, VECTOR > &  dwrc 
)

This function is used to calculate indicators based on DWR-type data containes. This means it is assumed that an additional SpaceTimeHandler is used for the calculation of the weights. See DWRDataContainter for details.

Template Parameters
<PROBLEM>The problem description
<STH>An additional SpaceTimeHandler
<EDC>An additional ElementDataContainer
<FDC>An additional FaceDataContainer
Parameters
pdeThe problem
dwrcThe data container for the error estimation. This object also stores the refinement indicators.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
template<typename PROBLEM >
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ComputeRefinementIndicators ( PROBLEM &  pde,
ResidualErrorContainer< VECTOR > &  dwrc 
)

This function is used to calculate indicators based on Residual-type data containes. See ResidualErrorContainer for details.

Template Parameters
<PROBLEM>The problem description
Parameters
pdeThe problem
dwrcThe data container for the error estimation. This object also stores the refinement indicators.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::DeleteDomainData ( std::string  name)
inline

This function is used to remove previously added domain data from the integrator.

Deleting data that is not present in the integrator will cause an exception.

Parameters
nameThe identifier for the data to be removed from the integrator.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::DeleteParamData ( std::string  name)
inline

This function is used to remove previously added parameter data from the integrator.

Deleting data that is not present in the integrator will cause an exception.

Parameters
nameThe identifier for the data to be removed from the integrator.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
const std::map< std::string, const VECTOR * > & DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::GetDomainData ( ) const
inlineprotected

This grants access to the domain data stored in the integrator.

Returns
The map with the stored data.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
INTEGRATORDATACONT & DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::GetIntegratorDataContainer ( ) const
inline
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
INTEGRATORDATACONT & DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::GetIntegratorDataContainerFunc ( ) const
inline
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
const std::map< std::string, const dealii::Vector< SCALAR > * > & DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::GetParamData ( ) const
inlineprotected

This grants access to the parameter data stored in the integrator.

Returns
The map with the stored data.
template<typename INTEGRATORDATACONT , typename VECTOR , typename SCALAR , int dim>
void DOpE::Integrator< INTEGRATORDATACONT, VECTOR, SCALAR, dim >::ReInit ( )

This Function should be called once after grid refinement, or changes in boundary values to recompute sparsity patterns, and constraint matrices.


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