GooseFEM 1.4.1.dev2+g78f16df
Loading...
Searching...
No Matches
ElementQuad4Axisymmetric.h
Go to the documentation of this file.
1
10#ifndef GOOSEFEM_ELEMENTQUAD4AXISYMMETRIC_H
11#define GOOSEFEM_ELEMENTQUAD4AXISYMMETRIC_H
12
13#include "config.h"
14
15namespace GooseFEM {
16namespace Element {
17namespace Quad4 {
18
33class QuadratureAxisymmetric : public QuadratureBaseCartesian<QuadratureAxisymmetric> {
34public:
35 QuadratureAxisymmetric() = default;
36
51 template <class T>
52 QuadratureAxisymmetric(const T& x) : QuadratureAxisymmetric(x, Gauss::xi(), Gauss::w())
53 {
54 }
55
72 template <class T, class X, class W>
73 QuadratureAxisymmetric(const T& x, const X& xi, const W& w)
74 {
75 m_x = x;
76 m_w = w;
77 m_xi = xi;
78 m_nip = w.size();
79 m_nelem = m_x.shape(0);
80
81 GOOSEFEM_ASSERT(m_x.shape(1) == s_nne);
82 GOOSEFEM_ASSERT(m_x.shape(2) == s_ndim);
83 GOOSEFEM_ASSERT(xt::has_shape(m_xi, {m_nip, s_ndim}));
84 GOOSEFEM_ASSERT(xt::has_shape(m_w, {m_nip}));
85
86 m_N = xt::empty<double>({m_nip, s_nne});
87 m_dNxi = xt::empty<double>({m_nip, s_nne, s_ndim});
88 m_B = xt::empty<double>({m_nelem, m_nip, s_nne, s_tdim, s_tdim, s_tdim});
89 m_vol = xt::empty<double>(this->shape_qscalar());
90
91 for (size_t q = 0; q < m_nip; ++q) {
92 m_N(q, 0) = 0.25 * (1.0 - m_xi(q, 0)) * (1.0 - m_xi(q, 1));
93 m_N(q, 1) = 0.25 * (1.0 + m_xi(q, 0)) * (1.0 - m_xi(q, 1));
94 m_N(q, 2) = 0.25 * (1.0 + m_xi(q, 0)) * (1.0 + m_xi(q, 1));
95 m_N(q, 3) = 0.25 * (1.0 - m_xi(q, 0)) * (1.0 + m_xi(q, 1));
96 }
97
98 for (size_t q = 0; q < m_nip; ++q) {
99 // - dN / dxi_0
100 m_dNxi(q, 0, 0) = -0.25 * (1.0 - m_xi(q, 1));
101 m_dNxi(q, 1, 0) = +0.25 * (1.0 - m_xi(q, 1));
102 m_dNxi(q, 2, 0) = +0.25 * (1.0 + m_xi(q, 1));
103 m_dNxi(q, 3, 0) = -0.25 * (1.0 + m_xi(q, 1));
104 // - dN / dxi_1
105 m_dNxi(q, 0, 1) = -0.25 * (1.0 - m_xi(q, 0));
106 m_dNxi(q, 1, 1) = -0.25 * (1.0 + m_xi(q, 0));
107 m_dNxi(q, 2, 1) = +0.25 * (1.0 + m_xi(q, 0));
108 m_dNxi(q, 3, 1) = +0.25 * (1.0 - m_xi(q, 0));
109 }
110
111 this->compute_dN_impl();
112 }
113
123 {
124 return m_B;
125 }
126
127private:
130
131 // qtensor(e, q, i, j) += B(e, q, m, i, j, k) * elemvec(e, q, m, k)
132 template <class T, class R>
133 void gradN_vector_impl(const T& elemvec, R& qtensor) const
134 {
135 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
136 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
137
138 qtensor.fill(0.0);
139
140#pragma omp parallel for
141 for (size_t e = 0; e < m_nelem; ++e) {
142
143 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
144
145 for (size_t q = 0; q < m_nip; ++q) {
146
147 auto B =
148 xt::adapt(&m_B(e, q, 0, 0, 0, 0), xt::xshape<s_nne, s_tdim, s_tdim, s_tdim>());
149 auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
150
151 // gradu(i,j) += B(m,i,j,k) * u(m,perm(k))
152 // (where perm(0) = 1, perm(2) = 0)
153 gradu(0, 0) = B(0, 0, 0, 0) * u(0, 1) + B(1, 0, 0, 0) * u(1, 1) +
154 B(2, 0, 0, 0) * u(2, 1) + B(3, 0, 0, 0) * u(3, 1);
155 gradu(1, 1) = B(0, 1, 1, 0) * u(0, 1) + B(1, 1, 1, 0) * u(1, 1) +
156 B(2, 1, 1, 0) * u(2, 1) + B(3, 1, 1, 0) * u(3, 1);
157 gradu(2, 2) = B(0, 2, 2, 2) * u(0, 0) + B(1, 2, 2, 2) * u(1, 0) +
158 B(2, 2, 2, 2) * u(2, 0) + B(3, 2, 2, 2) * u(3, 0);
159 gradu(0, 2) = B(0, 0, 2, 2) * u(0, 0) + B(1, 0, 2, 2) * u(1, 0) +
160 B(2, 0, 2, 2) * u(2, 0) + B(3, 0, 2, 2) * u(3, 0);
161 gradu(2, 0) = B(0, 2, 0, 0) * u(0, 1) + B(1, 2, 0, 0) * u(1, 1) +
162 B(2, 2, 0, 0) * u(2, 1) + B(3, 2, 0, 0) * u(3, 1);
163 }
164 }
165 }
166
167 template <class T, class R>
168 void gradN_vector_T_impl(const T& elemvec, R& qtensor) const
169 {
170 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
171 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
172
173 qtensor.fill(0.0);
174
175#pragma omp parallel for
176 for (size_t e = 0; e < m_nelem; ++e) {
177
178 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
179
180 for (size_t q = 0; q < m_nip; ++q) {
181
182 auto B =
183 xt::adapt(&m_B(e, q, 0, 0, 0, 0), xt::xshape<s_nne, s_tdim, s_tdim, s_tdim>());
184 auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
185
186 // gradu(j,i) += B(m,i,j,k) * u(m,perm(k))
187 // (where perm(0) = 1, perm(2) = 0)
188 gradu(0, 0) = B(0, 0, 0, 0) * u(0, 1) + B(1, 0, 0, 0) * u(1, 1) +
189 B(2, 0, 0, 0) * u(2, 1) + B(3, 0, 0, 0) * u(3, 1);
190 gradu(1, 1) = B(0, 1, 1, 0) * u(0, 1) + B(1, 1, 1, 0) * u(1, 1) +
191 B(2, 1, 1, 0) * u(2, 1) + B(3, 1, 1, 0) * u(3, 1);
192 gradu(2, 2) = B(0, 2, 2, 2) * u(0, 0) + B(1, 2, 2, 2) * u(1, 0) +
193 B(2, 2, 2, 2) * u(2, 0) + B(3, 2, 2, 2) * u(3, 0);
194 gradu(2, 0) = B(0, 0, 2, 2) * u(0, 0) + B(1, 0, 2, 2) * u(1, 0) +
195 B(2, 0, 2, 2) * u(2, 0) + B(3, 0, 2, 2) * u(3, 0);
196 gradu(0, 2) = B(0, 2, 0, 0) * u(0, 1) + B(1, 2, 0, 0) * u(1, 1) +
197 B(2, 2, 0, 0) * u(2, 1) + B(3, 2, 0, 0) * u(3, 1);
198 }
199 }
200 }
201
202 template <class T, class R>
203 void symGradN_vector_impl(const T& elemvec, R& qtensor) const
204 {
205 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
206 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
207
208 qtensor.fill(0.0);
209
210#pragma omp parallel for
211 for (size_t e = 0; e < m_nelem; ++e) {
212
213 auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
214
215 for (size_t q = 0; q < m_nip; ++q) {
216
217 auto B =
218 xt::adapt(&m_B(e, q, 0, 0, 0, 0), xt::xshape<s_nne, s_tdim, s_tdim, s_tdim>());
219 auto eps = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
220
221 // gradu(j,i) += B(m,i,j,k) * u(m,perm(k))
222 // eps(j,i) = 0.5 * (gradu(i,j) + gradu(j,i))
223 // (where perm(0) = 1, perm(2) = 0)
224 eps(0, 0) = B(0, 0, 0, 0) * u(0, 1) + B(1, 0, 0, 0) * u(1, 1) +
225 B(2, 0, 0, 0) * u(2, 1) + B(3, 0, 0, 0) * u(3, 1);
226 eps(1, 1) = B(0, 1, 1, 0) * u(0, 1) + B(1, 1, 1, 0) * u(1, 1) +
227 B(2, 1, 1, 0) * u(2, 1) + B(3, 1, 1, 0) * u(3, 1);
228 eps(2, 2) = B(0, 2, 2, 2) * u(0, 0) + B(1, 2, 2, 2) * u(1, 0) +
229 B(2, 2, 2, 2) * u(2, 0) + B(3, 2, 2, 2) * u(3, 0);
230 eps(2, 0) = 0.5 * (B(0, 0, 2, 2) * u(0, 0) + B(1, 0, 2, 2) * u(1, 0) +
231 B(2, 0, 2, 2) * u(2, 0) + B(3, 0, 2, 2) * u(3, 0) +
232 B(0, 2, 0, 0) * u(0, 1) + B(1, 2, 0, 0) * u(1, 1) +
233 B(2, 2, 0, 0) * u(2, 1) + B(3, 2, 0, 0) * u(3, 1));
234 eps(0, 2) = eps(2, 0);
235 }
236 }
237 }
238
239 // elemmat(e, q, m * ndim + i, n * ndim + i) +=
240 // N(e, q, m) * qscalar(e, q) * N(e, q, n) * dV(e, q)
241 template <class T, class R>
242 void int_N_scalar_NT_dV_impl(const T& qscalar, R& elemmat) const
243 {
244 GOOSEFEM_ASSERT(xt::has_shape(qscalar, this->shape_qscalar()));
245 GOOSEFEM_ASSERT(xt::has_shape(elemmat, this->shape_elemmat()));
246
247 elemmat.fill(0.0);
248
249#pragma omp parallel for
250 for (size_t e = 0; e < m_nelem; ++e) {
251
252 auto M = xt::adapt(&elemmat(e, 0, 0), xt::xshape<s_nne * s_ndim, s_nne * s_ndim>());
253
254 for (size_t q = 0; q < m_nip; ++q) {
255
256 auto N = xt::adapt(&m_N(q, 0), xt::xshape<s_nne>());
257 auto& vol = m_vol(e, q);
258 auto& rho = qscalar(e, q);
259
260 // M(m*ndim+i,n*ndim+i) += N(m) * scalar * N(n) * dV
261 for (size_t m = 0; m < s_nne; ++m) {
262 for (size_t n = 0; n < s_nne; ++n) {
263 M(m * s_ndim + 0, n * s_ndim + 0) += N(m) * rho * N(n) * vol;
264 M(m * s_ndim + 1, n * s_ndim + 1) += N(m) * rho * N(n) * vol;
265 }
266 }
267 }
268 }
269 }
270
271 // fm = ( Bm^T : qtensor ) dV
272 template <class T, class R>
273 void int_gradN_dot_tensor2_dV_impl(const T& qtensor, R& elemvec) const
274 {
275 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<2>()));
276 GOOSEFEM_ASSERT(xt::has_shape(elemvec, this->shape_elemvec()));
277
278 elemvec.fill(0.0);
279
280#pragma omp parallel for
281 for (size_t e = 0; e < m_nelem; ++e) {
282
283 auto f = xt::adapt(&elemvec(e, 0, 0), xt::xshape<s_nne, s_ndim>());
284
285 for (size_t q = 0; q < m_nip; ++q) {
286
287 auto B =
288 xt::adapt(&m_B(e, q, 0, 0, 0, 0), xt::xshape<s_nne, s_tdim, s_tdim, s_tdim>());
289 auto sig = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<s_tdim, s_tdim>());
290 auto& vol = m_vol(e, q);
291
292 // f(m,i) += B(m,i,j,perm(k)) * sig(i,j) * dV
293 // (where perm(0) = 1, perm(2) = 0)
294 for (size_t m = 0; m < s_nne; ++m) {
295 f(m, 0) += vol * (B(m, 2, 2, 2) * sig(2, 2) + B(m, 0, 2, 2) * sig(0, 2));
296 f(m, 1) += vol * (B(m, 0, 0, 0) * sig(0, 0) + B(m, 1, 1, 0) * sig(1, 1) +
297 B(m, 2, 0, 0) * sig(2, 0));
298 }
299 }
300 }
301 }
302
303 // Kmn = ( Bm^T : qtensor : Bn ) dV
304 template <class T, class R>
305 void int_gradN_dot_tensor4_dot_gradNT_dV_impl(const T& qtensor, R& elemmat) const
306 {
307 GOOSEFEM_ASSERT(xt::has_shape(qtensor, this->shape_qtensor<4>()));
308 GOOSEFEM_ASSERT(xt::has_shape(elemmat, this->shape_elemmat()));
309
310 elemmat.fill(0.0);
311
312#pragma omp parallel for
313 for (size_t e = 0; e < m_nelem; ++e) {
314
315 auto K = xt::adapt(&elemmat(e, 0, 0), xt::xshape<s_nne * s_ndim, s_nne * s_ndim>());
316
317 for (size_t q = 0; q < m_nip; ++q) {
318
319 auto B =
320 xt::adapt(&m_B(e, q, 0, 0, 0, 0), xt::xshape<s_nne, s_tdim, s_tdim, s_tdim>());
321 auto C = xt::adapt(
322 &qtensor(e, q, 0, 0, 0, 0), xt::xshape<s_tdim, s_tdim, s_tdim, s_tdim>()
323 );
324 auto& vol = m_vol(e, q);
325
326 // K(m*s_ndim+perm(c), n*s_ndim+perm(f)) = B(m,a,b,c) * C(a,b,d,e) * B(n,e,d,f) *
327 // vol; (where perm(0) = 1, perm(2) = 0)
328 for (size_t m = 0; m < s_nne; ++m) {
329 for (size_t n = 0; n < s_nne; ++n) {
330 K(m * s_ndim + 1, n * s_ndim + 1) +=
331 B(m, 0, 0, 0) * C(0, 0, 0, 0) * B(n, 0, 0, 0) * vol;
332 K(m * s_ndim + 1, n * s_ndim + 1) +=
333 B(m, 0, 0, 0) * C(0, 0, 1, 1) * B(n, 1, 1, 0) * vol;
334 K(m * s_ndim + 1, n * s_ndim + 0) +=
335 B(m, 0, 0, 0) * C(0, 0, 2, 2) * B(n, 2, 2, 2) * vol;
336 K(m * s_ndim + 1, n * s_ndim + 0) +=
337 B(m, 0, 0, 0) * C(0, 0, 2, 0) * B(n, 0, 2, 2) * vol;
338 K(m * s_ndim + 1, n * s_ndim + 1) +=
339 B(m, 0, 0, 0) * C(0, 0, 0, 2) * B(n, 2, 0, 0) * vol;
340
341 K(m * s_ndim + 1, n * s_ndim + 1) +=
342 B(m, 1, 1, 0) * C(1, 1, 0, 0) * B(n, 0, 0, 0) * vol;
343 K(m * s_ndim + 1, n * s_ndim + 1) +=
344 B(m, 1, 1, 0) * C(1, 1, 1, 1) * B(n, 1, 1, 0) * vol;
345 K(m * s_ndim + 1, n * s_ndim + 0) +=
346 B(m, 1, 1, 0) * C(1, 1, 2, 2) * B(n, 2, 2, 2) * vol;
347 K(m * s_ndim + 1, n * s_ndim + 0) +=
348 B(m, 1, 1, 0) * C(1, 1, 2, 0) * B(n, 0, 2, 2) * vol;
349 K(m * s_ndim + 1, n * s_ndim + 1) +=
350 B(m, 1, 1, 0) * C(1, 1, 0, 2) * B(n, 2, 0, 0) * vol;
351
352 K(m * s_ndim + 0, n * s_ndim + 1) +=
353 B(m, 2, 2, 2) * C(2, 2, 0, 0) * B(n, 0, 0, 0) * vol;
354 K(m * s_ndim + 0, n * s_ndim + 1) +=
355 B(m, 2, 2, 2) * C(2, 2, 1, 1) * B(n, 1, 1, 0) * vol;
356 K(m * s_ndim + 0, n * s_ndim + 0) +=
357 B(m, 2, 2, 2) * C(2, 2, 2, 2) * B(n, 2, 2, 2) * vol;
358 K(m * s_ndim + 0, n * s_ndim + 0) +=
359 B(m, 2, 2, 2) * C(2, 2, 2, 0) * B(n, 0, 2, 2) * vol;
360 K(m * s_ndim + 0, n * s_ndim + 1) +=
361 B(m, 2, 2, 2) * C(2, 2, 0, 2) * B(n, 2, 0, 0) * vol;
362
363 K(m * s_ndim + 0, n * s_ndim + 1) +=
364 B(m, 0, 2, 2) * C(0, 2, 0, 0) * B(n, 0, 0, 0) * vol;
365 K(m * s_ndim + 0, n * s_ndim + 1) +=
366 B(m, 0, 2, 2) * C(0, 2, 1, 1) * B(n, 1, 1, 0) * vol;
367 K(m * s_ndim + 0, n * s_ndim + 0) +=
368 B(m, 0, 2, 2) * C(0, 2, 2, 2) * B(n, 2, 2, 2) * vol;
369 K(m * s_ndim + 0, n * s_ndim + 0) +=
370 B(m, 0, 2, 2) * C(0, 2, 2, 0) * B(n, 0, 2, 2) * vol;
371 K(m * s_ndim + 0, n * s_ndim + 1) +=
372 B(m, 0, 2, 2) * C(0, 2, 0, 2) * B(n, 2, 0, 0) * vol;
373
374 K(m * s_ndim + 1, n * s_ndim + 1) +=
375 B(m, 2, 0, 0) * C(2, 0, 0, 0) * B(n, 0, 0, 0) * vol;
376 K(m * s_ndim + 1, n * s_ndim + 1) +=
377 B(m, 2, 0, 0) * C(2, 0, 1, 1) * B(n, 1, 1, 0) * vol;
378 K(m * s_ndim + 1, n * s_ndim + 0) +=
379 B(m, 2, 0, 0) * C(2, 0, 2, 2) * B(n, 2, 2, 2) * vol;
380 K(m * s_ndim + 1, n * s_ndim + 0) +=
381 B(m, 2, 0, 0) * C(2, 0, 2, 0) * B(n, 0, 2, 2) * vol;
382 K(m * s_ndim + 1, n * s_ndim + 1) +=
383 B(m, 2, 0, 0) * C(2, 0, 0, 2) * B(n, 2, 0, 0) * vol;
384 }
385 }
386 }
387 }
388 }
389
390 void compute_dN_impl()
391 {
392 // most components remain zero, and are not written
393 m_B.fill(0.0);
394
395#pragma omp parallel
396 {
397 array_type::tensor<double, 2> J = xt::empty<double>({2, 2});
398 array_type::tensor<double, 2> Jinv = xt::empty<double>({2, 2});
399
400#pragma omp for
401 for (size_t e = 0; e < m_nelem; ++e) {
402
403 auto x = xt::adapt(&m_x(e, 0, 0), xt::xshape<s_nne, s_ndim>());
404
405 for (size_t q = 0; q < m_nip; ++q) {
406
407 auto dNxi = xt::adapt(&m_dNxi(q, 0, 0), xt::xshape<s_nne, s_ndim>());
408 auto B = xt::adapt(
409 &m_B(e, q, 0, 0, 0, 0), xt::xshape<s_nne, s_tdim, s_tdim, s_tdim>()
410 );
411 auto N = xt::adapt(&m_N(q, 0), xt::xshape<s_nne>());
412
413 // J(i,j) += dNxi(m,i) * x(m,j);
414 J(0, 0) = dNxi(0, 0) * x(0, 0) + dNxi(1, 0) * x(1, 0) + dNxi(2, 0) * x(2, 0) +
415 dNxi(3, 0) * x(3, 0);
416 J(0, 1) = dNxi(0, 0) * x(0, 1) + dNxi(1, 0) * x(1, 1) + dNxi(2, 0) * x(2, 1) +
417 dNxi(3, 0) * x(3, 1);
418 J(1, 0) = dNxi(0, 1) * x(0, 0) + dNxi(1, 1) * x(1, 0) + dNxi(2, 1) * x(2, 0) +
419 dNxi(3, 1) * x(3, 0);
420 J(1, 1) = dNxi(0, 1) * x(0, 1) + dNxi(1, 1) * x(1, 1) + dNxi(2, 1) * x(2, 1) +
421 dNxi(3, 1) * x(3, 1);
422
423 double Jdet = detail::tensor<2>::inv(J, Jinv);
424
425 // radius for computation of volume
426 double rq = N(0) * x(0, 1) + N(1) * x(1, 1) + N(2) * x(2, 1) + N(3) * x(3, 1);
427
428 // dNx(m,i) += Jinv(i,j) * dNxi(m,j)
429 for (size_t m = 0; m < s_nne; ++m) {
430 // B(m, r, r, r) = dNdx(m,1)
431 B(m, 0, 0, 0) = Jinv(1, 0) * dNxi(m, 0) + Jinv(1, 1) * dNxi(m, 1);
432 // B(m, r, z, z) = dNdx(m,1)
433 B(m, 0, 2, 2) = Jinv(1, 0) * dNxi(m, 0) + Jinv(1, 1) * dNxi(m, 1);
434 // B(m, t, t, r)
435 B(m, 1, 1, 0) = 1.0 / rq * N(m);
436 // B(m, z, r, r) = dNdx(m,0)
437 B(m, 2, 0, 0) = Jinv(0, 0) * dNxi(m, 0) + Jinv(0, 1) * dNxi(m, 1);
438 // B(m, z, z, z) = dNdx(m,0)
439 B(m, 2, 2, 2) = Jinv(0, 0) * dNxi(m, 0) + Jinv(0, 1) * dNxi(m, 1);
440 }
441
442 m_vol(e, q) = m_w(q) * Jdet * 2.0 * M_PI * rq;
443 }
444 }
445 }
446 }
447
448 constexpr static size_t s_nne = 4;
449 constexpr static size_t s_ndim = 2;
450 constexpr static size_t s_tdim = 3;
451 size_t m_tdim = 3;
452 size_t m_nelem;
453 size_t m_nip;
454 array_type::tensor<double, 3> m_x;
455 array_type::tensor<double, 1> m_w;
456 array_type::tensor<double, 2> m_xi;
457 array_type::tensor<double, 2> m_N;
458 array_type::tensor<double, 3> m_dNxi;
459 array_type::tensor<double, 2> m_vol;
460 array_type::tensor<double, 6> m_B;
461};
462
463} // namespace Quad4
464} // namespace Element
465} // namespace GooseFEM
466
467#endif
QuadratureAxisymmetric(const T &x)
Constructor: use default Gauss integration.
const array_type::tensor< double, 6 > & B() const
Get the B-matrix (shape function gradients) (in global coordinates).
QuadratureAxisymmetric(const T &x, const X &xi, const W &w)
Constructor with custom 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.
xt::xtensor< T, N > tensor
Fixed (static) rank array.
Definition config.h:177
Toolbox to perform finite element computations.
Definition Allocate.h:14