DOpE
Public Member Functions | Data Fields | Protected Attributes
DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim > Class Template Reference

#include <pdeinterface.h>

Public Member Functions

 PDEInterface ()
 
virtual ~PDEInterface ()
 
virtual void ElementEquation (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void StrongElementResidual (const EDC< DH, VECTOR, dealdim > &, const EDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void ElementTimeEquation (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquation_U (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquation_UT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquation_UTT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquationExplicit (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquationExplicit_U (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquationExplicit_UT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquationExplicit_UTT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementTimeEquationExplicit_UU (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementEquation_U (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void StrongElementResidual_U (const EDC< DH, VECTOR, dealdim > &, const EDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void ElementEquation_UT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementEquation_UTT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementEquation_Q (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementEquation_QT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double scale_ico)
 
virtual void ElementEquation_QTT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementEquation_UU (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementEquation_QU (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementEquation_UQ (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementEquation_QQ (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void ElementRightHandSide (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ElementMatrix (const EDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void ElementTimeMatrix (const EDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &)
 
virtual void ElementTimeMatrix_T (const EDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &)
 
virtual void ElementTimeMatrixExplicit (const EDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &)
 
virtual void ElementTimeMatrixExplicit_T (const EDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &)
 
virtual void ElementMatrix_T (const EDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void ControlElementEquation (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ControlElementMatrix (const EDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double)
 
virtual void ControlBoundaryEquation (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void ControlBoundaryMatrix (const FDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double)
 
virtual void StrongElementResidual_Control (const EDC< DH, VECTOR, dealdim > &, const EDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void StrongFaceResidual_Control (const FDC< DH, VECTOR, dealdim > &, const FDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void StrongBoundaryResidual_Control (const FDC< DH, VECTOR, dealdim > &, const FDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void FaceEquation (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void StrongFaceResidual (const FDC< DH, VECTOR, dealdim > &, const FDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void FaceEquation_U (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void StrongFaceResidual_U (const FDC< DH, VECTOR, dealdim > &, const FDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void FaceEquation_UT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_UTT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_Q (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_QT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_QTT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_UU (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_QU (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_UQ (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceEquation_QQ (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void FaceRightHandSide (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void FaceMatrix (const FDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void FaceMatrix_T (const FDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void InterfaceMatrix (const FDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void InterfaceMatrix_T (const FDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void InterfaceEquation (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void InterfaceEquation_U (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void StrongBoundaryResidual (const FDC< DH, VECTOR, dealdim > &, const FDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void BoundaryEquation_U (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void StrongBoundaryResidual_U (const FDC< DH, VECTOR, dealdim > &, const FDC< DH, VECTOR, dealdim > &, double &, double)
 
virtual void BoundaryEquation_UT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_UTT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_Q (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_QT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_QTT (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_UU (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_QU (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_UQ (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryEquation_QQ (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double, double)
 
virtual void BoundaryRightHandSide (const FDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void BoundaryMatrix (const FDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void BoundaryMatrix_T (const FDC< DH, VECTOR, dealdim > &, dealii::FullMatrix< double > &, double, double)
 
virtual void Init_ElementEquation (const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale, double)
 
virtual void Init_ElementRhs_Q (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void Init_ElementRhs_QT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void Init_ElementRhs_QTT (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void Init_ElementRhs_QQ (const EDC< DH, VECTOR, dealdim > &, dealii::Vector< double > &, double)
 
virtual void Init_ElementRhs (const dealii::Function< dealdim > *init_values, const EDC< DH, VECTOR, dealdim > &edc, dealii::Vector< double > &local_vector, double scale)
 
virtual void Init_ElementMatrix (const EDC< DH, VECTOR, dealdim > &edc, dealii::FullMatrix< double > &local_entry_matrix, double scale, double)
 
virtual dealii::UpdateFlags GetUpdateFlags () const
 
virtual dealii::UpdateFlags GetFaceUpdateFlags () const
 
virtual bool HasFaces () const
 
virtual bool HasInterfaces () const
 
void SetProblemType (std::string type)
 
virtual unsigned int GetControlNBlocks () const
 
virtual unsigned int GetStateNBlocks () const
 
virtual std::vector< unsigned
int > & 
GetControlBlockComponent ()
 
virtual const std::vector
< unsigned int > & 
GetControlBlockComponent () const
 
virtual std::vector< unsigned
int > & 
GetStateBlockComponent ()
 
virtual const std::vector
< unsigned int > & 
GetStateBlockComponent () const
 
virtual void SetTime (double) const
 
unsigned int GetStateNComponents () const
 
template<typename ELEMENTITERATOR >
bool AtInterface (ELEMENTITERATOR &element, unsigned int face) const
 

Data Fields

boost::function1< void, double & > ResidualModifier
 
boost::function1< void,
dealii::Vector< double > & > 
VectorResidualModifier
 

Protected Attributes

std::string problem_type_
 

Detailed Description

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR, int dealdim>
class DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >

A template providing all evaluations of a PDE that may be used during the solution of a PDE or an optimization with a PDE constraint.

Whenever used below: u denotes the solution to the PDE q denotes a given control (i.e., parameter to the PDE) z denotes the adjoint solution dq denotes a given, fixed, direction in the control space du denotes the tangent solution to the PDE according to a given dq dz denotes an auxilliary adjoint used to caclulate second derivatives, e.g., the hessian.

and denote the basis functions in the state and control test space

Constructor & Destructor Documentation

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::PDEInterface ( )
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::~PDEInterface ( )
virtual

Member Function Documentation

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
template<typename ELEMENTITERATOR >
bool DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::AtInterface ( ELEMENTITERATOR &  element,
unsigned int  face 
) const
inline

Given a vector of active element iterators and a facenumber, checks if the face belongs to an 'interface' (i.e. the adjoining elements have different material ids).

Can be changed by the user to allow the use of dg-methods (i.e., then this function needs to return true for all input arguments.

ELEMENTITERATOR Class of the elementiterator.

Parameters
elementThe element in question.
faceLocal number of the face for which we ask if it is at the interface.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_Q ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_QQ ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_QT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_QTT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_QU ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_U ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_UQ ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_UT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_UTT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryEquation_UU ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryMatrix ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryMatrix_T ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::BoundaryRightHandSide ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ControlBoundaryEquation ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

This implements the scalar product on the boundary in the control space. The equation is used to calculate the representation of the cost functional gradient given the derivative of the cost functional

Parameters
fdcThe FaceDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ControlBoundaryMatrix ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double   
)
virtual

This implements the matrix corresponding to ControlBoundaryEquation

Parameters
fdcThe FaceDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_matrixThe matrix containing the integrals ordered according to the local number of the testfunction.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ControlElementEquation ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

This implements the scalar product in the control space. The equation is used to calculate the representation of the cost functional gradient given the derivative of the cost functional

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ControlElementMatrix ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double   
)
virtual

This implements the matrix corresponding to ControlElementEquation

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_matrixThe matrix containing the integrals ordered according to the local number of the testfunction.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

Assuming that the PDE is given in the form a(u;) = f(), this function implements all terms in a(u;) that are represented by intergrals over elements T. Hence, if a(u;) = a_T(u;) + ... then this function needs to implement a_T(u;) a_T may depend upon any spatial derivatives, but not on temporal derivatives.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_Q ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

Assuming that the element equation a_T depends not only on the state u, but also on a control q, this term implements the derivative of ElementEquation with respect to variations in q. This term implements the derivative of ElementEquation with respect to the u argument. I.e., if ElementEquation implements the term int_T a_T(u,q;) then this method implements a_T'_q(u,q;,z) .

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_QQ ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

Analog to ElementEquation_Q, but now considering second derivatives with respect to q, i.e., we calculate a_T''_{qq}(u,q;dq,,z) where dq is the given direction.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_QT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double  scale_ico 
)
virtual

Analog to ElementEqution_Q this term implements the derivative of ElementEquation with respect to the control argument. In contrast to ElementEqution_Q the test function is taken from the state space, while the argument for the control variation dq is fixed. ElementEquation implements the term int_T a_T(u,q;) then this method implements a_T'_q(u,q;dq,) .

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_QTT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

Analog to ElementEqution_Q, the only difference is that the argument z is exchanged by dz int_T a_T(u,q;) then this method implements a_T'_q(u,q;,dz) .

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_QU ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

Analog to ElementEquation_U and ElementEquation_Q, but now considering the mixed second derivatives with respect to u and q, i.e., we calculate a_T''_{qu}(u,q;dq,,z) where dq is a given variation for the control. This means the test function is taken in the state space.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_U ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

This term implements the derivative of ElementEquation with respect to the u argument. I.e., if ElementEquation implements the term int_T a_T(u;) then this method implements a_T'(u;,z) where denotes the direction to which the derivative is applied

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_UQ ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

Analog to ElementEquation_QU, but with different arguments, i.e., we calculate a_T''_{uq}(u,q;du,phi_q,z) where du is the given tangent direction. This means the test function is taken in the control space.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_UT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

This term implements the derivative of ElementEquation with respect to the u argument. I.e., if ElementEquation implements the term int_T a_T(u;) then this method implements a_T'(u;du,) . In contrast to ElementEquation_U the arguments are exchanged.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_UTT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

This term implements the derivative of ElementEquation with respect to the u argument. I.e., if ElementEquation implements the term int_T a_T(u;) then this method implements a_T'(u;phi,dz) .

This implements the same form as ElementEquation_U, but with exchanged functions, i.e., dz instead of z.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementEquation_UU ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

Analog to ElementEquation_U, but now considering second derivatives with respect to u, i.e., we calculate a_T''_{uu}(u,q;du,,z) where du is the given tangent direction

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementMatrix ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual

This implements the element integral used to calculate the stiffness matrix for the primal PDE. It corresponds to the evaluation of the matrix entries a_ij = a_T'(u;,)

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_matrixThe matrix containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementMatrix_T ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual

This implements the element integral used to calculate the stiffness matrix for the adjoint PDE. It corresponds to the evaluation of the matrix entries a_ji = a_T'(u;,)

By default, this function calls Element_Matrix and afterwards returns the transposed of the matrix.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_vectorThe matrix containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
scale_icoA special scaling parameter to be used in certain parts of the equation if they need to be treated differently in time-stepping schemes, see the PDF-documentation for more details.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementRightHandSide ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Implements the element integral corresponding to given volume data for the PDE.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquation ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Assuming that the discretization of temporal derivatives by a backward difference, i.e., u(t_i) 1/ t ( u(t_i) - u(t_{i-1}) in several cases the the temporal derivatives of the equation give rise to a spacial integral of the form

T(u(t_i); (t_i)) - T(u(t_{i-1}); (t_{i-1}))

This equation is used to implement the element contribution T(u,)

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquation_U ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Same as ElementTimeEquation, but here the derivative of T with respect to u is considered. Here, the derivative of T in u in a direction du for a fixed test function z is denoted as T'(u;du,z)

This equation is used to implement the element contribution T'(u;,z)

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquation_UT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Same as ElementTimeEquation_U, but we exchange the argument for the test function.

This equation is used to implement the element contribution T'(u;du,)

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquation_UTT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Same as ElementTimeEquation_UT, but we exchange the argument for the test function.

This equation is used to implement the element contribution T'(u;,dz)

Note that this is the same function as in ElementTimeEquation_U, but it is used with an other argument dz instead of z.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquationExplicit ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

In certain cases, the assumption of ElementTimeEquation are not meet, i.e., we can not use the same operator T at t_i and t_{i-1}. In these cases instead of ElementTimeEquation this funtion is used, where the user can implement the complete approximation of the temporal derivative.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquationExplicit_U ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Analog to ElementTimeEquationExplicit, this function is used to replace ElementTimeEquation_U if needed.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquationExplicit_UT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Analog to ElementTimeEquationExplicit, this function is used to replace ElementTimeEquation_UT if needed.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquationExplicit_UTT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

Analog to ElementTimeEquationExplicit, this function is used to replace ElementTimeEquation_UTT if needed.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeEquationExplicit_UU ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual

For nonlinear terms involved in the temporal derivative the second derivatives with respect to the state of the time derivative are implemented here.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_vectorThe vector containing the integrals ordered according to the local number of the testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeMatrix ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &   
)
virtual

This implements the element integral used to calculate the matrix for the primal PDE corresponding to the time derivatives given in ElementTimeEquation.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_matrixThe matrix containing the integrals ordered according to the local number of the testfunction.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeMatrix_T ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &   
)
virtual

The transposed of ElementTimeEquation.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_matrixThe matrix containing the integrals ordered according to the local number of the testfunction.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeMatrixExplicit ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &   
)
virtual

This implements the element integral used to calculate the matrix for the primal PDE corresponding to the time derivatives given in ElementTimeEquationExplicit.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_matrixThe matrix containing the integrals ordered according to the local number of the testfunction.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ElementTimeMatrixExplicit_T ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &   
)
virtual

The transposed of ElementTimeEquationExplicit.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
local_entry_matrixThe matrix containing the integrals ordered according to the local number of the testfunction.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual

The following Face... and Boundary... methods implement the analog terms as the corresponding Element... methods, except that now integrals on faces between elements or on the domain boundary are considered.

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_Q ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_QQ ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_QT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_QTT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_QU ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_U ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_UQ ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_UT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_UTT ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceEquation_UU ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceMatrix ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual

Documentation in optproblemcontainer.h.

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceMatrix_T ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::FaceRightHandSide ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
std::vector< unsigned int > & DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetControlBlockComponent ( )
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
const std::vector< unsigned int > & DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetControlBlockComponent ( ) const
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
unsigned int DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetControlNBlocks ( ) const
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
UpdateFlags DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetFaceUpdateFlags ( ) const
virtual

Returns the update flags needed by the integrator to decide which finite element informations need to be calculated on the next face (including on boundaries).

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
std::vector< unsigned int > & DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetStateBlockComponent ( )
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
const std::vector< unsigned int > & DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetStateBlockComponent ( ) const
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
unsigned int DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetStateNBlocks ( ) const
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
unsigned int DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetStateNComponents ( ) const
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
UpdateFlags DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::GetUpdateFlags ( ) const
virtual

Returns the update flags needed by the integrator to decide which finite element informations need to be calculated on the next element.

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
bool DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::HasFaces ( ) const
virtual

Should return true, if integration on interior faces is required; i.e., in dG implementations.

The default is false.

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
bool DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::HasInterfaces ( ) const
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::Init_ElementEquation ( const EDC< DH, VECTOR, dealdim > &  edc,
dealii::Vector< double > &  local_vector,
double  scale,
double   
)
inlinevirtual

The following functions Init_... implement the equation used for transfering the given initial values to the finite element space.

The initial data may depend on the control, but (obviously) not on the state itself, hence derivatives with respect to the control are considered.

Default is componentwise L2 projection.

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::Init_ElementMatrix ( const EDC< DH, VECTOR, dealdim > &  edc,
dealii::FullMatrix< double > &  local_entry_matrix,
double  scale,
double   
)
inlinevirtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::Init_ElementRhs ( const dealii::Function< dealdim > *  init_values,
const EDC< DH, VECTOR, dealdim > &  edc,
dealii::Vector< double > &  local_vector,
double  scale 
)
inlinevirtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::Init_ElementRhs_Q ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
inlinevirtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::Init_ElementRhs_QQ ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
inlinevirtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::Init_ElementRhs_QT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
inlinevirtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::Init_ElementRhs_QTT ( const EDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double   
)
inlinevirtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::InterfaceEquation ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::InterfaceEquation_U ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::Vector< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::InterfaceMatrix ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::InterfaceMatrix_T ( const FDC< DH, VECTOR, dealdim > &  ,
dealii::FullMatrix< double > &  ,
double  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::SetProblemType ( std::string  type)
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
virtual void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::SetTime ( double  ) const
inlinevirtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongBoundaryResidual ( const FDC< DH, VECTOR, dealdim > &  ,
const FDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongBoundaryResidual_Control ( const FDC< DH, VECTOR, dealdim > &  ,
const FDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual

Similar to the StongElementResidual, this function implements the strong boundary residual (i.e., jumps in conormal direction) for the gradient equation, i.e., j'(q) = 0, if present.

Parameters
fdcThe FaceDataContainer object which provides access to all information on the face, e.g., test-functions, mesh size,...
fdc_weightThe FaceDataContainer for the weight-function, e.g., the testfunction by which the residual needs to be multiplied
retThe value of the integral on the element of residual times testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongBoundaryResidual_U ( const FDC< DH, VECTOR, dealdim > &  ,
const FDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongElementResidual ( const EDC< DH, VECTOR, dealdim > &  ,
const EDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual

This function is used for error estimation and should implement the strong form of the residual on an element T.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
edc_wightThe ElementDataContainer for the weight-function, e.g., the testfunction by which the residual needs to be multiplied
retThe value of the integral on the element of residual times testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongElementResidual_Control ( const EDC< DH, VECTOR, dealdim > &  ,
const EDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual

Similar to the StongElementResidual, this function implements the strong element residual for the gradient equation, i.e., j'(q) = 0.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
edc_weightThe ElementDataContainer for the weight-function, e.g., the testfunction by which the residual needs to be multiplied
retThe value of the integral on the element of residual times testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongElementResidual_U ( const EDC< DH, VECTOR, dealdim > &  ,
const EDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual

Similar to the StongElementResidual, this function implements the strong element residual for the adjoint equation.

Parameters
edcThe ElementDataContainer object which provides access to all information on the element, e.g., test-functions, mesh size,...
edc_weightThe ElementDataContainer for the weight-function, e.g., the testfunction by which the residual needs to be multiplied
retThe value of the integral on the element of residual times testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongFaceResidual ( const FDC< DH, VECTOR, dealdim > &  ,
const FDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongFaceResidual_Control ( const FDC< DH, VECTOR, dealdim > &  ,
const FDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual

Similar to the StongElementResidual, this function implements the strong face residual (i.e., jumps in conormal direction) for the gradient equation, i.e., j'(q) = 0, if present.

Parameters
fdcThe FaceDataContainer object which provides access to all information on the face, e.g., test-functions, mesh size,...
fdc_weightThe FaceDataContainer for the weight-function, e.g., the testfunction by which the residual needs to be multiplied
retThe value of the integral on the element of residual times testfunction.
scaleA scaling parameter to be used in all equations.
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
void DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::StrongFaceResidual_U ( const FDC< DH, VECTOR, dealdim > &  ,
const FDC< DH, VECTOR, dealdim > &  ,
double &  ,
double   
)
virtual

Field Documentation

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
std::string DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::problem_type_
protected
template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
boost::function1<void, double&> DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::ResidualModifier

This function is set by the error estimators in order to allow the calculation of squared norms of the residual as needed for Residual Error estimators as well as the residual itself as needed by the DWR estimators.

template<template< template< int, int > class DH, typename VECTOR, int dealdim > class EDC, template< template< int, int > class DH, typename VECTOR, int dealdim > class FDC, template< int, int > class DH, typename VECTOR , int dealdim>
boost::function1<void, dealii::Vector<double>&> DOpE::PDEInterface< EDC, FDC, DH, VECTOR, dealdim >::VectorResidualModifier

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