13 #include "util/point.h"
91 BaseElem(
size_t order,
size_t element_type);
118 virtual double elemSize(
const std::vector<util::Point3> &nodes) = 0;
136 virtual std::vector<double>
146 virtual std::vector<std::vector<double>>
148 const std::vector<util::Point3> &nodes);
167 virtual std::vector<fe::QuadData>
186 virtual std::vector<fe::QuadData>
198 std::string
printStr(
int nt = 0,
int lvl = 0)
const;
207 void print(
int nt = 0,
int lvl = 0)
const { std::cout <<
printStr(nt, lvl); };
226 virtual std::vector<std::vector<double>>
239 const std::vector<util::Point3> &nodes);
251 const std::vector<util::Point3> &nodes,
252 std::vector<std::vector<double>> *J) = 0;
259 virtual void init() = 0;
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