GooseFEM 1.4.1.dev2+g78f16df
Loading...
Searching...
No Matches
ElementQuad4.h
Go to the documentation of this file.
1
10#ifndef GOOSEFEM_ELEMENTQUAD4_H
11#define GOOSEFEM_ELEMENTQUAD4_H
12
13#include "Element.h"
14#include "config.h"
15#include "detail.h"
16
17namespace GooseFEM {
18namespace Element {
19
23namespace Quad4 {
24
36namespace Gauss {
37
45inline size_t nip()
46{
47 return 4;
48}
49
56{
57 size_t nip = 4;
58 size_t ndim = 2;
59
60 array_type::tensor<double, 2> xi = xt::empty<double>({nip, ndim});
61
62 xi(0, 0) = -1.0 / std::sqrt(3.0);
63 xi(0, 1) = -1.0 / std::sqrt(3.0);
64 xi(1, 0) = +1.0 / std::sqrt(3.0);
65 xi(1, 1) = -1.0 / std::sqrt(3.0);
66 xi(2, 0) = +1.0 / std::sqrt(3.0);
67 xi(2, 1) = +1.0 / std::sqrt(3.0);
68 xi(3, 0) = -1.0 / std::sqrt(3.0);
69 xi(3, 1) = +1.0 / std::sqrt(3.0);
70
71 return xi;
72}
73
80{
81 size_t nip = 4;
82
83 array_type::tensor<double, 1> w = xt::empty<double>({nip});
84
85 w(0) = 1.0;
86 w(1) = 1.0;
87 w(2) = 1.0;
88 w(3) = 1.0;
89
90 return w;
91}
92
93} // namespace Gauss
94
103namespace Nodal {
104
112inline size_t nip()
113{
114 return 4;
115}
116
123{
124 size_t nip = 4;
125 size_t ndim = 2;
126
127 array_type::tensor<double, 2> xi = xt::empty<double>({nip, ndim});
128
129 xi(0, 0) = -1.0;
130 xi(0, 1) = -1.0;
131
132 xi(1, 0) = +1.0;
133 xi(1, 1) = -1.0;
134
135 xi(2, 0) = +1.0;
136 xi(2, 1) = +1.0;
137
138 xi(3, 0) = -1.0;
139 xi(3, 1) = +1.0;
140
141 return xi;
142}
143
150{
151 size_t nip = 4;
152
153 array_type::tensor<double, 1> w = xt::empty<double>({nip});
154
155 w(0) = 1.0;
156 w(1) = 1.0;
157 w(2) = 1.0;
158 w(3) = 1.0;
159
160 return w;
161}
162
163} // namespace Nodal
164
174namespace MidPoint {
175
183inline size_t nip()
184{
185 return 1;
186}
187
194{
195 size_t nip = 1;
196 size_t ndim = 2;
197
198 array_type::tensor<double, 2> xi = xt::empty<double>({nip, ndim});
199
200 xi(0, 0) = 0.0;
201 xi(0, 1) = 0.0;
202
203 return xi;
204}
205
212{
213 size_t nip = 1;
214
215 array_type::tensor<double, 1> w = xt::empty<double>({nip});
216
217 w(0) = 1.0;
218
219 return w;
220}
221
222} // namespace MidPoint
223
237class Quadrature : public QuadratureBaseCartesian<Quadrature> {
238public:
239 Quadrature() = default;
240
255 template <class T>
256 Quadrature(const T& x) : Quadrature(x, Gauss::xi(), Gauss::w())
257 {
258 }
259
276 template <class T, class X, class W>
277 Quadrature(const T& x, const X& xi, const W& w)
278 {
279 m_x = x;
280 m_w = w;
281 m_xi = xi;
282 m_nip = w.size();
283 m_nelem = m_x.shape(0);
284 m_N = xt::empty<double>({m_nip, s_nne});
285 m_dNxi = xt::empty<double>({m_nip, s_nne, s_ndim});
286
287 for (size_t q = 0; q < m_nip; ++q) {
288 m_N(q, 0) = 0.25 * (1.0 - xi(q, 0)) * (1.0 - xi(q, 1));
289 m_N(q, 1) = 0.25 * (1.0 + xi(q, 0)) * (1.0 - xi(q, 1));
290 m_N(q, 2) = 0.25 * (1.0 + xi(q, 0)) * (1.0 + xi(q, 1));
291 m_N(q, 3) = 0.25 * (1.0 - xi(q, 0)) * (1.0 + xi(q, 1));
292 }
293
294 for (size_t q = 0; q < m_nip; ++q) {
295 // - dN / dxi_0
296 m_dNxi(q, 0, 0) = -0.25 * (1.0 - xi(q, 1));
297 m_dNxi(q, 1, 0) = +0.25 * (1.0 - xi(q, 1));
298 m_dNxi(q, 2, 0) = +0.25 * (1.0 + xi(q, 1));
299 m_dNxi(q, 3, 0) = -0.25 * (1.0 + xi(q, 1));
300 // - dN / dxi_1
301 m_dNxi(q, 0, 1) = -0.25 * (1.0 - xi(q, 0));
302 m_dNxi(q, 1, 1) = -0.25 * (1.0 + xi(q, 0));
303 m_dNxi(q, 2, 1) = +0.25 * (1.0 + xi(q, 0));
304 m_dNxi(q, 3, 1) = +0.25 * (1.0 - xi(q, 0));
305 }
306
307 GOOSEFEM_ASSERT(m_x.shape(1) == s_nne);
308 GOOSEFEM_ASSERT(m_x.shape(2) == s_ndim);
309 GOOSEFEM_ASSERT(xt::has_shape(m_xi, {m_nip, s_ndim}));
310 GOOSEFEM_ASSERT(xt::has_shape(m_w, {m_nip}));
311 GOOSEFEM_ASSERT(xt::has_shape(m_N, {m_nip, s_nne}));
312 GOOSEFEM_ASSERT(xt::has_shape(m_dNxi, {m_nip, s_nne, s_ndim}));
313
314 m_dNx = xt::empty<double>({m_nelem, m_nip, s_nne, s_ndim});
315 m_vol = xt::empty<double>(this->shape_qscalar());
316
317 this->compute_dN_impl();
318 }
319
320private:
323
324 template <class T, class R>
325 void interpQuad_vector_impl(const T& elemvec, R& qvector) const
326 {
327 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
328 GOOSEFEM_ASSERT(xt::has_shape(qvector, this->shape_qvector()));
329
330 qvector.fill(0.0);
331
332#pragma omp parallel for
333 for (size_t e = 0; e < m_nelem; ++e) {
334
335 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
336
337 for (size_t q = 0; q < m_nip; ++q) {
338
339 auto N = xt::adapt(&m_N(q, 0), xt::xshape<s_nne>());
340 auto ui = xt::adapt(&qvector(e, q, 0), xt::xshape<s_ndim>());
341
342 ui(0) = N(0) * u(0, 0) + N(1) * u(1, 0) + N(2) * u(2, 0) + N(3) * u(3, 0);
343 ui(1) = N(0) * u(0, 1) + N(1) * u(1, 1) + N(2) * u(2, 1) + N(3) * u(3, 1);
344 }
345 }
346 }
347
348 template <class T, class R>
349 void gradN_vector_impl(const T& elemvec, R& qtensor) const
350 {
351 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
352 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
353
354#pragma omp parallel for
355 for (size_t e = 0; e < m_nelem; ++e) {
356
357 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
358
359 for (size_t q = 0; q < m_nip; ++q) {
360
361 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
362 auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_ndim, s_ndim>());
363
364 // gradu(i,j) += dNx(m,i) * u(m,j)
365 gradu(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
366 dNx(3, 0) * u(3, 0);
367 gradu(0, 1) = dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
368 dNx(3, 0) * u(3, 1);
369 gradu(1, 0) = dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) + dNx(2, 1) * u(2, 0) +
370 dNx(3, 1) * u(3, 0);
371 gradu(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
372 dNx(3, 1) * u(3, 1);
373 }
374 }
375 }
376
377 template <class T, class R>
378 void gradN_vector_T_impl(const T& elemvec, R& qtensor) const
379 {
380 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
381 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
382
383#pragma omp parallel for
384 for (size_t e = 0; e < m_nelem; ++e) {
385
386 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
387
388 for (size_t q = 0; q < m_nip; ++q) {
389
390 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
391 auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_ndim, s_ndim>());
392
393 // gradu(j,i) += dNx(m,i) * u(m,j)
394 gradu(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
395 dNx(3, 0) * u(3, 0);
396 gradu(1, 0) = dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
397 dNx(3, 0) * u(3, 1);
398 gradu(0, 1) = dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) + dNx(2, 1) * u(2, 0) +
399 dNx(3, 1) * u(3, 0);
400 gradu(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
401 dNx(3, 1) * u(3, 1);
402 }
403 }
404 }
405
406 template <class T, class R>
407 void symGradN_vector_impl(const T& elemvec, R& qtensor) const
408 {
409 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
410 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
411
412#pragma omp parallel for
413 for (size_t e = 0; e < m_nelem; ++e) {
414
415 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
416
417 for (size_t q = 0; q < m_nip; ++q) {
418
419 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
420 auto eps = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_ndim, s_ndim>());
421
422 // gradu(i,j) += dNx(m,i) * u(m,j)
423 // eps(j,i) = 0.5 * (gradu(i,j) + gradu(j,i))
424 eps(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
425 dNx(3, 0) * u(3, 0);
426 eps(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
427 dNx(3, 1) * u(3, 1);
428 eps(0, 1) = 0.5 * (dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
429 dNx(3, 0) * u(3, 1) + dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) +
430 dNx(2, 1) * u(2, 0) + dNx(3, 1) * u(3, 0));
431 eps(1, 0) = eps(0, 1);
432 }
433 }
434 }
435
436 template <class T, class R>
437 void int_N_scalar_NT_dV_impl(const T& qscalar, R& elemmat) const
438 {
439 GOOSEFEM_ASSERT(xt::has_shape(qscalar, this->shape_qscalar()));
440 GOOSEFEM_ASSERT(xt::has_shape(elemmat, this->shape_elemmat()));
441
442 elemmat.fill(0.0);
443
444#pragma omp parallel for
445 for (size_t e = 0; e < m_nelem; ++e) {
446
447 auto M = xt::adapt(&elemmat(e, 0, 0), xt::xshape<s_nne * s_ndim, s_nne * s_ndim>());
448
449 for (size_t q = 0; q < m_nip; ++q) {
450
451 auto N = xt::adapt(&m_N(q, 0), xt::xshape<s_nne>());
452 auto& vol = m_vol(e, q);
453 auto& rho = qscalar(e, q);
454
455 // M(m*ndim+i,n*ndim+i) += N(m) * scalar * N(n) * dV
456 for (size_t m = 0; m < s_nne; ++m) {
457 for (size_t n = 0; n < s_nne; ++n) {
458 M(m * s_ndim + 0, n * s_ndim + 0) += N(m) * rho * N(n) * vol;
459 M(m * s_ndim + 1, n * s_ndim + 1) += N(m) * rho * N(n) * vol;
460 }
461 }
462 }
463 }
464 }
465
466 template <class T, class R>
467 void int_gradN_dot_tensor2_dV_impl(const T& qtensor, R& elemvec) const
468 {
469 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
470 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
471
472 elemvec.fill(0.0);
473
474#pragma omp parallel for
475 for (size_t e = 0; e < m_nelem; ++e) {
476
477 auto f = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
478
479 for (size_t q = 0; q < m_nip; ++q) {
480
481 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
482 auto sig = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_ndim, s_ndim>());
483 auto& vol = m_vol(e, q);
484
485 for (size_t m = 0; m < s_nne; ++m) {
486 f(m, 0) += (dNx(m, 0) * sig(0, 0) + dNx(m, 1) * sig(1, 0)) * vol;
487 f(m, 1) += (dNx(m, 0) * sig(0, 1) + dNx(m, 1) * sig(1, 1)) * vol;
488 }
489 }
490 }
491 }
492
493 void compute_dN_impl()
494 {
495#pragma omp parallel
496 {
497 array_type::tensor<double, 2> J = xt::empty<double>({2, 2});
498 array_type::tensor<double, 2> Jinv = xt::empty<double>({2, 2});
499
500#pragma omp for
501 for (size_t e = 0; e < m_nelem; ++e) {
502
503 auto x = xt::adapt(&m_x(e, 0, 0), xt::xshape<s_nne, s_ndim>());
504
505 for (size_t q = 0; q < m_nip; ++q) {
506
507 auto dNxi = xt::adapt(&m_dNxi(q, 0, 0), xt::xshape<s_nne, s_ndim>());
508 auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<s_nne, s_ndim>());
509
510 // J(i,j) += dNxi(m,i) * x(m,j);
511 J(0, 0) = dNxi(0, 0) * x(0, 0) + dNxi(1, 0) * x(1, 0) + dNxi(2, 0) * x(2, 0) +
512 dNxi(3, 0) * x(3, 0);
513 J(0, 1) = dNxi(0, 0) * x(0, 1) + dNxi(1, 0) * x(1, 1) + dNxi(2, 0) * x(2, 1) +
514 dNxi(3, 0) * x(3, 1);
515 J(1, 0) = dNxi(0, 1) * x(0, 0) + dNxi(1, 1) * x(1, 0) + dNxi(2, 1) * x(2, 0) +
516 dNxi(3, 1) * x(3, 0);
517 J(1, 1) = dNxi(0, 1) * x(0, 1) + dNxi(1, 1) * x(1, 1) + dNxi(2, 1) * x(2, 1) +
518 dNxi(3, 1) * x(3, 1);
519
520 double Jdet = detail::tensor<2>::inv(J, Jinv);
521
522 // dNx(m,i) += Jinv(i,j) * dNxi(m,j);
523 for (size_t m = 0; m < s_nne; ++m) {
524 dNx(m, 0) = Jinv(0, 0) * dNxi(m, 0) + Jinv(0, 1) * dNxi(m, 1);
525 dNx(m, 1) = Jinv(1, 0) * dNxi(m, 0) + Jinv(1, 1) * dNxi(m, 1);
526 }
527
528 m_vol(e, q) = m_w(q) * Jdet;
529 }
530 }
531 }
532 }
533
534 constexpr static size_t s_nne = 4;
535 constexpr static size_t s_ndim = 2;
536 constexpr static size_t s_tdim = 2;
537 size_t m_tdim = 2;
538 size_t m_nelem;
539 size_t m_nip;
540 array_type::tensor<double, 3> m_x;
541 array_type::tensor<double, 1> m_w;
542 array_type::tensor<double, 2> m_xi;
543 array_type::tensor<double, 2> m_N;
544 array_type::tensor<double, 3> m_dNxi;
545 array_type::tensor<double, 4> m_dNx;
546 array_type::tensor<double, 2> m_vol;
547};
548
549} // namespace Quad4
550} // namespace Element
551} // namespace GooseFEM
552
553#endif
Convenience methods for integration point data.
Interpolation and quadrature.
Quadrature(const T &x, const X &xi, const W &w)
Constructor with custom integration.
Quadrature(const T &x)
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_qvector() const -> std::array< size_t, 3 >
Get the shape of a "qvector" (a "qtensor" of rank 1)
Definition Element.h:363
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
array_type::tensor< double, 1 > w()
Integration point weights.
size_t nip()
Number of integration points:
array_type::tensor< double, 2 > xi()
Integration point coordinates (local coordinates).
size_t nip()
Number of integration points::
array_type::tensor< double, 2 > xi()
Integration point coordinates (local coordinates).
array_type::tensor< double, 1 > w()
Integration point weights.
size_t nip()
Number of integration points::
array_type::tensor< double, 1 > w()
Integration point weights.
array_type::tensor< double, 2 > xi()
Integration point coordinates (local coordinates).
xt::xtensor< T, N > tensor
Fixed (static) rank array.
Definition config.h:177
Toolbox to perform finite element computations.
Definition Allocate.h:14