NLMech  0.1.0
quadElem.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_QUADELEM_H
10 #define FE_QUADELEM_H
11 
12 #include "baseElem.h" // base class BaseElem
13 
14 namespace fe {
15 
62 class QuadElem : public BaseElem {
63 
64 public:
69  explicit QuadElem(size_t order);
70 
89  double elemSize(const std::vector<util::Point3> &nodes) override;
90 
136  std::vector<fe::QuadData>
137  getQuadDatas(const std::vector<util::Point3> &nodes) override;
138 
154  std::vector<fe::QuadData>
155  getQuadPoints(const std::vector<util::Point3> &nodes) override;
156 
166  std::string printStr(int nt = 0, int lvl = 0) const;
167 
175  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
176 
177 private:
184  std::vector<double> getShapes(const util::Point3 &p) override;
185 
193  std::vector<std::vector<double>> getDerShapes(const util::Point3 &p) override;
194 
203  double getJacobian(const util::Point3 &p,
204  const std::vector<util::Point3> &nodes,
205  std::vector<std::vector<double>> *J) override;
206 
210  void init() override;
211 
212 };
213 
214 } // namespace fe
215 
216 #endif // FE_QUADELEM_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 bi-linear quadrangle element.
Definition: quadElem.h:62
std::vector< std::vector< double > > getDerShapes(const util::Point3 &p) override
Returns the values of derivative of shape function at point p on reference element.
Definition: quadElem.cpp:97
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: quadElem.cpp:67
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: quadElem.cpp:303
std::vector< double > getShapes(const util::Point3 &p) override
Returns the values of shape function at point p on reference element.
Definition: quadElem.cpp:87
void init() override
Compute the quadrature points for quadrangle element.
Definition: quadElem.cpp:153
double elemSize(const std::vector< util::Point3 > &nodes) override
Returns the area of element.
Definition: quadElem.cpp:19
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: quadElem.cpp:119
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: quadElem.cpp:27
QuadElem(size_t order)
Constructor for quadrangle element.
Definition: quadElem.cpp:13
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: quadElem.h:175
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