PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
lineElem.cpp
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#include "lineElem.h"
12#include "util/feElementDefs.h" // global definition of elements
13#include "util/function.h"
14#include <iostream>
15
17 : fe::BaseElem(order, util::vtk_type_line) {
18
19 // compute quad data
20 this->init();
21}
22
23double fe::LineElem::elemSize(const std::vector<util::Point> &nodes) {
24 return (nodes[1].d_x - nodes[0].d_x);
25}
26
27std::vector<double>
29 const std::vector<util::Point> &nodes) {
30 return getShapes(mapPointToRefElem(p, nodes));
31}
32
33std::vector<std::vector<double>>
35 const std::vector<util::Point> &nodes) {
36 // get derivatives of shape function in reference triangle
37 auto ders_ref = getDerShapes(mapPointToRefElem(p, nodes));
38
39 // get Jacobian
40 auto detJ = getJacobian(p, nodes, nullptr);
41
42 // modify derivatives of shape function
43 ders_ref[0][0] = ders_ref[0][0] / detJ;
44 ders_ref[1][0] = ders_ref[1][0] / detJ;
45
46 return ders_ref;
47}
48
49std::vector<fe::QuadData>
50fe::LineElem::getQuadDatas(const std::vector<util::Point> &nodes) {
51
52 // copy quad data associated to reference element
53 auto qds = d_quads;
54
55 // modify data
56 for (auto &qd : qds) {
57
58 // get Jacobian and determinant
59 qd.d_detJ = getJacobian(qd.d_p, nodes, &(qd.d_J));
60
61 // transform quad weight
62 qd.d_w *= qd.d_detJ;
63
64 // map point to line
65 qd.d_p.d_x = qd.d_shapes[0] * nodes[0].d_x + qd.d_shapes[1] * nodes[1]
66 .d_x;
67
68 // modify derivative of shape function
69 for (size_t i=0; i<2; i++)
70 qd.d_derShapes[i][0] = qd.d_derShapes[i][0] / qd.d_detJ;
71 }
72
73 return qds;
74}
75
76std::vector<fe::QuadData>
77fe::LineElem::getQuadPoints(const std::vector<util::Point> &nodes) {
78
79 // copy quad data associated to reference element
80 auto qds = d_quads;
81
82 // modify data
83 for (auto &qd : qds) {
84
85 // transform quad weight
86 qd.d_w *= getJacobian(qd.d_p, nodes, nullptr);
87
88 // map point to line
89 qd.d_p.d_x = qd.d_shapes[0] * nodes[0].d_x + qd.d_shapes[1] * nodes[1]
90 .d_x;
91 }
92
93 return qds;
94}
95
96std::vector<double> fe::LineElem::getShapes(const util::Point &p) {
97
98 // N1 = (1 - xi)/2
99 // N2 = (1 + xi)/2
100 return std::vector<double>{0.5 * (1. - p.d_x), 0.5 * (1. + p.d_x)};
101}
102
103std::vector<std::vector<double>>
105
106 // N1 = (1 - xi)/2
107 // --> d N1/d xi = -1/2
108 //
109 // N2 = (1 + xi)/2
110 // --> d N2/d xi = 1/2
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});
114 return r;
115}
116
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);
122
123 if (util::isLess(xi, -1. - 1.0E-8) ||
124 util::isGreater(xi, 1. + 1.0E-8) ) {
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";
129 exit(1);
130 }
131
132 if (util::isLess(xi, -1.))
133 xi = -1.;
134 if (util::isGreater(xi, 1.))
135 xi = 1.;
136
137 return {xi, 0., 0.};
138}
139
141 const std::vector<util::Point> &nodes,
142 std::vector<std::vector<double>> *J) {
143
144 if (J != nullptr) {
145 J->resize(1);
146 (*J)[0] = std::vector<double>{nodes[1].d_x - nodes[0].d_x};
147
148 return (*J)[0][0];
149 }
150
151 return (nodes[1].d_x - nodes[0].d_x);
152}
153
155
156 //
157 // compute quad data for reference line element with vertex at p1 = -1 and
158 // p2 = 1
159 //
160 // Shape functions are
161 // N1 = (1 - xi)/2
162 // N2 = (1 + xi)/2
163
164 if (!d_quads.empty())
165 return;
166
167 // no point in zeroth order
168 if (d_quadOrder == 0)
169 d_quads.resize(0);
170
171 // 1x1 identity matrix
172 std::vector<std::vector<double>> ident_mat;
173 ident_mat.push_back(std::vector<double>{1.});
174
175 //
176 // first order quad points
177 //
178 if (d_quadOrder == 1) {
179 d_quads.clear();
180 // 1-d points are: {0} and weights are: {2}
181 fe::QuadData qd;
182 qd.d_w = 2.;
183 qd.d_p = util::Point();
184 qd.d_shapes = getShapes(qd.d_p);
185 qd.d_derShapes = getDerShapes(qd.d_p);
186 qd.d_J = ident_mat;
187 qd.d_detJ = 1.;
188 d_quads.push_back(qd);
189 }
190
191 //
192 // second order quad points
193 //
194 if (d_quadOrder == 2) {
195 d_quads.clear();
196 // 1-d points are: {-1/sqrt{3], 1/sqrt{3}} and weights are: {1,1}
197 fe::QuadData qd;
198 qd.d_w = 1.;
199 qd.d_p = util::Point(-1. / std::sqrt(3.), 0., 0.);
200 qd.d_shapes = getShapes(qd.d_p);
201 qd.d_derShapes = getDerShapes(qd.d_p);
202 qd.d_J = ident_mat;
203 qd.d_detJ = 1.;
204 d_quads.push_back(qd);
205
206 qd.d_w = 1.;
207 qd.d_p = util::Point(1. / std::sqrt(3.), 0., 0.);
208 qd.d_shapes = getShapes(qd.d_p);
209 qd.d_derShapes = getDerShapes(qd.d_p);
210 qd.d_J = ident_mat;
211 qd.d_detJ = 1.;
212 d_quads.push_back(qd);
213 }
214
215 //
216 // third order quad points for triangle
217 //
218 if (d_quadOrder == 3) {
219 d_quads.clear();
220 // 1-d points are: {-sqrt{3}/sqrt{5}, 0, sqrt{3}/sqrt{5}}
221 // weights are: {5/9, 8/9, 5/9}
222 fe::QuadData qd;
223 qd.d_w = 5. / 9.;
224 qd.d_p = util::Point(-std::sqrt(3.) / std::sqrt(5.), 0., 0.);
225 qd.d_shapes = getShapes(qd.d_p);
226 qd.d_derShapes = getDerShapes(qd.d_p);
227 qd.d_J = ident_mat;
228 qd.d_detJ = 1.;
229 d_quads.push_back(qd);
230
231 qd.d_w = 8. / 9.;
232 qd.d_p = util::Point();
233 qd.d_shapes = getShapes(qd.d_p);
234 qd.d_derShapes = getDerShapes(qd.d_p);
235 qd.d_J = ident_mat;
236 qd.d_detJ = 1.;
237 d_quads.push_back(qd);
238
239 qd.d_w = 5. / 9.;
240 qd.d_p = util::Point(std::sqrt(3.) / std::sqrt(5.), 0., 0.);
241 qd.d_shapes = getShapes(qd.d_p);
242 qd.d_derShapes = getDerShapes(qd.d_p);
243 qd.d_J = ident_mat;
244 qd.d_detJ = 1.;
245 d_quads.push_back(qd);
246 }
247
248 //
249 // fourth order quad points for triangle
250 //
251 if (d_quadOrder == 4) {
252 d_quads.clear();
253 fe::QuadData qd;
254 qd.d_w = 0.6521451548625461;
255 qd.d_p = util::Point(-0.3399810435848563, 0., 0.);
256 qd.d_shapes = getShapes(qd.d_p);
257 qd.d_derShapes = getDerShapes(qd.d_p);
258 qd.d_J = ident_mat;
259 qd.d_detJ = 1.;
260 d_quads.push_back(qd);
261
262 qd.d_w = 0.6521451548625461;
263 qd.d_p = util::Point(0.3399810435848563, 0., 0.);
264 qd.d_shapes = getShapes(qd.d_p);
265 qd.d_derShapes = getDerShapes(qd.d_p);
266 qd.d_J = ident_mat;
267 qd.d_detJ = 1.;
268 d_quads.push_back(qd);
269
270 qd.d_w = 0.3478548451374538;
271 qd.d_p = util::Point(-0.8611363115940526, 0., 0.);
272 qd.d_shapes = getShapes(qd.d_p);
273 qd.d_derShapes = getDerShapes(qd.d_p);
274 qd.d_J = ident_mat;
275 qd.d_detJ = 1.;
276 d_quads.push_back(qd);
277
278 qd.d_w = 0.3478548451374538;
279 qd.d_p = util::Point(0.8611363115940526, 0., 0.);
280 qd.d_shapes = getShapes(qd.d_p);
281 qd.d_derShapes = getDerShapes(qd.d_p);
282 qd.d_J = ident_mat;
283 qd.d_detJ = 1.;
284 d_quads.push_back(qd);
285 }
286
287 //
288 // fifth order quad points for triangle
289 //
290 if (d_quadOrder == 5) {
291 d_quads.clear();
292 fe::QuadData qd;
293 qd.d_w = 0.5688888888888889;
294 qd.d_p = util::Point();
295 qd.d_shapes = getShapes(qd.d_p);
296 qd.d_derShapes = getDerShapes(qd.d_p);
297 qd.d_J = ident_mat;
298 qd.d_detJ = 1.;
299 d_quads.push_back(qd);
300
301 qd.d_w = 0.4786286704993665;
302 qd.d_p = util::Point(-0.5384693101056831, 0., 0.);
303 qd.d_shapes = getShapes(qd.d_p);
304 qd.d_derShapes = getDerShapes(qd.d_p);
305 qd.d_J = ident_mat;
306 qd.d_detJ = 1.;
307 d_quads.push_back(qd);
308
309 qd.d_w = 0.4786286704993665;
310 qd.d_p = util::Point(0.5384693101056831, 0., 0.);
311 qd.d_shapes = getShapes(qd.d_p);
312 qd.d_derShapes = getDerShapes(qd.d_p);
313 qd.d_J = ident_mat;
314 qd.d_detJ = 1.;
315 d_quads.push_back(qd);
316
317 qd.d_w = 0.2369268850561891;
318 qd.d_p = util::Point(-0.9061798459386640, 0., 0.);
319 qd.d_shapes = getShapes(qd.d_p);
320 qd.d_derShapes = getDerShapes(qd.d_p);
321 qd.d_J = ident_mat;
322 qd.d_detJ = 1.;
323 d_quads.push_back(qd);
324
325 qd.d_w = 0.2369268850561891;
326 qd.d_p = util::Point(0.9061798459386640, 0., 0.);
327 qd.d_shapes = getShapes(qd.d_p);
328 qd.d_derShapes = getDerShapes(qd.d_p);
329 qd.d_J = ident_mat;
330 qd.d_detJ = 1.;
331 d_quads.push_back(qd);
332 }
333}
A base class which provides methods to map points to/from reference element and to compute quadrature...
Definition baseElem.h:84
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
LineElem(size_t order)
Constructor for line element.
Definition lineElem.cpp:16
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
Collection of methods useful in simulation.
bool isGreater(const double &a, const double &b)
Returns true if a > b.
Definition function.cpp:15
bool isLess(const double &a, const double &b)
Returns true if a < b.
Definition function.cpp:20
A struct to store the quadrature data. List of data are.
Definition quadData.h:23
std::vector< double > d_shapes
Value of shape functions at quad point p.
Definition quadData.h:37
std::vector< std::vector< double > > d_derShapes
Value of derivative of shape functions at quad point p.
Definition quadData.h:50
double d_w
Quadrature weight.
Definition quadData.h:26
util::Point d_p
Quadrature point in 1-d, 2-d or 3-d.
Definition quadData.h:29
std::vector< std::vector< double > > d_J
Jacobian of the map from reference element to the element.
Definition quadData.h:58
double d_detJ
Determinant of the Jacobian of the map from reference element to the element.
Definition quadData.h:64
A structure to represent 3d vectors.
Definition point.h:30
double d_x
the x coordinate
Definition point.h:33