PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
fe::QuadElem Class Reference

A class for mapping and quadrature related operations for bi-linear quadrangle element. More...

#include <quadElem.h>

Inheritance diagram for fe::QuadElem:
Collaboration diagram for fe::QuadElem:

Public Member Functions

 QuadElem (size_t order)
 Constructor for quadrangle element.
 
double elemSize (const std::vector< util::Point > &nodes) override
 Returns the area of element.
 
std::vector< fe::QuadDatagetQuadDatas (const std::vector< util::Point > &nodes) override
 Get vector of quadrature data.
 
std::vector< fe::QuadDatagetQuadPoints (const std::vector< util::Point > &nodes) override
 Get vector of quadrature data.
 
- Public Member Functions inherited from fe::BaseElem
 BaseElem (size_t order, size_t element_type)
 Constructor.
 
size_t getElemType ()
 Get element type.
 
size_t getQuadOrder ()
 Get order of quadrature approximation.
 
size_t getNumQuadPoints ()
 Get number of quadrature points in the data.
 
virtual std::vector< double > getShapes (const util::Point &p, const std::vector< util::Point > &nodes)
 Returns the values of shape function at point p.
 
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.
 

Private Member Functions

std::vector< double > getShapes (const util::Point &p) override
 Returns the values of shape function at point p on reference element.
 
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.
 
double getJacobian (const util::Point &p, const std::vector< util::Point > &nodes, std::vector< std::vector< double > > *J) override
 Computes the Jacobian of map \( \Phi: T^0 \to T \).
 
void init () override
 Compute the quadrature points for quadrangle element.
 

Additional Inherited Members

- Protected Member Functions inherited from fe::BaseElem
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.
 
- Protected Attributes inherited from fe::BaseElem
size_t d_quadOrder
 Order of quadrature point integration approximation.
 
size_t d_numQuadPts
 Number of quadrature points for order d_quadOrder.
 
size_t d_elemType
 Element type.
 
std::vector< fe::QuadDatad_quads
 Quadrature data collection.
 

Detailed Description

A class for mapping and quadrature related operations for bi-linear quadrangle element.

The reference quadrangle element \( T^0 \) is given by vertices \( (-1, -1), \, (1,-1), \, (1,1), \, (-1,1) \).

  1. The shape functions at point \( (\xi, \eta ) \in T^0 \) are

    \[N^0_1(\xi, \eta) = \frac{(1- \xi)(1 - \eta)}{4}, \quad N^0_2(\xi, \eta) = \frac{(1+ \xi)(1 - \eta)}{4}, \]

    \[ N^0_3(\xi, \eta) = \frac{(1+ \xi)(1 + \eta)}{4},\quad N^0_4(\xi, \eta) = \frac{(1- \xi)(1 + \eta)}{4}. \]

  2. Derivative of shape functions at point \( (\xi, \eta ) \in T^0 \) are as follows

    \[\frac{d N^0_1(\xi, \eta)}{d\xi} = \frac{-(1 - \eta)}{4}, \quad \frac{d N^0_1(\xi, \eta)}{d\eta} = \frac{-(1 - \xi)}{4}, \]

    \[\frac{d N^0_2(\xi, \eta)}{d\xi} = \frac{(1 - \eta)}{4}, \quad \frac{d N^0_2 (\xi, \eta)}{d\eta} = \frac{-(1 + \xi)}{4}, \]

    \[\frac{d N^0_3(\xi, \eta)}{d\xi} = \frac{(1 + \eta)}{4}, \quad \frac{d N^0_3 (\xi, \eta)}{d\eta} = \frac{(1 + \xi)}{4}, \]

    \[\frac{d N^0_4(\xi, \eta)}{d\xi} = \frac{-(1 + \eta)}{4}, \quad \frac{d N^0_4 (\xi, \eta)}{d\eta} = \frac{(1 - \xi)}{4}. \]

  3. Map \( \Phi: T^0 \to T\) is given by

    \[ x(\xi, \eta) = \sum_{i=1}^4 N^0_i(\xi, \eta) v^i_x, \quad y(\xi, \eta) = \sum_{i=1}^4 N^0_i(\xi, \eta) v^i_y \]

    where \( v^1, v^2, v^3, v^4\) are vertices of element \( T \).
  4. Jacobian of the map \( \Phi: T^0 \to T\) is given by

    \[ J = \left[ { \begin{array}{cc} \frac{dx}{d\xi} &\frac{dy}{d\xi} \\ \frac{dx}{d\eta} & \frac{dy}{d\eta} \\ \end{array} } \right] \]

    and determinant of Jacobian is

    \[ det(J) = \frac{dx}{d\xi} \times \frac{dy}{d\eta} - \frac{dy}{d\xi}\times \frac{dx}{d\eta}. \]

Definition at line 64 of file quadElem.h.

Constructor & Destructor Documentation

◆ QuadElem()

fe::QuadElem::QuadElem ( size_t  order)
explicit

Constructor for quadrangle element.

Parameters
orderOrder of quadrature point approximation

Definition at line 14 of file quadElem.cpp.

16
17 // compute quad data
18 this->init();
19}
A base class which provides methods to map points to/from reference element and to compute quadrature...
Definition baseElem.h:84
void init() override
Compute the quadrature points for quadrangle element.
Definition quadElem.cpp:162
static const int vtk_type_quad
Integer flag for quad element.

References init().

Here is the call graph for this function:

Member Function Documentation

◆ elemSize()

double fe::QuadElem::elemSize ( const std::vector< util::Point > &  nodes)
overridevirtual

Returns the area of element.

If quadrangle \( T \) is given by points \( v^1, v^2, v^3, v^4 \) then the area is

\[ area(T) = \frac{(-v^1_1 + v^2_1 + v^3_1 - v^4_1) (-v^1_2 - v^2_2 + v^3_2 + v^4_2) - (-v^1_1 - v^2_1 + v^3_1 + v^4_1) (-v^1_2 + v^2_2 + v^3_2 - v^4_2)}{4}, \]

where \( v^i_1, v^i_2 \) are the x and y component of point \( v^i \).

Note that area and Jacobian of map \( \Phi: T^0 \to T \) are related as

\[ area(T) = area(T^0) \times det(J(\xi = 0, \eta = 0)), \]

where \( area(T^0) = 4 \) and \( J(\xi = 0, \eta = 0) \) is the Jacobian of map at point \( (0,0) \in T^0 \).

Parameters
nodesVertices of element
Returns
vector Vector of shape functions at point p

Implements fe::BaseElem.

Definition at line 21 of file quadElem.cpp.

21 {
22 return 0.25 *
23 ((-nodes[0].d_x + nodes[1].d_x + nodes[2].d_x - nodes[3].d_x) *
24 (-nodes[0].d_y - nodes[1].d_y + nodes[2].d_y + nodes[3].d_y) -
25 (-nodes[0].d_x - nodes[1].d_x + nodes[2].d_x + nodes[3].d_x) *
26 (-nodes[0].d_y + nodes[1].d_y + nodes[2].d_y - nodes[3].d_y));
27}

◆ getDerShapes()

std::vector< std::vector< double > > fe::QuadElem::getDerShapes ( const util::Point p)
overrideprivatevirtual

Returns the values of derivative of shape function at point p on reference element.

Parameters
pLocation of point
Returns
vector Vector of derivative of shape functions

Implements fe::BaseElem.

Definition at line 105 of file quadElem.cpp.

105 {
106
107 // N1 = (1 - xi)(1 - eta)/4
108 // --> d N1/d xi = -(1 - eta)/4, d N1/d eta = -(1 - xi)/4
109 //
110 // N2 = (1 + xi)(1 - eta)/4
111 // --> d N2/d xi = (1 - eta)/4, d N2/d eta = -(1 + xi)/4
112 //
113 // N3 = (1 + xi)(1 + eta)/4
114 // --> d N3/d xi = (1 + eta)/4, d N3/d eta = (1 + xi)/4
115 //
116 // N4 = (1 - xi)(1 + eta)/4
117 // --> d N4/d xi = -(1 + eta)/4, d N4/d eta = (1 - xi)/4
118 std::vector<std::vector<double>> r;
119 r.push_back(std::vector<double>{-0.25 * (1. - p.d_y), -0.25 * (1. - p.d_x)});
120 r.push_back(std::vector<double>{0.25 * (1. - p.d_y), -0.25 * (1. + p.d_x)});
121 r.push_back(std::vector<double>{0.25 * (1. + p.d_y), 0.25 * (1. + p.d_x)});
122 r.push_back(std::vector<double>{-0.25 * (1. + p.d_y), 0.25 * (1. - p.d_x)});
123
124 return r;
125}
double d_y
the y coordinate
Definition point.h:36
double d_x
the x coordinate
Definition point.h:33

References util::Point::d_x, and util::Point::d_y.

◆ getJacobian()

double fe::QuadElem::getJacobian ( const util::Point p,
const std::vector< util::Point > &  nodes,
std::vector< std::vector< double > > *  J 
)
overrideprivatevirtual

Computes the Jacobian of map \( \Phi: T^0 \to T \).

Parameters
pLocation of point in reference element
nodesVertices of element
JMatrix to store the Jacobian (if not nullptr)
Returns
det(J) Determinant of the Jacobain

Implements fe::BaseElem.

Definition at line 127 of file quadElem.cpp.

129 {
130
131 auto der_shapes = getDerShapes(p);
132 if (J != nullptr) {
133 J->resize(2);
134 (*J)[0] = std::vector<double>{
135 der_shapes[0][0] * nodes[0].d_x + der_shapes[1][0] * nodes[1].d_x +
136 der_shapes[2][0] * nodes[2].d_x + der_shapes[3][0] * nodes[3].d_x,
137 der_shapes[0][0] * nodes[0].d_y + der_shapes[1][0] * nodes[1].d_y +
138 der_shapes[2][0] * nodes[2].d_y + der_shapes[3][0] * nodes[3].d_y};
139 (*J)[1] = std::vector<double>{
140 der_shapes[0][1] * nodes[0].d_x + der_shapes[1][1] * nodes[1].d_x +
141 der_shapes[2][1] * nodes[2].d_x + der_shapes[3][1] * nodes[3].d_x,
142 der_shapes[0][1] * nodes[0].d_y + der_shapes[1][1] * nodes[1].d_y +
143 der_shapes[2][1] * nodes[2].d_y + der_shapes[3][1] * nodes[3].d_y};
144
145 return (*J)[0][0] * (*J)[1][1] - (*J)[0][1] * (*J)[1][0];
146 }
147
148 return (der_shapes[0][0] * nodes[0].d_x + der_shapes[1][0] * nodes[1].d_x +
149 der_shapes[2][0] * nodes[2].d_x + der_shapes[3][0] * nodes[3].d_x) *
150 (der_shapes[0][1] * nodes[0].d_y +
151 der_shapes[1][1] * nodes[1].d_y +
152 der_shapes[2][1] * nodes[2].d_y +
153 der_shapes[3][1] * nodes[3].d_y) -
154 (der_shapes[0][0] * nodes[0].d_y + der_shapes[1][0] * nodes[1].d_y +
155 der_shapes[2][0] * nodes[2].d_y + der_shapes[3][0] * nodes[3].d_y) *
156 (der_shapes[0][1] * nodes[0].d_x +
157 der_shapes[1][1] * nodes[1].d_x +
158 der_shapes[2][1] * nodes[2].d_x +
159 der_shapes[3][1] * nodes[3].d_x);
160}
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

◆ getQuadDatas()

std::vector< fe::QuadData > fe::QuadElem::getQuadDatas ( const std::vector< util::Point > &  nodes)
overridevirtual

Get vector of quadrature data.

Given element vertices, this method returns the list of quadrature point and essential quantities at quadrature points. Here, order of quadrature approximation is set in the constructor. List of data

  • quad point
  • quad weight
  • shape function evaluated at quad point
  • derivative of shape function evaluated at quad point
  • Jacobian matrix
  • determinant of the Jacobian

Let \( T \) is the given element with vertices \( v^1, v^2, v^3, v^4 \) and let \( T^0 \) is the reference element.

  1. To compute quadrature point, we first compute the quadrature points on \( T^0 \), and then map the point on \( T^0 \) to the point on \( T \) using the map \( \Phi: T^0 \to T\).
  2. To compute the quadrature weight, we compute the quadrature weight associated to the quadrature point in reference triangle \( T^0 \). Suppose \( w^0_q \) is the quadrature weight associated to quadrature point \( (\xi_q, \eta_q) \in T^0 \), then the quadrature point \( w_q \) associated to the mapped point \( (x(\xi_q, \eta_q), y(\xi_q, \eta_q)) \in T \) is given by

    \[ w_q = w^0_q * det(J(\xi_q, \eta_q)) \]

    where \( det(J) \) is the determinant of the Jacobian \( J \) of map \( \Phi: T^0 \to T \). \( J(\xi_q, \eta_q)\) is the value of \( J \) at point \( (\xi_q, \eta_q) \in T^0 \).
  3. We compute shape functions \( N_1, N_2, N_3, N_4\) associated to \( T \) at the quadrature point \( (x(\xi_q, \eta_q), y(\xi_q,\eta_q)) \) using formula

    \[ N_i(x(\xi_q, \eta_q), y(\xi_q, \eta_q)) = N^0_i(\xi_q, \eta_q). \]

  4. To compute the derivatives of shape functions such as \( \frac{\partial N_i}{\partial x}, \frac{\partial N_i}{\partial y}\) associated to \( T \), we use the relation between derivatives of shape function in \( T \) and \( T^0 \) described in fe::TriElem::getDerShapes.
Parameters
nodesVector of vertices of an element
Returns
vector Vector of QuadData

Implements fe::BaseElem.

Definition at line 32 of file quadElem.cpp.

32 {
33
34 // copy quad data associated to reference element
35 auto qds = d_quads;
36
37 // modify data
38 for (auto &qd : qds) {
39
40 // get Jacobian and determinant
41 qd.d_detJ = getJacobian(qd.d_p, nodes, &(qd.d_J));
42
43 // transform quad weight
44 qd.d_w *= qd.d_detJ;
45
46 // map point to triangle
47 qd.d_p.d_x = qd.d_shapes[0] * nodes[0].d_x + qd.d_shapes[1] * nodes[1].d_x +
48 qd.d_shapes[2] * nodes[2].d_x + qd.d_shapes[3] * nodes[3].d_x;
49 qd.d_p.d_y = qd.d_shapes[0] * nodes[0].d_y + qd.d_shapes[1] * nodes[1].d_y +
50 qd.d_shapes[2] * nodes[2].d_y + qd.d_shapes[3] * nodes[3].d_y;
51
52 // derivatives of shape function
53 std::vector<std::vector<double>> ders;
54
55 for (size_t i=3; i<4; i++) {
56 // partial N_i/ partial x
57 auto d1 = (qd.d_derShapes[i][0] * qd.d_J[1][1] - qd.d_derShapes[i][1] *
58 qd.d_J[0][1]) / qd.d_detJ;
59 // partial N_i/ partial y
60 auto d2 = (-qd.d_derShapes[i][0] * qd.d_J[1][0] + qd.d_derShapes[i][1] *
61 qd.d_J[0][0]) / qd.d_detJ;
62
63 ders.push_back(std::vector<double>{d1, d2});
64 }
65 qd.d_derShapes = ders;
66 }
67
68 return qds;
69}
std::vector< fe::QuadData > d_quads
Quadrature data collection.
Definition baseElem.h:252
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

◆ getQuadPoints()

std::vector< fe::QuadData > fe::QuadElem::getQuadPoints ( const std::vector< util::Point > &  nodes)
overridevirtual

Get vector of quadrature data.

Given element vertices, this method returns the list of quadrature point and essential quantities at quadrature points. Here, order of quadrature approximation is set in the constructor. List of data

  • quad point
  • quad weight
  • shape function evaluated at quad point

This function is lite version of QuadElem::getQuadDatas.

Parameters
nodesVector of vertices of an element
Returns
vector Vector of QuadData

Implements fe::BaseElem.

Definition at line 72 of file quadElem.cpp.

72 {
73
74 // copy quad data associated to reference element
75 auto qds = d_quads;
76
77 // modify data
78 for (auto &qd : qds) {
79
80 // transform quad weight
81 qd.d_w *= getJacobian(qd.d_p, nodes, nullptr);
82
83 // map point to triangle
84 qd.d_p.d_x = qd.d_shapes[0] * nodes[0].d_x + qd.d_shapes[1] * nodes[1].d_x +
85 qd.d_shapes[2] * nodes[2].d_x + qd.d_shapes[3] * nodes[3].d_x;
86 qd.d_p.d_y = qd.d_shapes[0] * nodes[0].d_y + qd.d_shapes[1] * nodes[1].d_y +
87 qd.d_shapes[2] * nodes[2].d_y + qd.d_shapes[3] * nodes[3].d_y;
88 }
89
90 return qds;
91}

◆ getShapes()

std::vector< double > fe::QuadElem::getShapes ( const util::Point p)
overrideprivatevirtual

Returns the values of shape function at point p on reference element.

Parameters
pLocation of point
Returns
vector Vector of shape functions at point p

Implements fe::BaseElem.

Definition at line 93 of file quadElem.cpp.

93 {
94
95 // N1 = (1 - xi)(1 - eta)/4
96 // N2 = (1 + xi)(1 - eta)/4
97 // N3 = (1 + xi)(1 + eta)/4
98 // N4 = (1 - xi)(1 + eta)/4
99 return std::vector<double>{
100 0.25 * (1. - p.d_x) * (1. - p.d_y), 0.25 * (1. + p.d_x) * (1. - p.d_y),
101 0.25 * (1. + p.d_x) * (1. + p.d_y), 0.25 * (1. - p.d_x) * (1. + p.d_y)};
102}

References util::Point::d_x, and util::Point::d_y.

◆ init()

void fe::QuadElem::init ( )
overrideprivatevirtual

Compute the quadrature points for quadrangle element.

Implements fe::BaseElem.

Definition at line 162 of file quadElem.cpp.

162 {
163
164 //
165 // compute quad data for reference quadrangle with vertex at
166 // p1 = (-1,-1), p2 = (1,-1), p3 = (1,1), p4 = (-1,1)
167 //
168 // Shape functions are
169 // N1 = (1 - xi)(1 - eta)/4
170 // N2 = (1 + xi)(1 - eta)/4
171 // N3 = (1 + xi)(1 + eta)/4
172 // N4 = (1 - xi)(1 + eta)/4
173 //
174 //
175 // Let [-1,1] is the 1-d reference element and {x1, x2, x3,.., xN} are N
176 // quad points for 1-d domain and {w1, w2, w3,..., wN} are respective
177 // weights. Then, the Nth order quad points in Quadrangle [-1,1]x[-1,1] is
178 // simply given by N^2 points and
179 //
180 // (i,j) point is (xi, xj) and weight is wi \times wj
181 //
182
183 if (!d_quads.empty())
184 return;
185
186 // no point in zeroth order
187 if (d_quadOrder == 0)
188 d_quads.resize(0);
189
190 // 2x2 identity matrix
191 std::vector<std::vector<double>> ident_mat;
192 ident_mat.push_back(std::vector<double>{1., 0.});
193 ident_mat.push_back(std::vector<double>{0., 1.});
194
195 //
196 // first order quad points
197 //
198 if (d_quadOrder == 1) {
199 d_quads.clear();
200 // 1-d points are: {0} and weights are: {2}
201 int npts = 1;
202 std::vector<double> x = std::vector<double>(1, 0.);
203 std::vector<double> w = std::vector<double>(1, 2.);
204 for (size_t i = 0; i < npts; i++)
205 for (size_t j = 0; j < npts; j++) {
206
207 fe::QuadData qd;
208 qd.d_w = w[i] * w[j];
209 qd.d_p = util::Point(x[i], x[j], 0.);
210 qd.d_shapes = getShapes(qd.d_p);
212 qd.d_J = ident_mat;
213 qd.d_detJ = 1.;
214 d_quads.push_back(qd);
215 }
216 }
217
218 //
219 // second order quad points
220 //
221 if (d_quadOrder == 2) {
222 d_quads.clear();
223 // 1-d points are: {-1/sqrt{3], 1/sqrt{3}} and weights are: {1,1}
224 int npts = 2;
225 std::vector<double> x =
226 std::vector<double>{-1. / std::sqrt(3.), 1. / std::sqrt(3.)};
227 std::vector<double> w = std::vector<double>{1., 1.};
228 for (size_t i = 0; i < npts; i++)
229 for (size_t j = 0; j < npts; j++) {
230
231 fe::QuadData qd;
232 qd.d_w = w[i] * w[j];
233 qd.d_p = util::Point(x[i], x[j], 0.);
234 qd.d_shapes = getShapes(qd.d_p);
236 qd.d_J = ident_mat;
237 qd.d_detJ = 1.;
238 d_quads.push_back(qd);
239 }
240 }
241
242 //
243 // third order quad points
244 //
245 if (d_quadOrder == 3) {
246 d_quads.clear();
247 // 1-d points are: {-sqrt{3}/sqrt{5}, 0, sqrt{3}/sqrt{5}}
248 // weights are: {5/9, 8/9, 5/9}
249 int npts = 3;
250 std::vector<double> x = std::vector<double>{
251 -std::sqrt(3.) / std::sqrt(5.), 0., std::sqrt(3.) / std::sqrt(5.)};
252 std::vector<double> w = std::vector<double>{5. / 9., 8. / 9., 5. / 9.};
253 for (size_t i = 0; i < npts; i++)
254 for (size_t j = 0; j < npts; j++) {
255
256 fe::QuadData qd;
257 qd.d_w = w[i] * w[j];
258 qd.d_p = util::Point(x[i], x[j], 0.);
259 qd.d_shapes = getShapes(qd.d_p);
261 qd.d_J = ident_mat;
262 qd.d_detJ = 1.;
263 d_quads.push_back(qd);
264 }
265 }
266
267 //
268 // fourth order quad points
269 //
270 if (d_quadOrder == 4) {
271 d_quads.clear();
272 int npts = 4;
273 std::vector<double> x =
274 std::vector<double>{-0.3399810435848563, 0.3399810435848563,
275 -0.8611363115940526, 0.8611363115940526};
276 std::vector<double> w =
277 std::vector<double>{0.6521451548625461, 0.6521451548625461,
278 0.3478548451374538, 0.3478548451374538};
279 for (size_t i = 0; i < npts; i++)
280 for (size_t j = 0; j < npts; j++) {
281
282 fe::QuadData qd;
283 qd.d_w = w[i] * w[j];
284 qd.d_p = util::Point(x[i], x[j], 0.);
285 qd.d_shapes = getShapes(qd.d_p);
287 qd.d_J = ident_mat;
288 qd.d_detJ = 1.;
289 d_quads.push_back(qd);
290 }
291 }
292
293 //
294 // fifth order quad points
295 //
296 if (d_quadOrder == 5) {
297 d_quads.clear();
298 int npts = 5;
299 std::vector<double> x =
300 std::vector<double>{0., -0.5384693101056831, 0.5384693101056831,
301 -0.9061798459386640, 0.9061798459386640};
302 std::vector<double> w = std::vector<double>{
303 0.5688888888888889, 0.4786286704993665, 0.4786286704993665,
304 0.2369268850561891, 0.2369268850561891};
305 for (size_t i = 0; i < npts; i++)
306 for (size_t j = 0; j < npts; j++) {
307
308 fe::QuadData qd;
309 qd.d_w = w[i] * w[j];
310 qd.d_p = util::Point(x[i], x[j], 0.);
311 qd.d_shapes = getShapes(qd.d_p);
313 qd.d_J = ident_mat;
314 qd.d_detJ = 1.;
315 d_quads.push_back(qd);
316 }
317 }
318}
size_t d_quadOrder
Order of quadrature point integration approximation.
Definition baseElem.h:243
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
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

References fe::QuadData::d_derShapes, fe::QuadData::d_detJ, fe::QuadData::d_J, fe::QuadData::d_p, fe::QuadData::d_shapes, and fe::QuadData::d_w.

Referenced by QuadElem().

Here is the caller graph for this function:

The documentation for this class was generated from the following files: