Programming tasks to Scientific Computing I
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp > Class Template Reference

Provides an easy interface to Finite Element operators of the form $ \mbox{div}(A(x)\nabla u)$, where $A$ is an ASYMMETRIC coefficient matrix. The corresponding matrix assembly yields $ \left(\int_\Omega \nabla\phi_i\cdot A(x)\nabla\phi_j dx\right)_{ij} $ for FE basis functions $ \phi_i,\phi_j $. More...

#include <unityTriangleIntegratorShellFE.h>

Public Member Functions

 UnitTriangleFELinWeightedStiffIntegrator (const ConfiguratorType &Config)
 
void getCoeffMatrix (const typename ConfiguratorType::ElementType &El, int QuadPoint, Matrix22 &Matrix) const
 This function has to be provided in the implementation (derived class) of the interface. More...
 
void prepareLocalMatrix (const typename ConfiguratorType::ElementType &El, LocalMatrixType &LocalMatrix) const
 Computes the numerical quadrature of the bilinear form and saves the values locally. More...
 
- Public Member Functions inherited from shellFE::MatrixValuedIntegratorBase< ConfiguratorType, UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp > >
 MatrixValuedIntegratorBase (const ConfiguratorType &conf, const RealType factor=1.0)
 
void assembleTripletListCached () const
 
void assemble (SparseMatrixType &Dest) const
 
void assembleDirichlet (SparseMatrixType &Dest, const MaskType &boundaryMask) const
 

Protected Types

typedef ConfiguratorType::RealType RealType
 
typedef ConfiguratorType::DomVecType DomVecType
 
typedef ConfiguratorType::Matrix22 Matrix22
 
typedef ConfiguratorType::LocalMatrixType LocalMatrixType
 

Protected Member Functions

Imp & asImp ()
 
const Imp & asImp () const
 
- Protected Member Functions inherited from shellFE::MatrixValuedIntegratorBase< ConfiguratorType, UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp > >
void assembleTripletList (std::vector< TripletType > &tripletList, const RealType Factor) const
 
UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp > & asImp ()
 
const UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp > & asImp () const
 

Protected Attributes

const ConfiguratorType_config
 
- Protected Attributes inherited from shellFE::MatrixValuedIntegratorBase< ConfiguratorType, UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp > >
const ConfiguratorType_config
 
std::vector< TripletType_tripletList
 
RealType _factor
 
bool _tripletListAssembled
 

Additional Inherited Members

- Public Types inherited from shellFE::MatrixValuedIntegratorBase< ConfiguratorType, UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp > >
typedef ConfiguratorType::RealType RealType
 
typedef ConfiguratorType::LocalMatrixType LocalMatrixType
 
typedef ConfiguratorType::TripletType TripletType
 
typedef ConfiguratorType::ElementType ElementType
 
typedef ConfiguratorType::MaskType MaskType
 

Detailed Description

template<typename ConfiguratorType, typename Imp>
class shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >

Provides an easy interface to Finite Element operators of the form $ \mbox{div}(A(x)\nabla u)$, where $A$ is an ASYMMETRIC coefficient matrix. The corresponding matrix assembly yields $ \left(\int_\Omega \nabla\phi_i\cdot A(x)\nabla\phi_j dx\right)_{ij} $ for FE basis functions $ \phi_i,\phi_j $.

Definition at line 212 of file unityTriangleIntegratorShellFE.h.

Member Typedef Documentation

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::DomVecType shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::DomVecType
protected

Definition at line 216 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::LocalMatrixType shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::LocalMatrixType
protected

Definition at line 218 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::Matrix22 shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::Matrix22
protected

Definition at line 217 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::RealType shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::RealType
protected

Definition at line 215 of file unityTriangleIntegratorShellFE.h.

Constructor & Destructor Documentation

template<typename ConfiguratorType, typename Imp>
shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::UnitTriangleFELinWeightedStiffIntegrator ( const ConfiguratorType Config)
inline

Definition at line 222 of file unityTriangleIntegratorShellFE.h.

222  :
223  MatrixValuedIntegratorBase< ConfiguratorType, UnitTriangleFELinWeightedStiffIntegrator<ConfiguratorType, Imp> > ( Config ),
224  _config ( Config ) {}

Member Function Documentation

template<typename ConfiguratorType, typename Imp>
Imp& shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::asImp ( )
inlineprotected

Definition at line 267 of file unityTriangleIntegratorShellFE.h.

267 { return static_cast<Imp&> ( *this ); }
template<typename ConfiguratorType, typename Imp>
const Imp& shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::asImp ( ) const
inlineprotected

Definition at line 268 of file unityTriangleIntegratorShellFE.h.

268 { return static_cast<const Imp&> ( *this ); }
template<typename ConfiguratorType, typename Imp>
void shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::getCoeffMatrix ( const typename ConfiguratorType::ElementType El,
int  QuadPoint,
Matrix22 Matrix 
) const
inline

This function has to be provided in the implementation (derived class) of the interface.

Definition at line 228 of file unityTriangleIntegratorShellFE.h.

229  {
230  this->asImp().getCoeffMatrix ( El, QuadPoint, Matrix );
231  }
template<typename ConfiguratorType, typename Imp>
void shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::prepareLocalMatrix ( const typename ConfiguratorType::ElementType El,
LocalMatrixType LocalMatrix 
) const
inline

Computes the numerical quadrature of the bilinear form and saves the values locally.

Definition at line 234 of file unityTriangleIntegratorShellFE.h.

235  {
236 
237 
238  const int numDofs = _config.getNumLocalDofs ( El );
239 
240  for ( int i = 0; i < numDofs; ++i )
241  for ( int j = 0; j < numDofs; ++j )
242  LocalMatrix(i,j) = 0.;
243 
244  Matrix22 mat;
245  DomVecType matgrad1;
246 
247  const typename ConfiguratorType::BaseFuncSetType &bfs = _config.getBaseFunctionSet ( El );
248  const int numQuadPoints = bfs.numQuadPoints( );
249 
250  for ( int q = 0; q < numQuadPoints; ++q ) {
251  getCoeffMatrix ( El, q, mat );
252  for ( int i = 0; i < numDofs; ++i ) {
253  DomVecType grad_i = bfs.evaluateGradient( i, q );
254  matgrad1 = mat * grad_i;
255  for ( int j = 0; j < numDofs; ++j ) {
256  DomVecType grad_j = bfs.evaluateGradient( j, q );
257  LocalMatrix(j,i) += matgrad1.dot( grad_j ) * bfs.getWeight ( q );
258  }
259  }
260  }
261  for ( int i = 0; i < numDofs; ++i )
262  for ( int j = 0; j < numDofs; ++j )
263  LocalMatrix(i,j) *= 0.5; // 0.5 is the volume of the unit triangle
264  }
int numQuadPoints() const
Returns the number of quadrature points.
RealType getWeight(int QuadPoint) const
Returns the quadrature weight at QuadPointa quadrature point.
void evaluateGradient(int BaseFuncNum, const DomVecType &, DomVecType &Gradient) const
Evaluates the gradient of a base function at RefCoord.
int getNumLocalDofs(const ElementType &) const
Get The number of local degrees of freedom on an element.
void getCoeffMatrix(const typename ConfiguratorType::ElementType &El, int QuadPoint, Matrix22 &Matrix) const
This function has to be provided in the implementation (derived class) of the interface.
const BaseFuncSetType & getBaseFunctionSet(const ElementType &) const
Returns the base funcitons set for an element.

Member Data Documentation

template<typename ConfiguratorType, typename Imp>
const ConfiguratorType& shellFE::UnitTriangleFELinWeightedStiffIntegrator< ConfiguratorType, Imp >::_config
protected

Definition at line 219 of file unityTriangleIntegratorShellFE.h.


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