PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
triElem.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_TRIELEM_H
12#define FE_TRIELEM_H
13
14#include "baseElem.h" // base class BaseElem
15
16namespace fe {
17
91class TriElem : public BaseElem {
92
93public:
98 explicit TriElem(size_t order);
99
116 double elemSize(const std::vector<util::Point> &nodes) override;
117
129 std::vector<double>
130 getShapes(const util::Point &p,
131 const std::vector<util::Point> &nodes) override;
132
180 std::vector<std::vector<double>>
181 getDerShapes(const util::Point &p,
182 const std::vector<util::Point> &nodes) override;
183
227 std::vector<fe::QuadData>
228 getQuadDatas(const std::vector<util::Point> &nodes) override;
229
245 std::vector<fe::QuadData>
246 getQuadPoints(const std::vector<util::Point> &nodes) override;
247
248private:
255 std::vector<double> getShapes(const util::Point &p) override;
256
264 std::vector<std::vector<double>> getDerShapes(const util::Point &p) override;
265
298 const std::vector<util::Point> &nodes) override;
299
308 double getJacobian(const util::Point &p,
309 const std::vector<util::Point> &nodes,
310 std::vector<std::vector<double>> *J) override;
311
315 void init() override;
316};
317
318} // namespace fe
319
320#endif // FE_TRIELEM_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 triangle element.
Definition triElem.h:91
double elemSize(const std::vector< util::Point > &nodes) override
Returns the area of element.
Definition triElem.cpp:23
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition triElem.cpp:98
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 triElem.cpp:170
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition triElem.cpp:58
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 triElem.cpp:35
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 triElem.cpp:29
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 triElem.cpp:139
void init() override
Compute the quadrature points for triangle element.
Definition triElem.cpp:188
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