LifeV
FESpace< MeshType, MapType > Class Template Reference

FESpace - Short description here please! More...

#include <FESpace.hpp>

+ Collaboration diagram for FESpace< MeshType, MapType >:

Protected Types

enum  spaceData {
  P0 = 1, P1 = 10, P1_HIGH = 15, P1Bubble = 16,
  P2 = 20, P2_HIGH = 25, P2Bubble = 26
}
 Set space Map (useful for switch syntax with strings) More...
 

Protected Attributes

std::map< std::string, spaceDataM_spaceMap
 
meshPtr_Type M_mesh
 reference to the mesh More...
 
const ReferenceFEM_refFE
 Reference FE for the velocity. More...
 
const QuadratureRuleM_Qr
 Quadrature rule for volumic elementary computations. More...
 
const QuadratureRuleM_bdQr
 Quadrature rule for surface elementary computations. More...
 
UInt M_fieldDim
 dimension of the field variable ( scalar/vector field) More...
 
std::shared_ptr< DOFM_dof
 A shared pointer to the DOF object. More...
 
UInt M_dim
 The number of total dofs. More...
 
std::shared_ptr< CurrentFEM_fe
 Current FE. More...
 
std::shared_ptr< CurrentFEManifoldM_feBd
 
mapPtr_Type M_map
 Map. More...
 

Public Types

typedef std::function< Real(Real const &, Real const &, Real const &, Real const &, ID const &) > function_Type
 
typedef MeshType mesh_Type
 
typedef std::shared_ptr< mesh_TypemeshPtr_Type
 
typedef MapType map_Type
 
typedef std::shared_ptr< map_TypemapPtr_Type
 
typedef map_Type::commPtr_Type commPtr_Type
 

Constructor & Destructor

 LIFEV_DEPRECATED (FESpace(MeshPartitioner< MeshType > &mesh, const ReferenceFE &refFE, const QuadratureRule &Qr, const QuadratureRule &bdQr, const Int fDim, const commPtr_Type &commptr))
 
 LIFEV_DEPRECATED (FESpace(MeshPartitioner< MeshType > &mesh, const std::string &space, const Int fDim, const commPtr_Type &commptr))
 
 FESpace (meshPtr_Type mesh, const ReferenceFE &refFE, const QuadratureRule &Qr, const QuadratureRule &bdQr, const Int fDim, const commPtr_Type &commptr)
 
 FESpace (meshPtr_Type mesh, const std::string &space, const Int fDim, const commPtr_Type &commptr)
 
virtual ~FESpace ()
 Do nothing destructor. More...
 

Methods

template<typename vector_type >
void interpolate (const function_Type &fct, vector_type &vect, const Real time=0.)
 Interpolate a given function nodally onto a vector. More...
 
template<typename vector_type >
void interpolateBC (BCHandler &BCh, vector_type &vect, const Real time)
 
template<typename ReturnType , typename vector_type >
void interpolate (const FEFunction< MeshType, MapType, ReturnType > *fEFunction, vector_type &vector, const Real time=0.)
 Interpolation method for FEFunctions. More...
 
template<typename vector_type >
void l2ScalarProduct (const function_Type &fct, vector_type &vec, const Real t)
 calculate L2 velocity error for given exact velocity function More...
 
template<typename vector_type >
Real l20Error (const function_Type &fexact, const vector_type &vec, const Real time, Real *relError=0)
 
template<typename vector_type >
Real l2Error (const function_Type &fexact, const vector_type &vec, const Real time, Real *relError=0)
 
template<typename function , typename vector_type >
Real h1Error (const function &fexact, const vector_type &vec, const Real time, Real *relError=0)
 
template<typename vector_type >
Real l2Norm (const vector_type &vec)
 
template<typename function >
Real l2NormFunction (const function &f, const Real time=0)
 
template<typename vector_type >
Real h1Norm (const vector_type &vec)
 
template<typename vector_type >
Real l2ErrorWeighted (const function_Type &exactSolution, const vector_type &solution, const function_Type &weight, const Real time)
 Method to computes the L2 error when using a weight function. More...
 
template<typename point_type , typename vector_type >
Real feInterpolateValue (const ID &elementID, const vector_type &solutionVector, const point_type &pt, const UInt &component=0) const
 This method computes the interpolate value of a given FE function in a given point. More...
 
template<typename point_type , typename vector_type >
Real feInterpolateValueLocal (const ID &elementID, const vector_type &solutionVector, const point_type &pt) const
 This method computes the interpolate value of a given FE function in a given point. More...
 
template<typename point_type , typename vector_type >
Real feInterpolateGradient (const ID &elementID, const vector_type &solutionVector, const point_type &pt, const UInt &gradientElement, const UInt &component=0) const
 This method computes the interpolated gradient of a given FE function in a given point. More...
 
template<typename point_type , typename vector_type >
Real feInterpolateGradientLocal (const ID &elementID, const vector_type &solutionVector, const point_type &pt, const UInt &gradientElement) const
 This method computes the interpolated gradient of a given FE function in a given point. More...
 
template<typename vector_type >
vector_type feToFEInterpolate (const FESpace< mesh_Type, map_Type > &originalSpace, const vector_type &originalVector, const MapEpetraType &outputMapType=Unique) const
 This method enables to pass a solution from one FESpace to the present one. More...
 
template<typename vector_type >
vector_type gradientRecovery (const vector_type &solution, const UInt &component) const
 This method reconstruct a gradient of a solution in the present FE space. More...
 
template<typename vector_type >
vector_type recoveryFunction (const vector_type &solution) const
 This method reconstruct a gradient of a solution in the present FE space. More...
 
template<typename vector_type >
vector_type laplacianRecovery (const vector_type &solution) const
 Reconstruction of the laplacian using gradientRecovery procedures. More...
 
UInt polynomialDegree () const
 Return the polynomial degree of the finite element used. More...
 

Set Methods

void setQuadRule (const QuadratureRule &Qr)
 Method to set replace the quadrule. More...
 
void setBdQuadRule (const QuadratureRule &bdQr)
 

Get Methods

const meshPtr_Typemesh () const
 
meshPtr_Typemesh ()
 
const map_Typemap () const
 Returns map. More...
 
map_Typemap ()
 
const mapPtr_TypemapPtr () const
 
const DOFdof () const
 Returns the velocity dof. More...
 
DOFdof ()
 
const std::shared_ptr< DOF > & dofPtr () const
 
const CurrentFEfe () const
 Returns the current FE. More...
 
CurrentFEfe ()
 
CurrentFEManifoldfeBd ()
 Returns the current boundary FE. More...
 
const ReferenceFErefFE () const
 Returns the res FE. More...
 
const QuadratureRuleqr () const
 Returns the volumic quadratic rule. More...
 
const QuadratureRulebdQr () const
 Returns the surfasic quadratic rule. More...
 
const UIntdim () const
 Returns FE space dimension. More...
 
const UIntfieldDim () const
 

Private Methods

 FESpace (const FESpace &fespace)
 copy constructor More...
 
void createMap (const commPtr_Type &commptr)
 Creates the map for interprocessor communication. More...
 
void resetBoundaryFE ()
 Resets boundary data if necessary. More...
 
void setSpace (const std::string &space, UInt dimension)
 Set space. More...
 
template<typename vector_type >
vector_type interpolateGeneric (const FESpace< mesh_Type, map_Type > &OriginalSpace, const vector_type &OriginalVector) const
 This is a generic function called by feToFEInterpolate method. More...
 
template<typename vector_type >
vector_type linearInterpolate (const FESpace< mesh_Type, map_Type > &originalSpace, const vector_type &originalVector) const
 This is a specialized function called by feToFEInterpolate method. More...
 
template<typename vector_type >
vector_type P2Interpolate (const FESpace< mesh_Type, map_Type > &original_space, const vector_type &original_vector) const
 This is a specialized function called by feToFEInterpolate method. More...
 
template<typename vector_type >
vector_type RT0ToP0Interpolate (const FESpace< mesh_Type, map_Type > &original_space, const vector_type &original_vector) const
 This is a specialized function called by FESpace::feToFEInterpolate method for RT0 to P0 interpolation. More...
 

Detailed Description

template<typename MeshType, typename MapType>
class LifeV::FESpace< MeshType, MapType >

FESpace - Short description here please!

Author
Gilles Fourestey

Class representing the FE space, i.e. the reference FE and the geometric mapping.

Definition at line 78 of file FESpace.hpp.

Member Typedef Documentation

◆ function_Type

typedef std::function< Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > function_Type

Definition at line 87 of file FESpace.hpp.

◆ mesh_Type

Definition at line 88 of file FESpace.hpp.

◆ meshPtr_Type

typedef std::shared_ptr<mesh_Type> meshPtr_Type

Definition at line 89 of file FESpace.hpp.

◆ map_Type

typedef MapType map_Type

Definition at line 90 of file FESpace.hpp.

◆ mapPtr_Type

typedef std::shared_ptr<map_Type> mapPtr_Type

Definition at line 91 of file FESpace.hpp.

◆ commPtr_Type

typedef map_Type::commPtr_Type commPtr_Type

Definition at line 92 of file FESpace.hpp.

Member Enumeration Documentation

◆ spaceData

enum spaceData
protected

Set space Map (useful for switch syntax with strings)

Enumerator
P0 
P1 
P1_HIGH 
P1Bubble 
P2 
P2_HIGH 
P2Bubble 

Definition at line 559 of file FESpace.hpp.

Constructor & Destructor Documentation

◆ FESpace() [1/3]

FESpace ( meshPtr_Type  mesh,
const ReferenceFE refFE,
const QuadratureRule Qr,
const QuadratureRule bdQr,
const Int  fDim,
const commPtr_Type commptr 
)

Definition at line 656 of file FESpace.hpp.

◆ FESpace() [2/3]

FESpace ( meshPtr_Type  mesh,
const std::string &  space,
const Int  fDim,
const commPtr_Type commptr 
)

Definition at line 680 of file FESpace.hpp.

◆ ~FESpace()

virtual ~FESpace ( )
inlinevirtual

Do nothing destructor.

Definition at line 141 of file FESpace.hpp.

◆ FESpace() [3/3]

FESpace ( const FESpace< MeshType, MapType > &  fespace)
protected

copy constructor

Member Function Documentation

◆ LIFEV_DEPRECATED() [1/2]

LIFEV_DEPRECATED ( FESpace< MeshType, MapType >(MeshPartitioner< MeshType > &mesh, const ReferenceFE &refFE, const QuadratureRule &Qr, const QuadratureRule &bdQr, const Int fDim, const commPtr_Type &commptr)  )
Parameters
data_fileGetPot data file
refFE_ureference FE for the velocity
refFE_preference FE for the pressure
Qr_uelement quadrature rule for the velocity
bdQr_usurface quadrature rule for the velocity
Qr_pelement quadrature rule for the pressure
bdQr_psurface quadrature rule for the pressure
BCh_fluidboundary conditions for the fluid
ord_bdforder of the bdf time advancing scheme and incremental pressure approach (default: Backward Euler)

◆ LIFEV_DEPRECATED() [2/2]

LIFEV_DEPRECATED ( FESpace< MeshType, MapType >(MeshPartitioner< MeshType > &mesh, const std::string &space, const Int fDim, const commPtr_Type &commptr)  )

◆ interpolate() [1/2]

void interpolate ( const function_Type fct,
vector_type vect,
const Real  time = 0. 
)

Interpolate a given function nodally onto a vector.

Definition at line 720 of file FESpace.hpp.

◆ interpolateBC()

void interpolateBC ( BCHandler BCh,
vector_type vect,
const Real  time 
)

Definition at line 838 of file FESpace.hpp.

◆ interpolate() [2/2]

void interpolate ( const FEFunction< MeshType, MapType, ReturnType > *  fEFunction,
vector_type vector,
const Real  time = 0. 
)

Interpolation method for FEFunctions.

Parameters
fEFunctionPointer to an FEFunction
vectorInterpolated function
timeTime in the interpolation

Definition at line 777 of file FESpace.hpp.

◆ l2ScalarProduct()

void l2ScalarProduct ( const function_Type fct,
vector_type vec,
const Real  t 
)

calculate L2 velocity error for given exact velocity function

Parameters
pexactthe exact velocity as a function
timethe time
relErrorReal* to store the relative error in

Definition at line 891 of file FESpace.hpp.

◆ l20Error()

Real l20Error ( const function_Type fexact,
const vector_type vec,
const Real  time,
Real relError = 0 
)

Definition at line 932 of file FESpace.hpp.

◆ l2Error()

Real l2Error ( const function_Type fexact,
const vector_type vec,
const Real  time,
Real relError = 0 
)

Definition at line 986 of file FESpace.hpp.

◆ h1Error()

Real h1Error ( const function &  fexact,
const vector_type vec,
const Real  time,
Real relError = 0 
)

Definition at line 1130 of file FESpace.hpp.

◆ l2Norm()

Real l2Norm ( const vector_type vec)

Definition at line 1185 of file FESpace.hpp.

◆ l2NormFunction()

Real l2NormFunction ( const function &  f,
const Real  time = 0 
)

Definition at line 1039 of file FESpace.hpp.

◆ h1Norm()

Real h1Norm ( const vector_type vec)

Definition at line 1215 of file FESpace.hpp.

◆ l2ErrorWeighted()

Real l2ErrorWeighted ( const function_Type exactSolution,
const vector_type solution,
const function_Type weight,
const Real  time 
)

Method to computes the L2 error when using a weight function.

The scope of this method is to compute $ \left( \int w (u_{exact} - u_h) \right)^{1/2} $. The usual L2 error norm can be retrieved by using $ w=1 $

Definition at line 1067 of file FESpace.hpp.

◆ feInterpolateValue()

Real feInterpolateValue ( const ID elementID,
const vector_type solutionVector,
const point_type &  pt,
const UInt component = 0 
) const

This method computes the interpolate value of a given FE function in a given point.

The user of this function has to provide the element and the vector of the DOF values. Given these informations and a point P, the function compute the value: value = sum_{dof i} v(i) phi_i(P)

Note that the point P can be outside of the element considered (this is NOT checked).

Warning: this method has been only checked in 3D.

Parameters
elementIDThe ID of the element considered. The ID is the local ID (not the global one, see in the MeshEntity class).
solutionVectorThe vector containing the values in all nodes (not only of the considered element).
ptThe point where to perform the interpolation. Note that pt must allow for STL-type accessor [] and must have the size() method (typically a std::vector<Real>).
componentThe component for which the interpolation has to be performed (usefull only for the vectorial FE, set to 0 if the FE is scalar).
Returns
The interpolated value in the point.
Author
Samuel Quinodoz
Date
13-01-2010

Definition at line 1248 of file FESpace.hpp.

◆ feInterpolateValueLocal()

Real feInterpolateValueLocal ( const ID elementID,
const vector_type solutionVector,
const point_type &  pt 
) const

This method computes the interpolate value of a given FE function in a given point.

This method is the same as feInterpolateValue, but for the definition of the solution vector. Here, the solution is given only for the degrees of freedom of the given element. The parameter solutionVector is typically a std::vector<Real> containing the values in the dofs.

It is not possible to specify the component, the user has to provide directly the values for the wanted component in the solutionVector.

Author
Samuel Quinodoz
Date
13-01-2010

Definition at line 1311 of file FESpace.hpp.

◆ feInterpolateGradient()

Real feInterpolateGradient ( const ID elementID,
const vector_type solutionVector,
const point_type &  pt,
const UInt gradientElement,
const UInt component = 0 
) const

This method computes the interpolated gradient of a given FE function in a given point.

This method is the same as feInterpolateValue, but it computes the gradient insted of the value. Therefor, the desired element of the gradient has to be expressed using the parameter gradientElement.

For example, if gradientElement=i and component=j, then the results corresponds to the value of partial u_j / partial x_i.

Warning: This method has not been tested deeply, buggy behaviour is possible.

Author
Samuel Quinodoz
Date
13-01-2010

Definition at line 1361 of file FESpace.hpp.

◆ feInterpolateGradientLocal()

Real feInterpolateGradientLocal ( const ID elementID,
const vector_type solutionVector,
const point_type &  pt,
const UInt gradientElement 
) const

This method computes the interpolated gradient of a given FE function in a given point.

This method is the same as feInterpolateGradient, but it works with local values only, just as feInterpolateValueLocal.

Warning: This method has not been tested deeply, buggy behaviour is possible.

Author
Samuel Quinodoz
Date
13-01-2010

Definition at line 1431 of file FESpace.hpp.

◆ feToFEInterpolate()

vector_type feToFEInterpolate ( const FESpace< mesh_Type, map_Type > &  originalSpace,
const vector_type originalVector,
const MapEpetraType outputMapType = Unique 
) const

This method enables to pass a solution from one FESpace to the present one.

This method interpolates the values of a FE solution (given by the vector and the FESpace) in the dofs of the present FESpace: we note $ v_i $ the component of the vector given in argument and $ phi_i $ the basis functions of the FESpace given in argument. This defines a function on the domain given by: $ f(x) = \sum_i v_i \phi_i(x) $.

We search then the function $ g(x) $ belonging to the present FESpace such that $ g(x_j) = f(x_j) $ for all $ x_j $ dofs of the present FESpace. This function returns the vector $ w_j $ such that $ g(x) = sum_j w_j psi_j(x) $ where $ psi_j $ are the basis functions of the present space.

Warning: It is NOT true in general that $ w_j = f(x_j) $ (it is for lagrangian FEs). For example, if we want to pass the solution from the P_1 to the P_1 with bubbles, then the value associated to the bubble function will be always 0 even if the value of f is not 0 at the location of the dof of the bubble!

Parameters
originalSpaceThe space where the solution is defined originally.
originalVectorThe vector of the solution in the original space.
outputMapTypeThe map type (default: Unique) of the the returned vector.
Returns
The vector in the current FESpace having map type outputMapType.
Author
Samuel Quinodoz
Date
13-01-2010

Definition at line 1485 of file FESpace.hpp.

◆ gradientRecovery()

vector_type gradientRecovery ( const vector_type solution,
const UInt component 
) const

This method reconstruct a gradient of a solution in the present FE space.

The goal of this method is to build an approximation of the gradient of a given FE function in this FESpace. Typically, when one use P1 elements for approximating the solution of a given problem, the gradient is only piecewise constant. However, one could need continuous gradient. The solutions to this problem is either to use specific finite elements (like Hermite FE) or rely on a recovery procedure for the gradient.

This method implements a recovery procedure that performs a local average with weights corresponding to the areas of the elements:

\[ Gr(P) = \frac{\sum_{P \in T} |T| \nabla u(P)}{\sum_{P \in T} |T|} \]

See Zienkiewicz and Zhu (1987) for more details.

Results might be very wrong if you are not using lagrangian FE for tetrahedra

Definition at line 1572 of file FESpace.hpp.

◆ recoveryFunction()

vector_type recoveryFunction ( const vector_type solution) const

This method reconstruct a gradient of a solution in the present FE space.

The goal of this method is to build an approximation of the gradient of a given FE function in this FESpace. Typically, when one use P1 elements for approximating the solution of a given problem, the gradient is only piecewise constant. However, one could need continuous gradient. The solutions to this problem is either to use specific finite elements (like Hermite FE) or rely on a recovery procedure for the gradient.

This method implements a recovery procedure that performs a local average with weights corresponding to the areas of the elements:

\[ Gr(P) = \frac{\sum_{P \in T} |T| \nabla u(P)}{\sum_{P \in T} |T|} \]

See Zienkiewicz and Zhu (1987) for more details.

Results might be very wrong if you are not using lagrangian FE for tetrahedra

◆ laplacianRecovery()

vector_type laplacianRecovery ( const vector_type solution) const

Reconstruction of the laplacian using gradientRecovery procedures.

This method simply uses the FESpace::gradientRecovery method several times so that one can get a continuous approximation of the laplacian of the given solution.

Results might be very wrong if you are not using lagrangian FE for tetrahedra

Definition at line 1639 of file FESpace.hpp.

◆ polynomialDegree()

UInt polynomialDegree ( ) const

Return the polynomial degree of the finite element used.

Definition at line 2166 of file FESpace.hpp.

◆ setQuadRule()

void setQuadRule ( const QuadratureRule Qr)

Method to set replace the quadrule.

Parameters
QrThe new quadrule to be used in the FESpace

Definition at line 1814 of file FESpace.hpp.

◆ setBdQuadRule()

void setBdQuadRule ( const QuadratureRule bdQr)

Definition at line 1824 of file FESpace.hpp.

◆ mesh() [1/2]

const meshPtr_Type& mesh ( ) const
inline

Definition at line 401 of file FESpace.hpp.

◆ mesh() [2/2]

meshPtr_Type& mesh ( )
inline

Definition at line 405 of file FESpace.hpp.

◆ map() [1/2]

const map_Type& map ( ) const
inline

Returns map.

Definition at line 411 of file FESpace.hpp.

◆ map() [2/2]

map_Type& map ( )
inline

Definition at line 415 of file FESpace.hpp.

◆ mapPtr()

const mapPtr_Type& mapPtr ( ) const
inline

Definition at line 420 of file FESpace.hpp.

◆ dof() [1/2]

const DOF& dof ( ) const
inline

Returns the velocity dof.

Definition at line 426 of file FESpace.hpp.

◆ dof() [2/2]

DOF& dof ( )
inline

Definition at line 430 of file FESpace.hpp.

◆ dofPtr()

const std::shared_ptr<DOF>& dofPtr ( ) const
inline

Definition at line 434 of file FESpace.hpp.

◆ fe() [1/2]

const CurrentFE& fe ( ) const
inline

Returns the current FE.

Definition at line 440 of file FESpace.hpp.

◆ fe() [2/2]

CurrentFE& fe ( )
inline

Definition at line 444 of file FESpace.hpp.

◆ feBd()

CurrentFEManifold& feBd ( )
inline

Returns the current boundary FE.

Definition at line 450 of file FESpace.hpp.

◆ refFE()

const ReferenceFE& refFE ( ) const
inline

Returns the res FE.

Definition at line 456 of file FESpace.hpp.

◆ qr()

const QuadratureRule& qr ( ) const
inline

Returns the volumic quadratic rule.

Definition at line 462 of file FESpace.hpp.

◆ bdQr()

const QuadratureRule& bdQr ( ) const
inline

Returns the surfasic quadratic rule.

Definition at line 468 of file FESpace.hpp.

◆ dim()

const UInt& dim ( ) const
inline

Returns FE space dimension.

Definition at line 474 of file FESpace.hpp.

◆ fieldDim()

const UInt& fieldDim ( ) const
inline

Definition at line 478 of file FESpace.hpp.

◆ createMap()

void createMap ( const commPtr_Type commptr)
protected

Creates the map for interprocessor communication.

Definition at line 1851 of file FESpace.hpp.

◆ resetBoundaryFE()

void resetBoundaryFE ( )
protected

Resets boundary data if necessary.

Definition at line 1873 of file FESpace.hpp.

◆ setSpace()

void setSpace ( const std::string &  space,
UInt  dimension 
)
inlineprotected

Set space.

Definition at line 1665 of file FESpace.hpp.

◆ interpolateGeneric()

vector_type interpolateGeneric ( const FESpace< mesh_Type, map_Type > &  OriginalSpace,
const vector_type OriginalVector 
) const
protected

This is a generic function called by feToFEInterpolate method.

It allows to interpolate vectors between any two continuous and scalar FE spaces. It is used when other specialized functions are not provided

Parameters
originalSpaceThe space where the solution is defined originally.
originalVectorThe vector of the solution in the original space (must be a Repeated vector).
Returns
The Repeated vector in the current FESpace.
Author
Mauro Perego
Date
27-06-2011

Definition at line 1886 of file FESpace.hpp.

◆ linearInterpolate()

vector_type linearInterpolate ( const FESpace< mesh_Type, map_Type > &  originalSpace,
const vector_type originalVector 
) const
protected

This is a specialized function called by feToFEInterpolate method.

It allows to interpolate P1bubble, P2 vectors into P1 vectors, P1 into P1bubble vectors and Q2 vectors into Q1 vectors.

Parameters
originalSpaceThe space where the solution is defined originally.
originalVectorThe vector of the solution in the original space.
Returns
The Unique vector in the current FESpace.
Author
Mauro Perego
Date
27-06-2011

Definition at line 1961 of file FESpace.hpp.

◆ P2Interpolate()

vector_type P2Interpolate ( const FESpace< mesh_Type, map_Type > &  original_space,
const vector_type original_vector 
) const
protected

This is a specialized function called by feToFEInterpolate method.

It allows to interpolate P1, P1bubble vectors into P2 vectors, and Q1 vectors into Q2 vectors.

Parameters
originalSpaceThe space where the solution is defined originally.
originalVectorThe vector of the solution in the original space (must be a Repeated vector).
Returns
The Repeated vector in the current FESpace.
Author
Samuel Quinodoz
Date
13-01-2010

Definition at line 2000 of file FESpace.hpp.

◆ RT0ToP0Interpolate()

vector_type RT0ToP0Interpolate ( const FESpace< mesh_Type, map_Type > &  original_space,
const vector_type original_vector 
) const
protected

This is a specialized function called by FESpace::feToFEInterpolate method for RT0 to P0 interpolation.

Author
Alessio Fumagalli
Date
13-01-2010

Definition at line 2068 of file FESpace.hpp.

Field Documentation

◆ M_spaceMap

std::map<std::string, spaceData> M_spaceMap
protected

Definition at line 560 of file FESpace.hpp.

◆ M_mesh

meshPtr_Type M_mesh
protected

reference to the mesh

Definition at line 563 of file FESpace.hpp.

◆ M_refFE

const ReferenceFE* M_refFE
protected

Reference FE for the velocity.

Definition at line 566 of file FESpace.hpp.

◆ M_Qr

const QuadratureRule* M_Qr
protected

Quadrature rule for volumic elementary computations.

Definition at line 569 of file FESpace.hpp.

◆ M_bdQr

const QuadratureRule* M_bdQr
protected

Quadrature rule for surface elementary computations.

Definition at line 572 of file FESpace.hpp.

◆ M_fieldDim

UInt M_fieldDim
protected

dimension of the field variable ( scalar/vector field)

Definition at line 575 of file FESpace.hpp.

◆ M_dof

std::shared_ptr<DOF> M_dof
protected

A shared pointer to the DOF object.

Definition at line 578 of file FESpace.hpp.

◆ M_dim

UInt M_dim
protected

The number of total dofs.

Definition at line 581 of file FESpace.hpp.

◆ M_fe

std::shared_ptr<CurrentFE> M_fe
protected

Current FE.

Definition at line 584 of file FESpace.hpp.

◆ M_feBd

std::shared_ptr<CurrentFEManifold> M_feBd
protected

Definition at line 585 of file FESpace.hpp.

◆ M_map

mapPtr_Type M_map
protected

Map.

Definition at line 588 of file FESpace.hpp.


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