NLMech  0.1.0
baseElem.h
1 // Copyright (c) 2019 Prashant K. Jha
3 // Copyright (c) 2019 Patrick Diehl
4 //
5 // Distributed under the Boost Software License, Version 1.0. (See accompanying
6 // file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
8 
9 #ifndef FE_BASEELEM_H
10 #define FE_BASEELEM_H
11 
12 #include "quadData.h" // definition of QuadData
13 #include "util/point.h" // definition of Point3
14 
15 namespace fe {
16 
82 class BaseElem {
83 
84 public:
91  BaseElem(size_t order, size_t element_type);
92 
97  size_t getElemType() { return d_elemType; }
98 
103  size_t getQuadOrder() { return d_quadOrder; }
104 
109  size_t getNumQuadPoints() { return d_numQuadPts; }
110 
118  virtual double elemSize(const std::vector<util::Point3> &nodes) = 0;
119 
136  virtual std::vector<double>
137  getShapes(const util::Point3 &p, const std::vector<util::Point3> &nodes);
138 
146  virtual std::vector<std::vector<double>>
147  getDerShapes(const util::Point3 &p,
148  const std::vector<util::Point3> &nodes);
149 
167  virtual std::vector<fe::QuadData>
168  getQuadDatas(const std::vector<util::Point3> &nodes) = 0;
169 
186  virtual std::vector<fe::QuadData>
187  getQuadPoints(const std::vector<util::Point3> &nodes) = 0;
188 
198  std::string printStr(int nt = 0, int lvl = 0) const;
199 
207  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
208 
209 protected:
217  virtual std::vector<double> getShapes(const util::Point3 &p) = 0;
218 
226  virtual std::vector<std::vector<double>>
227  getDerShapes(const util::Point3 &p) = 0;
228 
237  virtual util::Point3
239  const std::vector<util::Point3> &nodes);
240 
250  virtual double getJacobian(const util::Point3 &p,
251  const std::vector<util::Point3> &nodes,
252  std::vector<std::vector<double>> *J) = 0;
253 
259  virtual void init() = 0;
260 
262  size_t d_quadOrder;
263 
265  size_t d_numQuadPts;
266 
268  size_t d_elemType;
269 
271  std::vector<fe::QuadData> d_quads;
272 };
273 
274 } // namespace fe
275 
276 #endif // FE_BASEELEM_H
A base class which provides methods to map points to/from reference element and to compute quadrature...
Definition: baseElem.h:82
size_t getQuadOrder()
Get order of quadrature approximation.
Definition: baseElem.h:103
virtual std::vector< std::vector< double > > getDerShapes(const util::Point3 &p)=0
Returns the values of derivative of shape function at point p on reference element.
virtual std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point3 > &nodes)=0
Get vector of quadrature data.
virtual util::Point3 mapPointToRefElem(const util::Point3 &p, const std::vector< util::Point3 > &nodes)
Maps point p in a given element to the reference element and returns the mapped point.
Definition: baseElem.cpp:37
size_t d_numQuadPts
Number of quadrature points for order d_quadOrder.
Definition: baseElem.h:265
size_t getElemType()
Get element type.
Definition: baseElem.h:97
virtual std::vector< double > getShapes(const util::Point3 &p, const std::vector< util::Point3 > &nodes)
Returns the values of shape function at point p.
Definition: baseElem.cpp:19
BaseElem(size_t order, size_t element_type)
Constructor.
Definition: baseElem.cpp:14
std::vector< fe::QuadData > d_quads
Quadrature data collection.
Definition: baseElem.h:271
virtual std::vector< double > getShapes(const util::Point3 &p)=0
Returns the values of shape function at point p on reference element.
size_t getNumQuadPoints()
Get number of quadrature points in the data.
Definition: baseElem.h:109
size_t d_elemType
Element type.
Definition: baseElem.h:268
virtual std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point3 > &nodes)=0
Get vector of quadrature data.
virtual double elemSize(const std::vector< util::Point3 > &nodes)=0
Returns the size of element (length in 1-d, area in 2-d, volume in 3-d element)
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: baseElem.h:207
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: baseElem.cpp:50
size_t d_quadOrder
Order of quadrature point integration approximation.
Definition: baseElem.h:262
virtual double getJacobian(const util::Point3 &p, const std::vector< util::Point3 > &nodes, std::vector< std::vector< double >> *J)=0
Computes Jacobian of map from reference element to given element .
virtual std::vector< std::vector< double > > getDerShapes(const util::Point3 &p, const std::vector< util::Point3 > &nodes)
Returns the values of derivative of shape function at point p.
Definition: baseElem.cpp:28
virtual void init()=0
Compute the quadrature points.
Definition: baseElem.cpp:44
Collection of methods and data related to finite element and mesh.
Definition: baseElem.h:15
A structure to represent 3d vectors.
Definition: point.h:29