NLMech  0.1.0
triElem.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_TRIELEM_H
10 #define FE_TRIELEM_H
11 
12 #include "baseElem.h" // base class BaseElem
13 
14 namespace fe {
15 
89 class TriElem : public BaseElem {
90 
91 public:
96  explicit TriElem(size_t order);
97 
114  double elemSize(const std::vector<util::Point3> &nodes) override;
115 
127  std::vector<double>
128  getShapes(const util::Point3 &p,
129  const std::vector<util::Point3> &nodes) override;
130 
178  std::vector<std::vector<double>>
179  getDerShapes(const util::Point3 &p,
180  const std::vector<util::Point3> &nodes) override;
181 
225  std::vector<fe::QuadData>
226  getQuadDatas(const std::vector<util::Point3> &nodes) override;
227 
243  std::vector<fe::QuadData>
244  getQuadPoints(const std::vector<util::Point3> &nodes) override;
245 
255  std::string printStr(int nt = 0, int lvl = 0) const;
256 
264  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
265 
266 private:
273  std::vector<double> getShapes(const util::Point3 &p) override;
274 
282  std::vector<std::vector<double>> getDerShapes(const util::Point3 &p) override;
283 
317  const std::vector<util::Point3> &nodes) override;
318 
327  double getJacobian(const util::Point3 &p,
328  const std::vector<util::Point3> &nodes,
329  std::vector<std::vector<double>> *J) override;
330 
334  void init() override;
335 };
336 
337 } // namespace fe
338 
339 #endif // FE_TRIELEM_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 triangle element.
Definition: triElem.h:89
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: triElem.cpp:27
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: triElem.cpp:161
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: triElem.h:264
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: triElem.cpp:32
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: triElem.cpp:54
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: triElem.cpp:131
double elemSize(const std::vector< util::Point3 > &nodes) override
Returns the area of element.
Definition: triElem.cpp:22
void init() override
Compute the quadrature points for triangle element.
Definition: triElem.cpp:178
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: triElem.cpp:407
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: triElem.cpp:94
TriElem(size_t order)
Constructor.
Definition: triElem.cpp:16
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