PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
lineElem.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#ifndef FE_LINEELEM_H
11#define FE_LINEELEM_H
12
13#include "baseElem.h" // base class BaseElem
14
15namespace fe {
16
49class LineElem : public BaseElem {
50
51public:
56 explicit LineElem(size_t order);
57
67 double elemSize(const std::vector<util::Point> &nodes) override;
68
81 std::vector<double>
82 getShapes(const util::Point &p,
83 const std::vector<util::Point> &nodes) override;
84
108 std::vector<std::vector<double>>
109 getDerShapes(const util::Point &p,
110 const std::vector<util::Point> &nodes) override;
111
153 std::vector<fe::QuadData>
154 getQuadDatas(const std::vector<util::Point> &nodes) override;
155
172 std::vector<fe::QuadData>
173 getQuadPoints(const std::vector<util::Point> &nodes) override;
174
175private:
182 std::vector<double> getShapes(const util::Point &p) override;
183
191 std::vector<std::vector<double>> getDerShapes(const util::Point &p) override;
192
212 const std::vector<util::Point> &nodes) override;
213
222 double getJacobian(const util::Point &p,
223 const std::vector<util::Point> &nodes,
224 std::vector<std::vector<double>> *J) override;
225
229 void init() override;
230
231};
232
233} // namespace fe
234
235#endif // FE_LINEELEM_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 2-node line element.
Definition lineElem.h:49
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 lineElem.cpp:28
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition lineElem.cpp:77
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
Definition lineElem.cpp:50
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 lineElem.cpp:118
void init() override
Compute the quadrature points for line element.
Definition lineElem.cpp:154
double getJacobian(const util::Point &p, const std::vector< util::Point > &nodes, std::vector< std::vector< double > > *J) override
Computes Jacobian of the map .
Definition lineElem.cpp:140
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 lineElem.cpp:34
double elemSize(const std::vector< util::Point > &nodes) override
Returns the length of element.
Definition lineElem.cpp:23
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