Programming tasks to Scientific Computing I
rhs.h
Go to the documentation of this file.
1 
5 #ifndef __FORCESSHELLFE_H
6 #define __FORCESSHELLFE_H
7 
8 #include <cmath>
9 
11 #include <shellHandler.h>
13 
14 namespace shellFE {
15 
17 template< typename ConfiguratorType >
19  public UnitTriangleFENonlinOpIntegratorShellFE <ConfiguratorType, createConstantRHS <ConfiguratorType> >
20 {
21 protected:
24 
25  const RealType _val;
26 
27 public:
28  createConstantRHS(const ConfiguratorType &conf, const RealType val)
30  _val(val)
31  {}
32 
33  RealType getNonlinearity ( const typename ConfiguratorType::ElementType &el, int /*quadPoint*/) const
34  {
35  return 2.0 * el.getAreaOfFlattenedTriangle() * _val;
36  }
37 };
38 
40 template< typename ConfiguratorType >
42  public UnitTriangleFENonlinOpIntegratorShellFE <ConfiguratorType, createNonConstantRHS <ConfiguratorType> >
43 {
44 protected:
47 
48  static const int dim = 3;
49 
50  std::unique_ptr<typename ConfiguratorType::VectorType> _dof;
51  std::unique_ptr<DiscreteVectorFunctionDefaultShellFE<ConfiguratorType>> _x;
52 
53 public:
56  {
57  int numGlobalDofs = this->_config.getNumGlobalDofs();
58  const typename ConfiguratorType::InitType &grid = this->_config.getInitializer();
59  _dof.reset(new typename ConfiguratorType::VectorType (dim * numGlobalDofs));
60 
61  for (int i = 0; i < numGlobalDofs; ++i) {
62  for (int c = 0; c < dim; ++c) {
63  (*_dof)[c*numGlobalDofs + i] = grid.getVertex(i)[c];
64  }
65  }
66 
67  _x.reset(new DiscreteVectorFunctionDefaultShellFE<ConfiguratorType>(this->_config, *_dof, dim));
68  }
69 
70  RealType getNonlinearity ( const typename ConfiguratorType::ElementType &el, int q) const
71  {
72  const typename ConfiguratorType::BaseFuncSetType &bfs = this->_config.getBaseFunctionSet(el);
74 
75  typename ConfiguratorType::TangentVecType cartCoord;
76  _x->evaluateAtQuadPoint(el, q, cartCoord);
77 
78  RealType nl = 2.0 * el.getAreaOfFlattenedTriangle() * evalRHS(cartCoord);
79 
80  return nl;
81  }
82 
83  RealType evalRHS(typename ConfiguratorType::TangentVecType &cartCoord) const
84  {
85  RealType nl = std::pow(aol::Sqr(-1. + 2. * cartCoord[0]) + aol::Sqr(-1. + 2. * cartCoord[1]), 1./3.) *
86  sin(2./3. * atan2(-1. + 2.*cartCoord[1], -1. + 2.*cartCoord[0]));
87 
88  return nl;
89  }
90 };
91 
93 template <typename ConfiguratorType>
95  const ShellHandler<ConfiguratorType> &shellHandler,
96  typename ConfiguratorType::RealType val,
97  typename ConfiguratorType::VectorType &rhs,
98  const bool collapseBoundaryValues = true)
99 {
100  const typename ConfiguratorType::MaskType boundaryMask(
101  shellHandler.getDirichletMask());
102  int numGlobalDofs = conf.getNumGlobalDofs();
103  rhs.setZero();
104 
106 
107  //bc rhs_force
108  if (collapseBoundaryValues) {
109  for (int i = 0; i < numGlobalDofs; ++i) {
110  if (boundaryMask[i]) {
111  rhs[i] = 0.0;
112  }
113  }
114  }
115 }
116 
118 template <typename ConfiguratorType>
120  const ShellHandler<ConfiguratorType> &shellHandler,
121  typename ConfiguratorType::RealType val,
122  typename ConfiguratorType::VectorType &rhs,
123  const bool collapseBoundaryValues = true)
124 {
125  const typename ConfiguratorType::MaskType boundaryMask(
126  shellHandler.getDirichletMask());
127  int numGlobalDofs = conf.getNumGlobalDofs();
128  rhs.setZero();
129 
131 
132  //bc rhs_force
133  if (collapseBoundaryValues) {
134  for (int i = 0; i < numGlobalDofs; ++i) {
135  if (boundaryMask[i]) {
136  rhs[i] = 0.0;
137  }
138  }
139  }
140 }
141 
142 } // end namespae
143 
144 #endif
DataTypeContainer::MaskType MaskType
std::unique_ptr< DiscreteVectorFunctionDefaultShellFE< ConfiguratorType > > _x
Definition: rhs.h:51
ConfiguratorType::TangentVecType TangentVecType
Definition: rhs.h:46
DataTypeContainer::VectorType VectorType
const MaskType & getDirichletMask() const
Definition: shellHandler.h:133
const Point3DType & getVertex(const int num) const
Definition: triangMesh.h:66
Integrator for , of some scalar valued function .
void createMaskedConstantRHS(const ConfiguratorType &conf, const ShellHandler< ConfiguratorType > &shellHandler, typename ConfiguratorType::RealType val, typename ConfiguratorType::VectorType &rhs, const bool collapseBoundaryValues=true)
Assembles constant rhs.
Definition: rhs.h:94
createNonConstantRHS(const ConfiguratorType &conf)
Definition: rhs.h:54
RealType getNonlinearity(const typename ConfiguratorType::ElementType &el, int) const
Definition: rhs.h:33
DataTypeContainer::DomVecType DomVecType
DataTypeContainer::RealType RealType
RealType evalRHS(typename ConfiguratorType::TangentVecType &cartCoord) const
Definition: rhs.h:83
ConfiguratorType::RealType RealType
Definition: rhs.h:22
T Sqr(const T a)
Definition: aol.h:134
Assembly operator for constant RHS.
Definition: rhs.h:18
const DomVecType & getRefCoord(int QuadPoint) const
Returns the coordinates of a quadrature point.
void createMaskedNonConstantRHS(const ConfiguratorType &conf, const ShellHandler< ConfiguratorType > &shellHandler, typename ConfiguratorType::RealType val, typename ConfiguratorType::VectorType &rhs, const bool collapseBoundaryValues=true)
Assembles sin-type rhs.
Definition: rhs.h:119
const RealType _val
Definition: rhs.h:25
Additional information about TriangleMeshes.
Definition: shellHandler.h:20
Helper class to evaluate a discrete vector-valued nodal function on a given mesh. ...
createConstantRHS(const ConfiguratorType &conf, const RealType val)
Definition: rhs.h:28
DataTypeContainer::TangentVecType TangentVecType
const InitType & getInitializer() const
Returns the mesh.
Definition: rhs.h:14
int getNumGlobalDofs() const
Returns the number of global degrees of freedom.
Assembly operator for non-constant RHS.
Definition: rhs.h:41
std::unique_ptr< typename ConfiguratorType::VectorType > _dof
Definition: rhs.h:50
ConfiguratorType::TangentVecType TangentVecType
Definition: rhs.h:23
ConfiguratorType::RealType RealType
Definition: rhs.h:45
Base classes for integration type operators defined using the unit triangle.
Configurator for Finite Elements.
Triangle which has a tangent space at each node.
const BaseFuncSetType & getBaseFunctionSet(const ElementType &) const
Returns the base funcitons set for an element.
RealType getNonlinearity(const typename ConfiguratorType::ElementType &el, int q) const
Definition: rhs.h:70