GooseFEM 1.4.1.dev2+g78f16df
Loading...
Searching...
No Matches
ElementQuad4Planar.h
Go to the documentation of this file.
1
12#ifndef GOOSEFEM_ELEMENTQUAD4PLANAR_H
13#define GOOSEFEM_ELEMENTQUAD4PLANAR_H
14
15#include "config.h"
16#include "detail.h"
17
18namespace GooseFEM {
19namespace Element {
20namespace Quad4 {
21
38class QuadraturePlanar : public QuadratureBaseCartesian<QuadraturePlanar> {
39public:
40 QuadraturePlanar() = default;
41
57 template <class T>
58 QuadraturePlanar(const T& x, double thick = 1.0)
59 : QuadraturePlanar(x, Gauss::xi(), Gauss::w(), thick)
60 {
61 }
62
80 template <class T, class X, class W>
81 QuadraturePlanar(const T& x, const X& xi, const W& w, double thick = 1.0)
82 {
83 m_x = x;
84 m_w = w;
85 m_xi = xi;
86 m_nip = w.size();
87 m_nelem = m_x.shape(0);
88 m_N = xt::empty<double>({m_nip, s_nne});
89 m_dNxi = xt::empty<double>({m_nip, s_nne, s_ndim});
90 m_thick = thick;
91
92 for (size_t q = 0; q < m_nip; ++q) {
93 m_N(q, 0) = 0.25 * (1.0 - xi(q, 0)) * (1.0 - xi(q, 1));
94 m_N(q, 1) = 0.25 * (1.0 + xi(q, 0)) * (1.0 - xi(q, 1));
95 m_N(q, 2) = 0.25 * (1.0 + xi(q, 0)) * (1.0 + xi(q, 1));
96 m_N(q, 3) = 0.25 * (1.0 - xi(q, 0)) * (1.0 + xi(q, 1));
97 }
98
99 for (size_t q = 0; q < m_nip; ++q) {
100 // - dN / dxi_0
101 m_dNxi(q, 0, 0) = -0.25 * (1.0 - xi(q, 1));
102 m_dNxi(q, 1, 0) = +0.25 * (1.0 - xi(q, 1));
103 m_dNxi(q, 2, 0) = +0.25 * (1.0 + xi(q, 1));
104 m_dNxi(q, 3, 0) = -0.25 * (1.0 + xi(q, 1));
105 // - dN / dxi_1
106 m_dNxi(q, 0, 1) = -0.25 * (1.0 - xi(q, 0));
107 m_dNxi(q, 1, 1) = -0.25 * (1.0 + xi(q, 0));
108 m_dNxi(q, 2, 1) = +0.25 * (1.0 + xi(q, 0));
109 m_dNxi(q, 3, 1) = +0.25 * (1.0 - xi(q, 0));
110 }
111
112 GOOSEFEM_ASSERT(m_x.shape(1) == s_nne);
113 GOOSEFEM_ASSERT(m_x.shape(2) == s_ndim);
114 GOOSEFEM_ASSERT(xt::has_shape(m_xi, {m_nip, s_ndim}));
115 GOOSEFEM_ASSERT(xt::has_shape(m_w, {m_nip}));
116 GOOSEFEM_ASSERT(xt::has_shape(m_N, {m_nip, s_nne}));
117 GOOSEFEM_ASSERT(xt::has_shape(m_dNxi, {m_nip, s_nne, s_ndim}));
118
119 m_dNx = xt::empty<double>({m_nelem, m_nip, s_nne, s_ndim});
120 m_vol = xt::empty<double>(this->shape_qscalar());
121
122 this->compute_dN_impl();
123 }
124
125private:
128
129 template <class T, class R>
130 void gradN_vector_impl(const T& elemvec, R& qtensor) const
131 {
132 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
133 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
134
135 qtensor.fill(0.0);
136
137#pragma omp parallel for
138 for (size_t e = 0; e < m_nelem; ++e) {
139
140 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
141
142 for (size_t q = 0; q < m_nip; ++q) {
143
144 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
145 auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
146
147 // gradu(i,j) += dNx(m,i) * u(m,j)
148 gradu(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
149 dNx(3, 0) * u(3, 0);
150 gradu(0, 1) = dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
151 dNx(3, 0) * u(3, 1);
152 gradu(1, 0) = dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) + dNx(2, 1) * u(2, 0) +
153 dNx(3, 1) * u(3, 0);
154 gradu(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
155 dNx(3, 1) * u(3, 1);
156 }
157 }
158 }
159
160 template <class T, class R>
161 void gradN_vector_T_impl(const T& elemvec, R& qtensor) const
162 {
163 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
164 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
165
166 qtensor.fill(0.0);
167
168#pragma omp parallel for
169 for (size_t e = 0; e < m_nelem; ++e) {
170
171 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
172
173 for (size_t q = 0; q < m_nip; ++q) {
174
175 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
176 auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
177
178 // gradu(j,i) += dNx(m,i) * u(m,j)
179 gradu(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
180 dNx(3, 0) * u(3, 0);
181 gradu(1, 0) = dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
182 dNx(3, 0) * u(3, 1);
183 gradu(0, 1) = dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) + dNx(2, 1) * u(2, 0) +
184 dNx(3, 1) * u(3, 0);
185 gradu(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
186 dNx(3, 1) * u(3, 1);
187 }
188 }
189 }
190
191 template <class T, class R>
192 void symGradN_vector_impl(const T& elemvec, R& qtensor) const
193 {
194 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
195 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
196
197 qtensor.fill(0.0);
198
199#pragma omp parallel for
200 for (size_t e = 0; e < m_nelem; ++e) {
201
202 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
203
204 for (size_t q = 0; q < m_nip; ++q) {
205
206 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
207 auto eps = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
208
209 // gradu(i,j) += dNx(m,i) * u(m,j)
210 // eps(j,i) = 0.5 * (gradu(i,j) + gradu(j,i))
211 eps(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
212 dNx(3, 0) * u(3, 0);
213 eps(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
214 dNx(3, 1) * u(3, 1);
215 eps(0, 1) = 0.5 * (dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
216 dNx(3, 0) * u(3, 1) + dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) +
217 dNx(2, 1) * u(2, 0) + dNx(3, 1) * u(3, 0));
218 eps(1, 0) = eps(0, 1);
219 }
220 }
221 }
222
223 template <class T, class R>
224 void int_N_scalar_NT_dV_impl(const T& qscalar, R& elemmat) const
225 {
226 GOOSEFEM_ASSERT(xt::has_shape(qscalar, this->shape_qscalar()));
227 GOOSEFEM_ASSERT(xt::has_shape(elemmat, this->shape_elemmat()));
228
229 elemmat.fill(0.0);
230
231#pragma omp parallel for
232 for (size_t e = 0; e < m_nelem; ++e) {
233
234 auto M = xt::adapt(&elemmat(e, 0, 0), xt::xshape<s_nne * s_ndim, s_nne * s_ndim>());
235
236 for (size_t q = 0; q < m_nip; ++q) {
237
238 auto N = xt::adapt(&m_N(q, 0), xt::xshape<s_nne>());
239 auto& vol = m_vol(e, q);
240 auto& rho = qscalar(e, q);
241
242 // M(m*ndim+i,n*ndim+i) += N(m) * scalar * N(n) * dV
243 for (size_t m = 0; m < s_nne; ++m) {
244 for (size_t n = 0; n < s_nne; ++n) {
245 M(m * s_ndim + 0, n * s_ndim + 0) += N(m) * rho * N(n) * vol;
246 M(m * s_ndim + 1, n * s_ndim + 1) += N(m) * rho * N(n) * vol;
247 }
248 }
249 }
250 }
251 }
252
253 template <class T, class R>
254 void int_gradN_dot_tensor2_dV_impl(const T& qtensor, R& elemvec) const
255 {
256 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
257 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
258
259 elemvec.fill(0.0);
260
261#pragma omp parallel for
262 for (size_t e = 0; e < m_nelem; ++e) {
263
264 auto f = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
265
266 for (size_t q = 0; q < m_nip; ++q) {
267
268 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
269 auto sig = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
270 auto& vol = m_vol(e, q);
271
272 for (size_t m = 0; m < s_nne; ++m) {
273 f(m, 0) += (dNx(m, 0) * sig(0, 0) + dNx(m, 1) * sig(1, 0)) * vol;
274 f(m, 1) += (dNx(m, 0) * sig(0, 1) + dNx(m, 1) * sig(1, 1)) * vol;
275 }
276 }
277 }
278 }
279
280 void compute_dN_impl()
281 {
282#pragma omp parallel
283 {
284 array_type::tensor<double, 2> J = xt::empty<double>({2, 2});
285 array_type::tensor<double, 2> Jinv = xt::empty<double>({2, 2});
286
287#pragma omp for
288 for (size_t e = 0; e < m_nelem; ++e) {
289
290 auto x = xt::adapt(&m_x(e, 0, 0), xt::xshape<s_nne, s_ndim>());
291
292 for (size_t q = 0; q < m_nip; ++q) {
293
294 auto dNxi = xt::adapt(&m_dNxi(q, 0, 0), xt::xshape<s_nne, s_ndim>());
295 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
296
297 // J(i,j) += dNxi(m,i) * x(m,j);
298 J(0, 0) = dNxi(0, 0) * x(0, 0) + dNxi(1, 0) * x(1, 0) + dNxi(2, 0) * x(2, 0) +
299 dNxi(3, 0) * x(3, 0);
300 J(0, 1) = dNxi(0, 0) * x(0, 1) + dNxi(1, 0) * x(1, 1) + dNxi(2, 0) * x(2, 1) +
301 dNxi(3, 0) * x(3, 1);
302 J(1, 0) = dNxi(0, 1) * x(0, 0) + dNxi(1, 1) * x(1, 0) + dNxi(2, 1) * x(2, 0) +
303 dNxi(3, 1) * x(3, 0);
304 J(1, 1) = dNxi(0, 1) * x(0, 1) + dNxi(1, 1) * x(1, 1) + dNxi(2, 1) * x(2, 1) +
305 dNxi(3, 1) * x(3, 1);
306
307 double Jdet = detail::tensor<2>::inv(J, Jinv);
308
309 // dNx(m,i) += Jinv(i,j) * dNxi(m,j);
310 for (size_t m = 0; m < s_nne; ++m) {
311 dNx(m, 0) = Jinv(0, 0) * dNxi(m, 0) + Jinv(0, 1) * dNxi(m, 1);
312 dNx(m, 1) = Jinv(1, 0) * dNxi(m, 0) + Jinv(1, 1) * dNxi(m, 1);
313 }
314
315 m_vol(e, q) = m_w(q) * Jdet * m_thick;
316 }
317 }
318 }
319 }
320
321 constexpr static size_t s_nne = 4;
322 constexpr static size_t s_ndim = 2;
323 constexpr static size_t s_tdim = 3;
324 size_t m_tdim = 3;
325 size_t m_nelem;
326 size_t m_nip;
327 array_type::tensor<double, 3> m_x;
328 array_type::tensor<double, 1> m_w;
329 array_type::tensor<double, 2> m_xi;
330 array_type::tensor<double, 2> m_N;
331 array_type::tensor<double, 3> m_dNxi;
332 array_type::tensor<double, 4> m_dNx;
333 array_type::tensor<double, 2> m_vol;
334 double m_thick;
335};
336
337} // namespace Quad4
338} // namespace Element
339} // namespace GooseFEM
340
341#endif
QuadraturePlanar(const T &x, const X &xi, const W &w, double thick=1.0)
Constructor with custom integration.
QuadraturePlanar(const T &x, double thick=1.0)
Constructor: use default Gauss integration.
CRTP base class for interpolation and quadrature for a generic element in Cartesian coordinates.
Definition Element.h:546
CRTP base class for quadrature.
Definition Element.h:151
auto shape_qscalar() const -> std::array< size_t, 2 >
Get the shape of a "qscalar" (a "qtensor" of rank 0)
Definition Element.h:354
auto shape_elemmat() const -> std::array< size_t, 3 >
Get the shape of an "elemmat".
Definition Element.h:273
auto shape_elemvec() const -> std::array< size_t, 3 >
Get the shape of an "elemvec".
Definition Element.h:252
auto AsTensor(const T &arg) const
Convert "qscalar" to "qtensor" of certain rank.
Definition Element.h:229
Basic configuration:
#define GOOSEFEM_ASSERT(expr)
All assertions are implementation as::
Definition config.h:97
@ Quad4
Quadrilateral: 4-noded element in 2-d.
Toolbox to perform finite element computations.
Definition Allocate.h:14