PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
quadElem.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_QUADELEM_H
12#define FE_QUADELEM_H
13
14#include "baseElem.h" // base class BaseElem
15
16namespace fe {
17
64class QuadElem : public BaseElem {
65
66public:
71 explicit QuadElem(size_t order);
72
91 double elemSize(const std::vector<util::Point> &nodes) override;
92
138 std::vector<fe::QuadData>
139 getQuadDatas(const std::vector<util::Point> &nodes) override;
140
156 std::vector<fe::QuadData>
157 getQuadPoints(const std::vector<util::Point> &nodes) override;
158
159private:
166 std::vector<double> getShapes(const util::Point &p) override;
167
175 std::vector<std::vector<double>> getDerShapes(const util::Point &p) override;
176
185 double getJacobian(const util::Point &p,
186 const std::vector<util::Point> &nodes,
187 std::vector<std::vector<double>> *J) override;
188
192 void init() override;
193
194};
195
196} // namespace fe
197
198#endif // FE_QUADELEM_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 bi-linear quadrangle element.
Definition quadElem.h:64
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 quadElem.cpp:127
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition quadElem.cpp:32
void init() override
Compute the quadrature points for quadrangle element.
Definition quadElem.cpp:162
std::vector< double > getShapes(const util::Point &p) override
Returns the values of shape function at point p on reference element.
Definition quadElem.cpp:93
std::vector< std::vector< double > > getDerShapes(const util::Point &p) override
Returns the values of derivative of shape function at point p on reference element.
Definition quadElem.cpp:105
double elemSize(const std::vector< util::Point > &nodes) override
Returns the area of element.
Definition quadElem.cpp:21
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition quadElem.cpp:72
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