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

Additional information about TriangleMeshes. More...

#include <shellHandler.h>

Public Types

typedef ConfiguratorType::RealType RealType
 
typedef ConfiguratorType::InitType MeshType
 
typedef ConfiguratorType::MaskType MaskType
 
typedef ConfiguratorType::TangentVecType TangentVecType
 
typedef ConfiguratorType::Point3DType Point3DType
 
typedef ConfiguratorType::VectorType VectorType
 

Public Member Functions

 ShellHandler (const ConfiguratorType &conf)
 
void generateChart_xA () const
 
VectorTypegetChartToUndeformedShell () const
 
void generateDirichletBoundaryMask (MaskType &mask, int &numBoundaryNodes, ShellBoundaryType boundaryType=ALLBOUNDARY, bool clampedBoundaryCondition=false) const
 
const MaskTypegetDirichletMask () const
 
const int & getNumBoundaryNodes () const
 

Public Attributes

const ConfiguratorType_conf
 
const MeshType_mesh
 
const int _numVertices
 
const int _numGlobalDofs
 
MaskType _DirichletMask
 
int _numBoundaryNodes
 
VectorType _xA
 

Detailed Description

template<typename ConfiguratorType>
class shellFE::ShellHandler< ConfiguratorType >

Additional information about TriangleMeshes.

Definition at line 20 of file shellHandler.h.

Member Typedef Documentation

template<typename ConfiguratorType >
typedef ConfiguratorType::MaskType shellFE::ShellHandler< ConfiguratorType >::MaskType

Definition at line 24 of file shellHandler.h.

template<typename ConfiguratorType >
typedef ConfiguratorType::InitType shellFE::ShellHandler< ConfiguratorType >::MeshType

Definition at line 23 of file shellHandler.h.

template<typename ConfiguratorType >
typedef ConfiguratorType::Point3DType shellFE::ShellHandler< ConfiguratorType >::Point3DType

Definition at line 26 of file shellHandler.h.

template<typename ConfiguratorType >
typedef ConfiguratorType::RealType shellFE::ShellHandler< ConfiguratorType >::RealType

Definition at line 22 of file shellHandler.h.

Definition at line 25 of file shellHandler.h.

template<typename ConfiguratorType >
typedef ConfiguratorType::VectorType shellFE::ShellHandler< ConfiguratorType >::VectorType

Definition at line 27 of file shellHandler.h.

Constructor & Destructor Documentation

template<typename ConfiguratorType >
shellFE::ShellHandler< ConfiguratorType >::ShellHandler ( const ConfiguratorType conf)
inline

Definition at line 38 of file shellHandler.h.

38  :
39  _conf ( conf ),
40  _mesh( conf.getInitializer() ),
43  _xA ( 3 * conf.getNumGlobalDofs() )
44  {
45  generateChart_xA ( );
47  }
int getNumVertices() const
Definition: triangMesh.h:42
const int _numGlobalDofs
Definition: shellHandler.h:31
const MeshType & _mesh
Definition: shellHandler.h:30
void generateDirichletBoundaryMask(MaskType &mask, int &numBoundaryNodes, ShellBoundaryType boundaryType=ALLBOUNDARY, bool clampedBoundaryCondition=false) const
Definition: shellHandler.h:64
void generateChart_xA() const
Definition: shellHandler.h:49
const InitType & getInitializer() const
Returns the mesh.
int getNumGlobalDofs() const
Returns the number of global degrees of freedom.
const ConfiguratorType & _conf
Definition: shellHandler.h:29

Member Function Documentation

template<typename ConfiguratorType >
void shellFE::ShellHandler< ConfiguratorType >::generateChart_xA ( ) const
inline

Definition at line 49 of file shellHandler.h.

50  {
51  for ( int nodeIdx=0; nodeIdx < _numVertices; ++nodeIdx ) {
52  const Point3DType& coords ( _mesh.getVertex(nodeIdx) );
53 
54  for( int comp=0; comp<3; ++comp )
55  _xA[nodeIdx + _numGlobalDofs * comp] = coords[comp];
56  }
57  }
const Point3DType & getVertex(const int num) const
Definition: triangMesh.h:66
const int _numGlobalDofs
Definition: shellHandler.h:31
const MeshType & _mesh
Definition: shellHandler.h:30
ConfiguratorType::Point3DType Point3DType
Definition: shellHandler.h:26
template<typename ConfiguratorType >
void shellFE::ShellHandler< ConfiguratorType >::generateDirichletBoundaryMask ( MaskType mask,
int &  numBoundaryNodes,
ShellBoundaryType  boundaryType = ALLBOUNDARY,
bool  clampedBoundaryCondition = false 
) const
inline

Definition at line 64 of file shellHandler.h.

66  {
67  mask.resize( _numGlobalDofs, false );
68  numBoundaryNodes = 0;
69 
70  switch ( boundaryType ){
71  case NOBOUNDARY : {
72  } break;
73 
74  case ALLBOUNDARY : {
75  for( int ElementIndex = 0; ElementIndex < _mesh.getNumTriangs(); ++ElementIndex ){
76  for( int i = 0; i < 3; i++ ){
77  if ( _mesh.getNeighbour( ElementIndex , i ) == -1 ){
78  mask[ _mesh.getTriangNodeIdx( ElementIndex, (i+1)%3) ] = true;
79  mask[ _mesh.getTriangNodeIdx( ElementIndex, (i+2)%3) ] = true;
80  }
81  }
82  }
83  } break;
84 
85  case PlateLeft : {
86  for ( int nodeIdx=0; nodeIdx < _mesh.getNumVertices(); ++nodeIdx ) {
87  const Point3DType& coords ( _mesh.getVertex(nodeIdx) );
88  if (coords [0] == 0. ){
89  ++numBoundaryNodes;
90  mask[nodeIdx] = true;
91  }
92  }
93  } break;
94 
95  case PlateLeftTop : {
96  for ( int nodeIdx=0; nodeIdx < _mesh.getNumVertices(); ++nodeIdx ) {
97  const Point3DType& coords ( _mesh.getVertex(nodeIdx) );
98  if (coords [0] == 0. || coords[1] == 1. ){
99  ++numBoundaryNodes;
100  mask[nodeIdx] = true;
101  }
102  }
103  } break;
104 
105  case PlateAll : {
106  for ( int nodeIdx=0; nodeIdx < _mesh.getNumVertices(); ++nodeIdx ) {
107  const Point3DType& coords ( _mesh.getVertex(nodeIdx) );
108  if (coords [0] == 0. || coords [0] == 1. || coords [1] == 0. || coords[1] == 1. ){
109  ++numBoundaryNodes;
110  mask[nodeIdx] = true;
111  }
112  }
113  } break;
114 
115 
116  default :
117  throw std::invalid_argument( aol::strprintf ( "Wrong boundary condition. In File %s at line %d.", __FILE__, __LINE__ ).c_str() );
118  break;
119  }
120 
121  //clamped boundary condition
122  if( (ConfiguratorType::_ShellFEType == C1Dofs) && (clampedBoundaryCondition) ){
123  for( int i=0; i<_numVertices; ++i ){
124  if( mask[i] ){
125  mask[ i + _numVertices] = true;
126  mask[ i + 2. * _numVertices] = true;
127  }
128  }
129  }
130 
131  }
const Point3DType & getVertex(const int num) const
Definition: triangMesh.h:66
int getNumVertices() const
Definition: triangMesh.h:42
int getNeighbour(const int elementID, const int acrossLocalEdge) const
Get neighbor on edge.
Definition: triangMesh.h:76
const int _numGlobalDofs
Definition: shellHandler.h:31
const MeshType & _mesh
Definition: shellHandler.h:30
ConfiguratorType::Point3DType Point3DType
Definition: shellHandler.h:26
int getNumTriangs() const
Definition: triangMesh.h:46
int getTriangNodeIdx(const int num, const int localNode) const
Definition: triangMesh.h:104
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
template<typename ConfiguratorType >
VectorType& shellFE::ShellHandler< ConfiguratorType >::getChartToUndeformedShell ( ) const
inline

Definition at line 59 of file shellHandler.h.

60  {
61  return _xA;
62  }
template<typename ConfiguratorType >
const MaskType& shellFE::ShellHandler< ConfiguratorType >::getDirichletMask ( ) const
inline

Definition at line 133 of file shellHandler.h.

133 { return _DirichletMask;}
template<typename ConfiguratorType >
const int& shellFE::ShellHandler< ConfiguratorType >::getNumBoundaryNodes ( ) const
inline

Definition at line 134 of file shellHandler.h.

134 { return _numBoundaryNodes; }

Member Data Documentation

template<typename ConfiguratorType >
const ConfiguratorType& shellFE::ShellHandler< ConfiguratorType >::_conf

Definition at line 29 of file shellHandler.h.

template<typename ConfiguratorType >
MaskType shellFE::ShellHandler< ConfiguratorType >::_DirichletMask

Definition at line 32 of file shellHandler.h.

template<typename ConfiguratorType >
const MeshType& shellFE::ShellHandler< ConfiguratorType >::_mesh

Definition at line 30 of file shellHandler.h.

template<typename ConfiguratorType >
int shellFE::ShellHandler< ConfiguratorType >::_numBoundaryNodes

Definition at line 33 of file shellHandler.h.

template<typename ConfiguratorType >
const int shellFE::ShellHandler< ConfiguratorType >::_numGlobalDofs

Definition at line 31 of file shellHandler.h.

template<typename ConfiguratorType >
const int shellFE::ShellHandler< ConfiguratorType >::_numVertices

Definition at line 31 of file shellHandler.h.

template<typename ConfiguratorType >
VectorType shellFE::ShellHandler< ConfiguratorType >::_xA
mutable

Definition at line 34 of file shellHandler.h.


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