24 return (nodes[1].d_x - nodes[0].d_x);
29 const std::vector<util::Point> &nodes) {
30 return getShapes(mapPointToRefElem(p, nodes));
33std::vector<std::vector<double>>
35 const std::vector<util::Point> &nodes) {
37 auto ders_ref = getDerShapes(mapPointToRefElem(p, nodes));
40 auto detJ = getJacobian(p, nodes,
nullptr);
43 ders_ref[0][0] = ders_ref[0][0] / detJ;
44 ders_ref[1][0] = ders_ref[1][0] / detJ;
49std::vector<fe::QuadData>
56 for (
auto &qd : qds) {
59 qd.d_detJ = getJacobian(qd.d_p, nodes, &(qd.d_J));
65 qd.d_p.d_x = qd.d_shapes[0] * nodes[0].d_x + qd.d_shapes[1] * nodes[1]
69 for (
size_t i=0; i<2; i++)
70 qd.d_derShapes[i][0] = qd.d_derShapes[i][0] / qd.d_detJ;
76std::vector<fe::QuadData>
83 for (
auto &qd : qds) {
86 qd.d_w *= getJacobian(qd.d_p, nodes,
nullptr);
89 qd.d_p.d_x = qd.d_shapes[0] * nodes[0].d_x + qd.d_shapes[1] * nodes[1]
100 return std::vector<double>{0.5 * (1. - p.
d_x), 0.5 * (1. + p.
d_x)};
103std::vector<std::vector<double>>
111 std::vector<std::vector<double>> r;
112 r.push_back(std::vector<double>{-0.5});
113 r.push_back(std::vector<double>{0.5});
119 const std::vector<util::Point> &nodes) {
120 auto xi = (2. * p.
d_x - nodes[0].d_x - nodes[1].d_x) /
121 (nodes[0].d_x - nodes[1].d_x);
125 std::cerr <<
"Error: Trying to map point p = " << p.
d_x
126 <<
" in given line to reference line.\n"
127 <<
"But the point p does not belong to line = {"
128 << nodes[0].d_x <<
", " << nodes[1].d_x <<
"}.\n";
141 const std::vector<util::Point> &nodes,
142 std::vector<std::vector<double>> *J) {
146 (*J)[0] = std::vector<double>{nodes[1].d_x - nodes[0].d_x};
151 return (nodes[1].d_x - nodes[0].d_x);
164 if (!d_quads.empty())
168 if (d_quadOrder == 0)
172 std::vector<std::vector<double>> ident_mat;
173 ident_mat.push_back(std::vector<double>{1.});
178 if (d_quadOrder == 1) {
188 d_quads.push_back(qd);
194 if (d_quadOrder == 2) {
204 d_quads.push_back(qd);
212 d_quads.push_back(qd);
218 if (d_quadOrder == 3) {
229 d_quads.push_back(qd);
237 d_quads.push_back(qd);
245 d_quads.push_back(qd);
251 if (d_quadOrder == 4) {
254 qd.
d_w = 0.6521451548625461;
260 d_quads.push_back(qd);
262 qd.
d_w = 0.6521451548625461;
268 d_quads.push_back(qd);
270 qd.
d_w = 0.3478548451374538;
276 d_quads.push_back(qd);
278 qd.
d_w = 0.3478548451374538;
284 d_quads.push_back(qd);
290 if (d_quadOrder == 5) {
293 qd.
d_w = 0.5688888888888889;
299 d_quads.push_back(qd);
301 qd.
d_w = 0.4786286704993665;
307 d_quads.push_back(qd);
309 qd.
d_w = 0.4786286704993665;
315 d_quads.push_back(qd);
317 qd.
d_w = 0.2369268850561891;
323 d_quads.push_back(qd);
325 qd.
d_w = 0.2369268850561891;
331 d_quads.push_back(qd);
A base class which provides methods to map points to/from reference element and to compute quadrature...
std::vector< double > getShapes(const util::Point &p, const std::vector< util::Point > &nodes) override
Returns the values of shape function at point p.
std::vector< fe::QuadData > getQuadPoints(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
std::vector< fe::QuadData > getQuadDatas(const std::vector< util::Point > &nodes) override
Get vector of quadrature data.
LineElem(size_t order)
Constructor for line element.
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.
void init() override
Compute the quadrature points for line element.
double getJacobian(const util::Point &p, const std::vector< util::Point > &nodes, std::vector< std::vector< double > > *J) override
Computes Jacobian of the map .
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.
double elemSize(const std::vector< util::Point > &nodes) override
Returns the length of element.
Collection of methods and data related to finite element and mesh.
Collection of methods useful in simulation.
bool isGreater(const double &a, const double &b)
Returns true if a > b.
bool isLess(const double &a, const double &b)
Returns true if a < b.
A struct to store the quadrature data. List of data are.
std::vector< double > d_shapes
Value of shape functions at quad point p.
std::vector< std::vector< double > > d_derShapes
Value of derivative of shape functions at quad point p.
double d_w
Quadrature weight.
util::Point d_p
Quadrature point in 1-d, 2-d or 3-d.
std::vector< std::vector< double > > d_J
Jacobian of the map from reference element to the element.
double d_detJ
Determinant of the Jacobian of the map from reference element to the element.
A structure to represent 3d vectors.
double d_x
the x coordinate