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

Integrator for $ (\int_\Omega s(x) w_i(x) da )_{i} $, of some scalar valued function $ s$. More...

#include <unityTriangleIntegratorShellFE.h>

Public Member Functions

 UnitTriangleFENonlinOpIntegratorShellFE (const ConfiguratorType &Conf)
 
virtual ~UnitTriangleFENonlinOpIntegratorShellFE ()
 
void assembleAdd (VectorType &Dest) const
 
RealType getNonlinearity (const typename ConfiguratorType::ElementType &El, int QuadPoint) const
 interface function, has to be provided in derived classes. More...
 

Protected Types

typedef ConfiguratorType::RealType RealType
 
typedef ConfiguratorType::VectorType VectorType
 
typedef ConfiguratorType::ElementType ElementType
 

Protected Member Functions

Imp & asImp ()
 
const Imp & asImp () const
 

Protected Attributes

const ConfiguratorType_config
 

Detailed Description

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

Integrator for $ (\int_\Omega s(x) w_i(x) da )_{i} $, of some scalar valued function $ s$.

Definition at line 16 of file unityTriangleIntegratorShellFE.h.

Member Typedef Documentation

template<typename ConfiguratorType , typename Imp >
typedef ConfiguratorType::ElementType shellFE::UnitTriangleFENonlinOpIntegratorShellFE< ConfiguratorType, Imp >::ElementType
protected

Definition at line 21 of file unityTriangleIntegratorShellFE.h.

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

Definition at line 19 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType , typename Imp >
typedef ConfiguratorType::VectorType shellFE::UnitTriangleFENonlinOpIntegratorShellFE< ConfiguratorType, Imp >::VectorType
protected

Definition at line 20 of file unityTriangleIntegratorShellFE.h.

Constructor & Destructor Documentation

template<typename ConfiguratorType , typename Imp >
shellFE::UnitTriangleFENonlinOpIntegratorShellFE< ConfiguratorType, Imp >::UnitTriangleFENonlinOpIntegratorShellFE ( const ConfiguratorType Conf)
inline

Definition at line 25 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType , typename Imp >
virtual shellFE::UnitTriangleFENonlinOpIntegratorShellFE< ConfiguratorType, Imp >::~UnitTriangleFENonlinOpIntegratorShellFE ( )
inlinevirtual

Definition at line 27 of file unityTriangleIntegratorShellFE.h.

27 {}

Member Function Documentation

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

Definition at line 64 of file unityTriangleIntegratorShellFE.h.

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

Definition at line 65 of file unityTriangleIntegratorShellFE.h.

65 { return static_cast<const Imp&> ( *this ); }
template<typename ConfiguratorType , typename Imp >
void shellFE::UnitTriangleFENonlinOpIntegratorShellFE< ConfiguratorType, Imp >::assembleAdd ( VectorType Dest) const
inline

Definition at line 29 of file unityTriangleIntegratorShellFE.h.

29  {
30  RealType *nl_cache = new RealType[ _config.maxNumQuadPoints() ];
31 
32  for ( int elementIdx = 0; elementIdx < _config.getInitializer().getNumTriangs(); ++elementIdx){
33  const ElementType& El ( _config.getInitializer().getTriang( elementIdx ) );
34  const int numLocalDofs = _config.getNumLocalDofs ( El );
35 
37  const int numQuadPoints = bfs.numQuadPoints( );
38 
39  for ( int q = 0; q < numQuadPoints; ++q )
40  nl_cache[q] = this->asImp().getNonlinearity ( El, q );
41 
42  RealType aux;
43 
44  for ( int dof = 0; dof < numLocalDofs; ++dof ) {
45  aux = 0.;
46 
47  for ( int q = 0; q < numQuadPoints; ++q )
48  aux += nl_cache[q] * bfs.evaluate ( dof, q ) * bfs.getWeight ( q );
49 
50  Dest[ _config.localToGlobal ( El, dof ) ] += 0.5 * aux;
51  }
52  }
53  delete[] nl_cache;
54  }
const TriangleType & getTriang(const int num) const
Definition: triangMesh.h:87
RealType evaluate(int BaseFuncNum, const DomVecType &RefCoord) const
Evaluates a base function at.
int numQuadPoints() const
Returns the number of quadrature points.
int maxNumQuadPoints() const
Returns the maximum number of quadrature points per element.
RealType getWeight(int QuadPoint) const
Returns the quadrature weight at QuadPointa quadrature point.
int getNumTriangs() const
Definition: triangMesh.h:46
int getNumLocalDofs(const ElementType &) const
Get The number of local degrees of freedom on an element.
double RealType
Definition: ex2.cpp:21
int localToGlobal(const ElementType &T, int localIndex) const
Returns global index of the dof with number localIndex.
const InitType & getInitializer() const
Returns the mesh.
const BaseFuncSetType & getBaseFunctionSet(const ElementType &) const
Returns the base funcitons set for an element.
template<typename ConfiguratorType , typename Imp >
RealType shellFE::UnitTriangleFENonlinOpIntegratorShellFE< ConfiguratorType, Imp >::getNonlinearity ( const typename ConfiguratorType::ElementType El,
int  QuadPoint 
) const
inline

interface function, has to be provided in derived classes.

Definition at line 57 of file unityTriangleIntegratorShellFE.h.

57  {
58  throw std::invalid_argument( aol::strprintf ( "Called the interface function. In File %s at line %d.", __FILE__, __LINE__ ).c_str() );
59  return this->asImp().getNonlinearity ( El, QuadPoint );
60  }
string strprintf(const char *format,...)
Give back formatted string, analogously to sprintf, but save the long way &#39;round with char arrays...
Definition: aol.cpp:5

Member Data Documentation

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

Definition at line 22 of file unityTriangleIntegratorShellFE.h.


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