PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
baseElem.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_BASEELEM_H
12#define FE_BASEELEM_H
13
14#include "quadData.h" // definition of QuadData
15#include "util/point.h" // definition of Point
16
17namespace fe {
18
84class BaseElem {
85
86public:
93 BaseElem(size_t order, size_t element_type);
94
99 size_t getElemType() { return d_elemType; }
100
105 size_t getQuadOrder() { return d_quadOrder; }
106
111 size_t getNumQuadPoints() { return d_numQuadPts; }
112
120 virtual double elemSize(const std::vector<util::Point> &nodes) = 0;
121
138 virtual std::vector<double>
139 getShapes(const util::Point &p, const std::vector<util::Point> &nodes);
140
148 virtual std::vector<std::vector<double>>
149 getDerShapes(const util::Point &p,
150 const std::vector<util::Point> &nodes);
151
169 virtual std::vector<fe::QuadData>
170 getQuadDatas(const std::vector<util::Point> &nodes) = 0;
171
188 virtual std::vector<fe::QuadData>
189 getQuadPoints(const std::vector<util::Point> &nodes) = 0;
190
191protected:
199 virtual std::vector<double> getShapes(const util::Point &p) = 0;
200
208 virtual std::vector<std::vector<double>>
210
220 const std::vector<util::Point> &nodes);
221
231 virtual double getJacobian(const util::Point &p,
232 const std::vector<util::Point> &nodes,
233 std::vector<std::vector<double>> *J) = 0;
234
240 virtual void init() = 0;
241
244
247
250
252 std::vector<fe::QuadData> d_quads;
253};
254
255} // namespace fe
256
257#endif // FE_BASEELEM_H
A base class which provides methods to map points to/from reference element and to compute quadrature...
Definition baseElem.h:84
size_t getQuadOrder()
Get order of quadrature approximation.
Definition baseElem.h:105
virtual util::Point mapPointToRefElem(const util::Point &p, const std::vector< util::Point > &nodes)
Maps point p in a given element to the reference element and returns the mapped point.
Definition baseElem.cpp:40
size_t d_numQuadPts
Number of quadrature points for order d_quadOrder.
Definition baseElem.h:246
virtual std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point > &nodes)=0
Get vector of quadrature data.
virtual std::vector< double > getShapes(const util::Point &p)=0
Returns the values of shape function at point p on reference element.
size_t getElemType()
Get element type.
Definition baseElem.h:99
std::vector< fe::QuadData > d_quads
Quadrature data collection.
Definition baseElem.h:252
virtual std::vector< std::vector< double > > getDerShapes(const util::Point &p, const std::vector< util::Point > &nodes)
Returns the values of derivative of shape function at point p.
Definition baseElem.cpp:30
virtual std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point > &nodes)=0
Get vector of quadrature data.
size_t getNumQuadPoints()
Get number of quadrature points in the data.
Definition baseElem.h:111
virtual std::vector< double > getShapes(const util::Point &p, const std::vector< util::Point > &nodes)
Returns the values of shape function at point p.
Definition baseElem.cpp:20
size_t d_elemType
Element type.
Definition baseElem.h:249
virtual std::vector< std::vector< double > > getDerShapes(const util::Point &p)=0
Returns the values of derivative of shape function at point p on reference element.
virtual double elemSize(const std::vector< util::Point > &nodes)=0
Returns the size of element (length in 1-d, area in 2-d, volume in 3-d element)
virtual double getJacobian(const util::Point &p, const std::vector< util::Point > &nodes, std::vector< std::vector< double > > *J)=0
Computes Jacobian of map from reference element to given element .
size_t d_quadOrder
Order of quadrature point integration approximation.
Definition baseElem.h:243
virtual void init()=0
Compute the quadrature points.
Definition baseElem.cpp:47
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