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

General interface for matrix valued integrators. More...

#include <unityTriangleIntegratorShellFE.h>

Public Types

typedef ConfiguratorType::RealType RealType
 
typedef ConfiguratorType::LocalMatrixType LocalMatrixType
 
typedef ConfiguratorType::TripletType TripletType
 
typedef ConfiguratorType::ElementType ElementType
 
typedef ConfiguratorType::MaskType MaskType
 

Public Member Functions

 MatrixValuedIntegratorBase (const ConfiguratorType &conf, const RealType factor=1.0)
 
void assembleTripletListCached () const
 
template<typename SparseMatrixType >
void assemble (SparseMatrixType &Dest) const
 
template<typename SparseMatrixType >
void assembleDirichlet (SparseMatrixType &Dest, const MaskType &boundaryMask) const
 

Protected Member Functions

void assembleTripletList (std::vector< TripletType > &tripletList, const RealType Factor) const
 
Imp & asImp ()
 
const Imp & asImp () const
 

Protected Attributes

const ConfiguratorType_config
 
std::vector< TripletType_tripletList
 
RealType _factor
 
bool _tripletListAssembled
 

Detailed Description

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

General interface for matrix valued integrators.

Definition at line 121 of file unityTriangleIntegratorShellFE.h.

Member Typedef Documentation

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::ElementType shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::ElementType

Definition at line 126 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::LocalMatrixType shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::LocalMatrixType

Definition at line 124 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::MaskType shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::MaskType

Definition at line 127 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::RealType shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::RealType

Definition at line 123 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
typedef ConfiguratorType::TripletType shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::TripletType

Definition at line 125 of file unityTriangleIntegratorShellFE.h.

Constructor & Destructor Documentation

template<typename ConfiguratorType, typename Imp>
shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::MatrixValuedIntegratorBase ( const ConfiguratorType conf,
const RealType  factor = 1.0 
)
inlineexplicit

Member Function Documentation

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

Definition at line 195 of file unityTriangleIntegratorShellFE.h.

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

Definition at line 196 of file unityTriangleIntegratorShellFE.h.

196 { return static_cast<const Imp&> ( *this ); }
template<typename ConfiguratorType, typename Imp>
template<typename SparseMatrixType >
void shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::assemble ( SparseMatrixType Dest) const
inline

Definition at line 165 of file unityTriangleIntegratorShellFE.h.

165  {
167  Dest.setFromTriplets( _tripletList.cbegin(), _tripletList.cend() );
168  }
template<typename ConfiguratorType, typename Imp>
template<typename SparseMatrixType >
void shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::assembleDirichlet ( SparseMatrixType Dest,
const MaskType boundaryMask 
) const
inline

Definition at line 171 of file unityTriangleIntegratorShellFE.h.

171  {
173 
174  std::vector<TripletType> tripletListMasked;
175  tripletListMasked.reserve(_config.getInitializer().getNumTriangs() * aol::Sqr( _config.getNumLocalDofs() ) );
176 
177  for( unsigned iter=0; iter < _tripletList.size(); ++iter ){
178  if( (boundaryMask[_tripletList[iter].row()]) || (boundaryMask[_tripletList[iter].col()]) ){
179  //Boundary node!
180  } else {
181  tripletListMasked.push_back( _tripletList[iter] );
182  }
183  }
184 
185  for ( int i = 0; i < _config.getNumGlobalDofs(); ++i ){
186  if ( boundaryMask[i] )
187  tripletListMasked.push_back( TripletType( i, i, 1.0 ) );
188  }
189 
190  Dest.setFromTriplets( tripletListMasked.begin(), tripletListMasked.end() );
191  }
T Sqr(const T a)
Definition: aol.h:134
int getNumTriangs() const
Definition: triangMesh.h:46
int getNumLocalDofs(const ElementType &) const
Get The number of local degrees of freedom on an element.
const InitType & getInitializer() const
Returns the mesh.
int getNumGlobalDofs() const
Returns the number of global degrees of freedom.
template<typename ConfiguratorType, typename Imp>
void shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::assembleTripletList ( std::vector< TripletType > &  tripletList,
const RealType  Factor 
) const
inlineprotected

Definition at line 134 of file unityTriangleIntegratorShellFE.h.

134  {
135  tripletList.reserve(_config.getInitializer().getNumTriangs() * aol::Sqr( _config.getNumLocalDofs() ) );
136  LocalMatrixType localMatrix;
137  int globalDofs[ ConfiguratorType::maxNumLocalDofs ];
138  for ( int elementIdx = 0; elementIdx < _config.getInitializer().getNumTriangs(); ++elementIdx){
139  const ElementType& El ( _config.getInitializer().getTriang( elementIdx ) );
140  this->asImp().prepareLocalMatrix ( El, localMatrix );
141  const int numLocalDofs = _config.getNumLocalDofs ( El );
142 
143  for ( int i = 0; i < numLocalDofs; ++i )
144  globalDofs[ i ] = _config.localToGlobal ( El, i );
145 
146  for ( int i = 0; i < numLocalDofs; ++i ) {
147  int glob_i = globalDofs[ i ];
148  for ( int j = 0; j < numLocalDofs; ++j ) {
149  int glob_j = globalDofs[ j ];
150  tripletList.push_back( TripletType( glob_i, glob_j, Factor * localMatrix(i,j) ) );
151  }
152  }
153  }
154  }
const TriangleType & getTriang(const int num) const
Definition: triangMesh.h:87
T Sqr(const T a)
Definition: aol.h:134
int getNumTriangs() const
Definition: triangMesh.h:46
int getNumLocalDofs(const ElementType &) const
Get The number of local degrees of freedom on an element.
int localToGlobal(const ElementType &T, int localIndex) const
Returns global index of the dof with number localIndex.
const InitType & getInitializer() const
Returns the mesh.
ConfiguratorType::LocalMatrixType LocalMatrixType
template<typename ConfiguratorType, typename Imp>
void shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::assembleTripletListCached ( ) const
inline

Definition at line 157 of file unityTriangleIntegratorShellFE.h.

157  {
158  if(!_tripletListAssembled) {
160  _tripletListAssembled = true;
161  }
162  }
void assembleTripletList(std::vector< TripletType > &tripletList, const RealType Factor) const

Member Data Documentation

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

Definition at line 198 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
RealType shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::_factor
protected

Definition at line 200 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
std::vector<TripletType> shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::_tripletList
mutableprotected

Definition at line 199 of file unityTriangleIntegratorShellFE.h.

template<typename ConfiguratorType, typename Imp>
bool shellFE::MatrixValuedIntegratorBase< ConfiguratorType, Imp >::_tripletListAssembled
mutableprotected

Definition at line 201 of file unityTriangleIntegratorShellFE.h.


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