NLMech  0.1.0
lineElem.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_LINEELEM_H
10 #define FE_LINEELEM_H
11 
12 #include "baseElem.h" // base class BaseElem
13 
14 namespace fe {
15 
48 class LineElem : public BaseElem {
49 
50 public:
55  explicit LineElem(size_t order);
56 
66  double elemSize(const std::vector<util::Point3> &nodes) override;
67 
80  std::vector<double>
81  getShapes(const util::Point3 &p,
82  const std::vector<util::Point3> &nodes) override;
83 
107  std::vector<std::vector<double>>
108  getDerShapes(const util::Point3 &p,
109  const std::vector<util::Point3> &nodes) override;
110 
152  std::vector<fe::QuadData>
153  getQuadDatas(const std::vector<util::Point3> &nodes) override;
154 
171  std::vector<fe::QuadData>
172  getQuadPoints(const std::vector<util::Point3> &nodes) override;
173 
183  std::string printStr(int nt = 0, int lvl = 0) const;
184 
192  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
193 
194 private:
201  std::vector<double> getShapes(const util::Point3 &p) override;
202 
210  std::vector<std::vector<double>> getDerShapes(const util::Point3 &p) override;
211 
232  const std::vector<util::Point3> &nodes) override;
233 
242  double getJacobian(const util::Point3 &p,
243  const std::vector<util::Point3> &nodes,
244  std::vector<std::vector<double>> *J) override;
245 
249  void init() override;
250 
251 };
252 
253 } // namespace fe
254 
255 #endif // FE_LINEELEM_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 2-node line element.
Definition: lineElem.h:48
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: lineElem.cpp:45
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: lineElem.cpp:25
LineElem(size_t order)
Constructor for line element.
Definition: lineElem.cpp:15
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: lineElem.h:192
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: lineElem.cpp:316
void init() override
Compute the quadrature points for line element.
Definition: lineElem.cpp:138
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: lineElem.cpp:30
double elemSize(const std::vector< util::Point3 > &nodes) override
Returns the length of element.
Definition: lineElem.cpp:21
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: lineElem.cpp:105
double getJacobian(const util::Point3 &p, const std::vector< util::Point3 > &nodes, std::vector< std::vector< double >> *J) override
Computes Jacobian of the map .
Definition: lineElem.cpp:125
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point3 > &nodes) override
Get vector of quadrature data.
Definition: lineElem.cpp:69
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