PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
tetElem.h
Go to the documentation of this file.
1/*
2 * -------------------------------------------
3 * Copyright (c) 2021 - 2024 Prashant K. Jha
4 * -------------------------------------------
5 * PeriDEM https://github.com/prashjha/PeriDEM
6 *
7 * Distributed under the Boost Software License, Version 1.0. (See accompanying
8 * file LICENSE)
9 */
10
11#ifndef FE_TETELEM_H
12#define FE_TETELEM_H
13
14#include "baseElem.h" // base class BaseElem
15
16namespace fe {
17
141class TetElem : public BaseElem {
142
143public:
148 explicit TetElem(size_t order);
149
171 double elemSize(const std::vector<util::Point> &nodes) override;
172
184 std::vector<double>
185 getShapes(const util::Point &p,
186 const std::vector<util::Point> &nodes) override;
187
247 std::vector<std::vector<double>>
248 getDerShapes(const util::Point &p,
249 const std::vector<util::Point> &nodes) override;
250
295 std::vector<fe::QuadData>
296 getQuadDatas(const std::vector<util::Point> &nodes) override;
297
313 std::vector<fe::QuadData>
314 getQuadPoints(const std::vector<util::Point> &nodes) override;
315
316private:
323 std::vector<double> getShapes(const util::Point &p) override;
324
332 std::vector<std::vector<double>> getDerShapes(const util::Point &p) override;
333
366 const std::vector<util::Point> &nodes) override;
367
376 double getJacobian(const util::Point &p,
377 const std::vector<util::Point> &nodes,
378 std::vector<std::vector<double>> *J) override;
379
383 void init() override;
384};
385
386} // namespace fe
387
388#endif // FE_TETELEM_H
A base class which provides methods to map points to/from reference element and to compute quadrature...
Definition baseElem.h:84
A class for mapping and quadrature related operations for linear tetrahedron element.
Definition tetElem.h:141
void init() override
Compute the quadrature points for triangle element.
Definition tetElem.cpp:270
std::vector< std::vector< double > > getDerShapes(const util::Point &p, const std::vector< util::Point > &nodes) override
Returns the values of derivative of shape function at point p.
Definition tetElem.cpp:95
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition tetElem.cpp:118
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition tetElem.cpp:158
double elemSize(const std::vector< util::Point > &nodes) override
Returns the volume of element.
Definition tetElem.cpp:81
double getJacobian(const util::Point &p, const std::vector< util::Point > &nodes, std::vector< std::vector< double > > *J) override
Computes the Jacobian of map .
Definition tetElem.cpp:237
util::Point mapPointToRefElem(const util::Point &p, const std::vector< util::Point > &nodes) override
Maps point p in a given element to the reference element.
Definition tetElem.cpp:203
std::vector< double > getShapes(const util::Point &p, const std::vector< util::Point > &nodes) override
Returns the values of shape function at point p.
Definition tetElem.cpp:90
Collection of methods and data related to finite element and mesh.
Definition baseElem.h:17
A structure to represent 3d vectors.
Definition point.h:30