NLMech  0.1.0
tetElem.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_TETELEM_H
10 #define FE_TETELEM_H
11 
12 #include "baseElem.h" // base class BaseElem
13 
14 namespace fe {
15 
139 class TetElem : public BaseElem {
140 
141 public:
146  explicit TetElem(size_t order);
147 
169  double elemSize(const std::vector<util::Point3> &nodes) override;
170 
182  std::vector<double>
183  getShapes(const util::Point3 &p,
184  const std::vector<util::Point3> &nodes) override;
185 
245  std::vector<std::vector<double>>
246  getDerShapes(const util::Point3 &p,
247  const std::vector<util::Point3> &nodes) override;
248 
293  std::vector<fe::QuadData>
294  getQuadDatas(const std::vector<util::Point3> &nodes) override;
295 
311  std::vector<fe::QuadData>
312  getQuadPoints(const std::vector<util::Point3> &nodes) override;
313 
323  std::string printStr(int nt = 0, int lvl = 0) const;
324 
332  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
333 
334 private:
341  std::vector<double> getShapes(const util::Point3 &p) override;
342 
350  std::vector<std::vector<double>> getDerShapes(const util::Point3 &p) override;
351 
385  const std::vector<util::Point3> &nodes) override;
386 
395  double getJacobian(const util::Point3 &p,
396  const std::vector<util::Point3> &nodes,
397  std::vector<std::vector<double>> *J) override;
398 
402  void init() override;
403 };
404 
405 } // namespace fe
406 
407 #endif // FE_TETELEM_H
A base class which provides methods to map points to/from reference element and to compute quadrature...
Definition: baseElem.h:82
A class for mapping and quadrature related operations for linear tetrahedron element.
Definition: tetElem.h:139
void init() override
Compute the quadrature points for triangle element.
Definition: tetElem.cpp:252
TetElem(size_t order)
Constructor.
Definition: tetElem.cpp:61
std::vector< double > getShapes(const util::Point3 &p, const std::vector< util::Point3 > &nodes) override
Returns the values of shape function at point p.
Definition: tetElem.cpp:82
util::Point3 mapPointToRefElem(const util::Point3 &p, const std::vector< util::Point3 > &nodes) override
Maps point p in a given element to the reference element.
Definition: tetElem.cpp:188
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: tetElem.cpp:394
double elemSize(const std::vector< util::Point3 > &nodes) override
Returns the volume of element.
Definition: tetElem.cpp:73
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: tetElem.cpp:146
std::vector< std::vector< double > > getDerShapes(const util::Point3 &p, const std::vector< util::Point3 > &nodes) override
Returns the values of derivative of shape function at point p.
Definition: tetElem.cpp:87
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: tetElem.cpp:108
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: tetElem.h:332
double getJacobian(const util::Point3 &p, const std::vector< util::Point3 > &nodes, std::vector< std::vector< double >> *J) override
Computes the Jacobian of map .
Definition: tetElem.cpp:219
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