MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
quadinterpolator.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "quadinterpolator.hpp"
13#include "qinterp/dispatch.hpp"
14#include "qspace.hpp"
15#include "../general/forall.hpp"
16#include "../linalg/dtensor.hpp"
17#include "../linalg/kernels.hpp"
18
19namespace mfem
20{
21
23 const IntegrationRule &ir):
24
25 fespace(&fes),
26 qspace(nullptr),
27 IntRule(&ir),
28 q_layout(QVectorLayout::byNODES),
29 use_tensor_products(UsesTensorBasis(fes))
30{
31 d_buffer.UseDevice(true);
32 if (fespace->GetNE() == 0) { return; }
33 const FiniteElement *fe = fespace->GetFE(0);
34 MFEM_VERIFY(dynamic_cast<const ScalarFiniteElement*>(fe) != NULL,
35 "Only scalar finite elements are supported");
36}
37
39 const QuadratureSpace &qs):
40
41 fespace(&fes),
42 qspace(&qs),
43 IntRule(nullptr),
44 q_layout(QVectorLayout::byNODES),
45 use_tensor_products(UsesTensorBasis(fes))
46{
47 d_buffer.UseDevice(true);
48 if (fespace->GetNE() == 0) { return; }
49 const FiniteElement *fe = fespace->GetFE(0);
50 MFEM_VERIFY(dynamic_cast<const ScalarFiniteElement*>(fe) != NULL,
51 "Only scalar finite elements are supported");
52}
53
54namespace internal
55{
56
57namespace quadrature_interpolator
58{
59
60// Compute kernel for 1D quadrature interpolation:
61// * non-tensor product version,
62// * assumes 'e_vec' is using ElementDofOrdering::NATIVE,
63// * assumes 'maps.mode == FULL'.
64static void Eval1D(const int NE,
65 const int vdim,
66 const QVectorLayout q_layout,
67 const GeometricFactors *geom,
68 const DofToQuad &maps,
69 const Vector &e_vec,
70 Vector &q_val,
71 Vector &q_der,
72 Vector &q_det,
73 const int eval_flags)
74{
75 using QI = QuadratureInterpolator;
76
77 const int nd = maps.ndof;
78 const int nq = maps.nqpt;
79 MFEM_ASSERT(maps.mode == DofToQuad::FULL, "internal error");
80 MFEM_ASSERT(!geom || geom->mesh->SpaceDimension() == 1, "");
81 MFEM_VERIFY(vdim == 1 || !(eval_flags & QI::DETERMINANTS), "");
82 MFEM_VERIFY(bool(geom) == bool(eval_flags & QI::PHYSICAL_DERIVATIVES),
83 "'geom' must be given (non-null) only when evaluating physical"
84 " derivatives");
85 const auto B = Reshape(maps.B.Read(), nq, nd);
86 const auto G = Reshape(maps.G.Read(), nq, nd);
87 const auto J = Reshape(geom ? geom->J.Read() : nullptr, nq, NE);
88 const auto E = Reshape(e_vec.Read(), nd, vdim, NE);
89 auto val = q_layout == QVectorLayout::byNODES ?
90 Reshape(q_val.Write(), nq, vdim, NE):
91 Reshape(q_val.Write(), vdim, nq, NE);
92 auto der = q_layout == QVectorLayout::byNODES ?
93 Reshape(q_der.Write(), nq, vdim, NE):
94 Reshape(q_der.Write(), vdim, nq, NE);
95 auto det = Reshape(q_det.Write(), nq, NE);
96 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
97 {
98 for (int q = 0; q < nq; ++q)
99 {
100 if (eval_flags & QI::VALUES)
101 {
102 for (int c = 0; c < vdim; c++)
103 {
104 real_t q_val = 0.0;
105 for (int d = 0; d < nd; ++d)
106 {
107 q_val += B(q,d)*E(d,c,e);
108 }
109 if (q_layout == QVectorLayout::byVDIM) { val(c,q,e) = q_val; }
110 if (q_layout == QVectorLayout::byNODES) { val(q,c,e) = q_val; }
111 }
112 }
113 if ((eval_flags & QI::DERIVATIVES) ||
114 (eval_flags & QI::PHYSICAL_DERIVATIVES) ||
115 (eval_flags & QI::DETERMINANTS))
116 {
117 for (int c = 0; c < vdim; c++)
118 {
119 real_t q_d = 0.0;
120 for (int d = 0; d < nd; ++d)
121 {
122 q_d += G(q,d)*E(d,c,e);
123 }
124 if (eval_flags & QI::PHYSICAL_DERIVATIVES)
125 {
126 q_d /= J(q,e);
127 }
128 if (eval_flags & QI::DERIVATIVES || eval_flags & QI::PHYSICAL_DERIVATIVES)
129 {
130 if (q_layout == QVectorLayout::byVDIM) { der(c,q,e) = q_d; }
131 if (q_layout == QVectorLayout::byNODES) { der(q,c,e) = q_d; }
132 }
133 if (vdim == 1 && (eval_flags & QI::DETERMINANTS))
134 {
135 det(q,e) = q_d;
136 }
137 }
138 }
139 }
140 });
141}
142
143// Template compute kernel for 2D quadrature interpolation:
144// * non-tensor product version,
145// * assumes 'e_vec' is using ElementDofOrdering::NATIVE,
146// * assumes 'maps.mode == FULL'.
147template<const int T_VDIM, const int T_ND, const int T_NQ>
148static void Eval2D(const int NE,
149 const int vdim,
150 const QVectorLayout q_layout,
151 const GeometricFactors *geom,
152 const DofToQuad &maps,
153 const Vector &e_vec,
154 Vector &q_val,
155 Vector &q_der,
156 Vector &q_det,
157 const int eval_flags)
158{
159 using QI = QuadratureInterpolator;
160
161 const int nd = maps.ndof;
162 const int nq = maps.nqpt;
163 const int ND = T_ND ? T_ND : nd;
164 const int NQ = T_NQ ? T_NQ : nq;
165 const int NMAX = NQ > ND ? NQ : ND;
166 const int VDIM = T_VDIM ? T_VDIM : vdim;
167 MFEM_ASSERT(maps.mode == DofToQuad::FULL, "internal error");
168 MFEM_ASSERT(!geom || geom->mesh->SpaceDimension() == 2, "");
169 MFEM_VERIFY(ND <= QI::MAX_ND2D, "");
170 MFEM_VERIFY(NQ <= QI::MAX_NQ2D, "");
171 MFEM_VERIFY(bool(geom) == bool(eval_flags & QI::PHYSICAL_DERIVATIVES),
172 "'geom' must be given (non-null) only when evaluating physical"
173 " derivatives");
174 const auto B = Reshape(maps.B.Read(), NQ, ND);
175 const auto G = Reshape(maps.G.Read(), NQ, 2, ND);
176 const auto J = Reshape(geom ? geom->J.Read() : nullptr, NQ, 2, 2, NE);
177 const auto E = Reshape(e_vec.Read(), ND, VDIM, NE);
178 auto val = q_layout == QVectorLayout::byNODES ?
179 Reshape(q_val.Write(), NQ, VDIM, NE):
180 Reshape(q_val.Write(), VDIM, NQ, NE);
181 auto der = q_layout == QVectorLayout::byNODES ?
182 Reshape(q_der.Write(), NQ, VDIM, 2, NE):
183 Reshape(q_der.Write(), VDIM, 2, NQ, NE);
184 auto det = Reshape(q_det.Write(), NQ, NE);
185 mfem::forall_2D(NE, NMAX, 1, [=] MFEM_HOST_DEVICE (int e)
186 {
187 const int ND = T_ND ? T_ND : nd;
188 const int NQ = T_NQ ? T_NQ : nq;
189 const int VDIM = T_VDIM ? T_VDIM : vdim;
190 constexpr int max_ND = T_ND ? T_ND : QI::MAX_ND2D;
191 constexpr int max_VDIM = T_VDIM ? T_VDIM : QI::MAX_VDIM2D;
192 MFEM_SHARED real_t s_E[max_VDIM*max_ND];
193 MFEM_FOREACH_THREAD(d, x, ND)
194 {
195 for (int c = 0; c < VDIM; c++)
196 {
197 s_E[c+d*VDIM] = E(d,c,e);
198 }
199 }
200 MFEM_SYNC_THREAD;
201
202 MFEM_FOREACH_THREAD(q, x, NQ)
203 {
204 if (eval_flags & QI::VALUES)
205 {
206 real_t ed[max_VDIM];
207 for (int c = 0; c < VDIM; c++) { ed[c] = 0.0; }
208 for (int d = 0; d < ND; ++d)
209 {
210 const real_t b = B(q,d);
211 for (int c = 0; c < VDIM; c++) { ed[c] += b*s_E[c+d*VDIM]; }
212 }
213 for (int c = 0; c < VDIM; c++)
214 {
215 if (q_layout == QVectorLayout::byVDIM) { val(c,q,e) = ed[c]; }
216 if (q_layout == QVectorLayout::byNODES) { val(q,c,e) = ed[c]; }
217 }
218 }
219 if ((eval_flags & QI::DERIVATIVES) ||
220 (eval_flags & QI::PHYSICAL_DERIVATIVES) ||
221 (eval_flags & QI::DETERMINANTS))
222 {
223 // use MAX_VDIM2D to avoid "subscript out of range" warnings
224 real_t D[QI::MAX_VDIM2D*2];
225 for (int i = 0; i < 2*VDIM; i++) { D[i] = 0.0; }
226 for (int d = 0; d < ND; ++d)
227 {
228 const real_t wx = G(q,0,d);
229 const real_t wy = G(q,1,d);
230 for (int c = 0; c < VDIM; c++)
231 {
232 real_t s_e = s_E[c+d*VDIM];
233 D[c+VDIM*0] += s_e * wx;
234 D[c+VDIM*1] += s_e * wy;
235 }
236 }
237 if (eval_flags & QI::DERIVATIVES)
238 {
239 for (int c = 0; c < VDIM; c++)
240 {
241 if (q_layout == QVectorLayout::byVDIM)
242 {
243 der(c,0,q,e) = D[c+VDIM*0];
244 der(c,1,q,e) = D[c+VDIM*1];
245 }
246 if (q_layout == QVectorLayout::byNODES)
247 {
248 der(q,c,0,e) = D[c+VDIM*0];
249 der(q,c,1,e) = D[c+VDIM*1];
250 }
251 }
252 }
253 if (eval_flags & QI::PHYSICAL_DERIVATIVES)
254 {
255 real_t Jloc[4], Jinv[4];
256 Jloc[0] = J(q,0,0,e);
257 Jloc[1] = J(q,1,0,e);
258 Jloc[2] = J(q,0,1,e);
259 Jloc[3] = J(q,1,1,e);
260 kernels::CalcInverse<2>(Jloc, Jinv);
261 for (int c = 0; c < VDIM; c++)
262 {
263 const real_t u = D[c+VDIM*0];
264 const real_t v = D[c+VDIM*1];
265 const real_t JiU = Jinv[0]*u + Jinv[1]*v;
266 const real_t JiV = Jinv[2]*u + Jinv[3]*v;
267 if (q_layout == QVectorLayout::byVDIM)
268 {
269 der(c,0,q,e) = JiU;
270 der(c,1,q,e) = JiV;
271 }
272 if (q_layout == QVectorLayout::byNODES)
273 {
274 der(q,c,0,e) = JiU;
275 der(q,c,1,e) = JiV;
276 }
277 }
278 }
279 if (eval_flags & QI::DETERMINANTS)
280 {
281 if (VDIM == 2) { det(q,e) = kernels::Det<2>(D); }
282 else
283 {
284 DeviceTensor<2> j(D, 3, 2);
285 const double E = j(0,0)*j(0,0) + j(1,0)*j(1,0) + j(2,0)*j(2,0);
286 const double F = j(0,0)*j(0,1) + j(1,0)*j(1,1) + j(2,0)*j(2,1);
287 const double G = j(0,1)*j(0,1) + j(1,1)*j(1,1) + j(2,1)*j(2,1);
288 det(q,e) = sqrt(E*G - F*F);
289 }
290 }
291 }
292 }
293 });
294}
295
296// Template compute kernel for 3D quadrature interpolation:
297// * non-tensor product version,
298// * assumes 'e_vec' is using ElementDofOrdering::NATIVE,
299// * assumes 'maps.mode == FULL'.
300template<const int T_VDIM, const int T_ND, const int T_NQ>
301static void Eval3D(const int NE,
302 const int vdim,
303 const QVectorLayout q_layout,
304 const GeometricFactors *geom,
305 const DofToQuad &maps,
306 const Vector &e_vec,
307 Vector &q_val,
308 Vector &q_der,
309 Vector &q_det,
310 const int eval_flags)
311{
312 using QI = QuadratureInterpolator;
313
314 const int nd = maps.ndof;
315 const int nq = maps.nqpt;
316 const int ND = T_ND ? T_ND : nd;
317 const int NQ = T_NQ ? T_NQ : nq;
318 const int NMAX = NQ > ND ? NQ : ND;
319 const int VDIM = T_VDIM ? T_VDIM : vdim;
320 MFEM_ASSERT(maps.mode == DofToQuad::FULL, "internal error");
321 MFEM_ASSERT(!geom || geom->mesh->SpaceDimension() == 3, "");
322 MFEM_VERIFY(ND <= QI::MAX_ND3D, "");
323 MFEM_VERIFY(NQ <= QI::MAX_NQ3D, "");
324 MFEM_VERIFY(VDIM == 3 || !(eval_flags & QI::DETERMINANTS), "");
325 MFEM_VERIFY(bool(geom) == bool(eval_flags & QI::PHYSICAL_DERIVATIVES),
326 "'geom' must be given (non-null) only when evaluating physical"
327 " derivatives");
328 const auto B = Reshape(maps.B.Read(), NQ, ND);
329 const auto G = Reshape(maps.G.Read(), NQ, 3, ND);
330 const auto J = Reshape(geom ? geom->J.Read() : nullptr, NQ, 3, 3, NE);
331 const auto E = Reshape(e_vec.Read(), ND, VDIM, NE);
332 auto val = q_layout == QVectorLayout::byNODES ?
333 Reshape(q_val.Write(), NQ, VDIM, NE):
334 Reshape(q_val.Write(), VDIM, NQ, NE);
335 auto der = q_layout == QVectorLayout::byNODES ?
336 Reshape(q_der.Write(), NQ, VDIM, 3, NE):
337 Reshape(q_der.Write(), VDIM, 3, NQ, NE);
338 auto det = Reshape(q_det.Write(), NQ, NE);
339 mfem::forall_2D(NE, NMAX, 1, [=] MFEM_HOST_DEVICE (int e)
340 {
341 const int ND = T_ND ? T_ND : nd;
342 const int NQ = T_NQ ? T_NQ : nq;
343 const int VDIM = T_VDIM ? T_VDIM : vdim;
344 constexpr int max_ND = T_ND ? T_ND : QI::MAX_ND3D;
345 constexpr int max_VDIM = T_VDIM ? T_VDIM : QI::MAX_VDIM3D;
346 MFEM_SHARED real_t s_E[max_VDIM*max_ND];
347 MFEM_FOREACH_THREAD(d, x, ND)
348 {
349 for (int c = 0; c < VDIM; c++)
350 {
351 s_E[c+d*VDIM] = E(d,c,e);
352 }
353 }
354 MFEM_SYNC_THREAD;
355
356 MFEM_FOREACH_THREAD(q, x, NQ)
357 {
358 if (eval_flags & QI::VALUES)
359 {
360 real_t ed[max_VDIM];
361 for (int c = 0; c < VDIM; c++) { ed[c] = 0.0; }
362 for (int d = 0; d < ND; ++d)
363 {
364 const real_t b = B(q,d);
365 for (int c = 0; c < VDIM; c++) { ed[c] += b*s_E[c+d*VDIM]; }
366 }
367 for (int c = 0; c < VDIM; c++)
368 {
369 if (q_layout == QVectorLayout::byVDIM) { val(c,q,e) = ed[c]; }
370 if (q_layout == QVectorLayout::byNODES) { val(q,c,e) = ed[c]; }
371 }
372 }
373 if ((eval_flags & QI::DERIVATIVES) ||
374 (eval_flags & QI::PHYSICAL_DERIVATIVES) ||
375 (eval_flags & QI::DETERMINANTS))
376 {
377 // use MAX_VDIM3D to avoid "subscript out of range" warnings
378 real_t D[QI::MAX_VDIM3D*3];
379 for (int i = 0; i < 3*VDIM; i++) { D[i] = 0.0; }
380 for (int d = 0; d < ND; ++d)
381 {
382 const real_t wx = G(q,0,d);
383 const real_t wy = G(q,1,d);
384 const real_t wz = G(q,2,d);
385 for (int c = 0; c < VDIM; c++)
386 {
387 real_t s_e = s_E[c+d*VDIM];
388 D[c+VDIM*0] += s_e * wx;
389 D[c+VDIM*1] += s_e * wy;
390 D[c+VDIM*2] += s_e * wz;
391 }
392 }
393 if (eval_flags & QI::DERIVATIVES)
394 {
395 for (int c = 0; c < VDIM; c++)
396 {
397 if (q_layout == QVectorLayout::byVDIM)
398 {
399 der(c,0,q,e) = D[c+VDIM*0];
400 der(c,1,q,e) = D[c+VDIM*1];
401 der(c,2,q,e) = D[c+VDIM*2];
402 }
403 if (q_layout == QVectorLayout::byNODES)
404 {
405 der(q,c,0,e) = D[c+VDIM*0];
406 der(q,c,1,e) = D[c+VDIM*1];
407 der(q,c,2,e) = D[c+VDIM*2];
408 }
409 }
410 }
411 if (eval_flags & QI::PHYSICAL_DERIVATIVES)
412 {
413 real_t Jloc[9], Jinv[9];
414 for (int col = 0; col < 3; col++)
415 {
416 for (int row = 0; row < 3; row++)
417 {
418 Jloc[row+3*col] = J(q,row,col,e);
419 }
420 }
421 kernels::CalcInverse<3>(Jloc, Jinv);
422 for (int c = 0; c < VDIM; c++)
423 {
424 const real_t u = D[c+VDIM*0];
425 const real_t v = D[c+VDIM*1];
426 const real_t w = D[c+VDIM*2];
427 const real_t JiU = Jinv[0]*u + Jinv[1]*v + Jinv[2]*w;
428 const real_t JiV = Jinv[3]*u + Jinv[4]*v + Jinv[5]*w;
429 const real_t JiW = Jinv[6]*u + Jinv[7]*v + Jinv[8]*w;
430 if (q_layout == QVectorLayout::byVDIM)
431 {
432 der(c,0,q,e) = JiU;
433 der(c,1,q,e) = JiV;
434 der(c,2,q,e) = JiW;
435 }
436 if (q_layout == QVectorLayout::byNODES)
437 {
438 der(q,c,0,e) = JiU;
439 der(q,c,1,e) = JiV;
440 der(q,c,2,e) = JiW;
441 }
442 }
443 }
444 if (VDIM == 3 && (eval_flags & QI::DETERMINANTS))
445 {
446 // The check (VDIM == 3) should eliminate this block when VDIM is
447 // known at compile time and (VDIM != 3).
448 det(q,e) = kernels::Det<3>(D);
449 }
450 }
451 }
452 });
453}
454
455} // namespace quadrature_interpolator
456
457} // namespace internal
458
460 unsigned eval_flags,
461 Vector &q_val,
462 Vector &q_der,
463 Vector &q_det) const
464{
465 using namespace internal::quadrature_interpolator;
466
467 const int ne = fespace->GetNE();
468 if (ne == 0) { return; }
469 const int vdim = fespace->GetVDim();
470 const FiniteElement *fe = fespace->GetFE(0);
471 const bool use_tensor_eval =
473 dynamic_cast<const TensorBasisElement*>(fe) != nullptr;
474 const IntegrationRule *ir =
476 const DofToQuad::Mode mode =
477 use_tensor_eval ? DofToQuad::TENSOR : DofToQuad::FULL;
478 const DofToQuad &maps = fe->GetDofToQuad(*ir, mode);
479 const int dim = maps.FE->GetDim();
480 const GeometricFactors *geom = nullptr;
481 if (eval_flags & PHYSICAL_DERIVATIVES)
482 {
483 const int jacobians = GeometricFactors::JACOBIANS;
484 geom = fespace->GetMesh()->GetGeometricFactors(*ir, jacobians);
485 }
486
487 MFEM_ASSERT(!(eval_flags & DETERMINANTS) || dim == vdim ||
488 (dim == 2 && vdim == 3), "Invalid dimensions for determinants.");
489 MFEM_ASSERT(fespace->GetMesh()->GetNumGeometries(
490 fespace->GetMesh()->Dimension()) == 1,
491 "mixed meshes are not supported");
492
493 if (use_tensor_eval)
494 {
495 // TODO: use fused kernels
497 {
498 if (eval_flags & VALUES)
499 {
500 TensorValues<QVectorLayout::byNODES>(ne, vdim, maps, e_vec, q_val);
501 }
502 if (eval_flags & DERIVATIVES)
503 {
504 TensorDerivatives<QVectorLayout::byNODES>(
505 ne, vdim, maps, e_vec, q_der);
506 }
507 if (eval_flags & PHYSICAL_DERIVATIVES)
508 {
509 TensorPhysDerivatives<QVectorLayout::byNODES>(
510 ne, vdim, maps, *geom, e_vec, q_der);
511 }
512 }
513
515 {
516 if (eval_flags & VALUES)
517 {
518 TensorValues<QVectorLayout::byVDIM>(ne, vdim, maps, e_vec, q_val);
519 }
520 if (eval_flags & DERIVATIVES)
521 {
522 TensorDerivatives<QVectorLayout::byVDIM>(
523 ne, vdim, maps, e_vec, q_der);
524 }
525 if (eval_flags & PHYSICAL_DERIVATIVES)
526 {
527 TensorPhysDerivatives<QVectorLayout::byVDIM>(
528 ne, vdim, maps, *geom, e_vec, q_der);
529 }
530 }
531 if (eval_flags & DETERMINANTS)
532 {
533 TensorDeterminants(ne, vdim, maps, e_vec, q_det, d_buffer);
534 }
535 }
536 else // use_tensor_eval == false
537 {
538 const int nd = maps.ndof;
539 const int nq = maps.nqpt;
540
541 void (*mult)(const int NE,
542 const int vdim,
544 const GeometricFactors *geom,
545 const DofToQuad &maps,
546 const Vector &e_vec,
547 Vector &q_val,
548 Vector &q_der,
549 Vector &q_det,
550 const int eval_flags) = NULL;
551
552 if (dim == 1)
553 {
554 mult = &Eval1D;
555 }
556 else if (vdim == 1) // dim == 2 || dim == 3
557 {
558 if (dim == 2)
559 {
560 switch (100*nd + nq)
561 {
562 // Q0
563 case 101: mult = &Eval2D<1,1,1>; break;
564 case 104: mult = &Eval2D<1,1,4>; break;
565 // Q1
566 case 404: mult = &Eval2D<1,4,4>; break;
567 case 409: mult = &Eval2D<1,4,9>; break;
568 // Q2
569 case 909: mult = &Eval2D<1,9,9>; break;
570 case 916: mult = &Eval2D<1,9,16>; break;
571 // Q3
572 case 1616: mult = &Eval2D<1,16,16>; break;
573 case 1625: mult = &Eval2D<1,16,25>; break;
574 case 1636: mult = &Eval2D<1,16,36>; break;
575 // Q4
576 case 2525: mult = &Eval2D<1,25,25>; break;
577 case 2536: mult = &Eval2D<1,25,36>; break;
578 case 2549: mult = &Eval2D<1,25,49>; break;
579 case 2564: mult = &Eval2D<1,25,64>; break;
580 }
581 if (nq >= 100 || !mult)
582 {
583 mult = &Eval2D<1,0,0>;
584 }
585 }
586 else if (dim == 3)
587 {
588 switch (1000*nd + nq)
589 {
590 // Q0
591 case 1001: mult = &Eval3D<1,1,1>; break;
592 case 1008: mult = &Eval3D<1,1,8>; break;
593 // Q1
594 case 8008: mult = &Eval3D<1,8,8>; break;
595 case 8027: mult = &Eval3D<1,8,27>; break;
596 // Q2
597 case 27027: mult = &Eval3D<1,27,27>; break;
598 case 27064: mult = &Eval3D<1,27,64>; break;
599 // Q3
600 case 64064: mult = &Eval3D<1,64,64>; break;
601 case 64125: mult = &Eval3D<1,64,125>; break;
602 case 64216: mult = &Eval3D<1,64,216>; break;
603 // Q4
604 case 125125: mult = &Eval3D<1,125,125>; break;
605 case 125216: mult = &Eval3D<1,125,216>; break;
606 }
607 if (nq >= 1000 || !mult)
608 {
609 mult = &Eval3D<1,0,0>;
610 }
611 }
612 }
613 else if (vdim == 3 && dim == 2)
614 {
615 switch (100*nd + nq)
616 {
617 // Q0
618 case 101: mult = &Eval2D<3,1,1>; break;
619 case 104: mult = &Eval2D<3,1,4>; break;
620 // Q1
621 case 404: mult = &Eval2D<3,4,4>; break;
622 case 409: mult = &Eval2D<3,4,9>; break;
623 // Q2
624 case 904: mult = &Eval2D<3,9,4>; break;
625 case 909: mult = &Eval2D<3,9,9>; break;
626 case 916: mult = &Eval2D<3,9,16>; break;
627 case 925: mult = &Eval2D<3,9,25>; break;
628 // Q3
629 case 1616: mult = &Eval2D<3,16,16>; break;
630 case 1625: mult = &Eval2D<3,16,25>; break;
631 case 1636: mult = &Eval2D<3,16,36>; break;
632 // Q4
633 case 2525: mult = &Eval2D<3,25,25>; break;
634 case 2536: mult = &Eval2D<3,25,36>; break;
635 case 2549: mult = &Eval2D<3,25,49>; break;
636 case 2564: mult = &Eval2D<3,25,64>; break;
637 default: mult = &Eval2D<3,0,0>;
638 }
639 }
640 else if (vdim == dim)
641 {
642 if (dim == 2)
643 {
644 switch (100*nd + nq)
645 {
646 // Q1
647 case 404: mult = &Eval2D<2,4,4>; break;
648 case 409: mult = &Eval2D<2,4,9>; break;
649 // Q2
650 case 909: mult = &Eval2D<2,9,9>; break;
651 case 916: mult = &Eval2D<2,9,16>; break;
652 // Q3
653 case 1616: mult = &Eval2D<2,16,16>; break;
654 case 1625: mult = &Eval2D<2,16,25>; break;
655 case 1636: mult = &Eval2D<2,16,36>; break;
656 // Q4
657 case 2525: mult = &Eval2D<2,25,25>; break;
658 case 2536: mult = &Eval2D<2,25,36>; break;
659 case 2549: mult = &Eval2D<2,25,49>; break;
660 case 2564: mult = &Eval2D<2,25,64>; break;
661 }
662 if (nq >= 100 || !mult) { mult = &Eval2D<2,0,0>; }
663 }
664 else if (dim == 3)
665 {
666 switch (1000*nd + nq)
667 {
668 // Q1
669 case 8008: mult = &Eval3D<3,8,8>; break;
670 case 8027: mult = &Eval3D<3,8,27>; break;
671 // Q2
672 case 27027: mult = &Eval3D<3,27,27>; break;
673 case 27064: mult = &Eval3D<3,27,64>; break;
674 case 27125: mult = &Eval3D<3,27,125>; break;
675 // Q3
676 case 64064: mult = &Eval3D<3,64,64>; break;
677 case 64125: mult = &Eval3D<3,64,125>; break;
678 case 64216: mult = &Eval3D<3,64,216>; break;
679 // Q4
680 case 125125: mult = &Eval3D<3,125,125>; break;
681 case 125216: mult = &Eval3D<3,125,216>; break;
682 }
683 if (nq >= 1000 || !mult) { mult = &Eval3D<3,0,0>; }
684 }
685 }
686 if (mult)
687 {
688 mult(ne,vdim,q_layout,geom,maps,e_vec,q_val,q_der,q_det,eval_flags);
689 }
690 else { MFEM_ABORT("case not supported yet"); }
691 }
692}
693
695 const Vector &q_val,
696 const Vector &q_der,
697 Vector &e_vec) const
698{
699 MFEM_CONTRACT_VAR(eval_flags);
700 MFEM_CONTRACT_VAR(q_val);
701 MFEM_CONTRACT_VAR(q_der);
702 MFEM_CONTRACT_VAR(e_vec);
703 MFEM_ABORT("this method is not implemented yet");
704}
705
707 Vector &q_val) const
708{
709 Vector empty;
710 Mult(e_vec, VALUES, q_val, empty, empty);
711}
712
714 Vector &q_der) const
715{
716 Vector empty;
717 Mult(e_vec, DERIVATIVES, empty, q_der, empty);
718}
719
721 Vector &q_der) const
722{
723 Vector empty;
724 Mult(e_vec, PHYSICAL_DERIVATIVES, empty, q_der, empty);
725}
726
728 Vector &q_det) const
729{
730 Vector empty;
731 Mult(e_vec, DETERMINANTS, empty, empty, q_det);
732}
733
734} // namespace mfem
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Definition fe_base.hpp:137
Mode mode
Describes the contents of the B, Bt, G, and Gt arrays, see Mode.
Definition fe_base.hpp:170
Array< real_t > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition fe_base.hpp:210
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Definition fe_base.hpp:150
@ FULL
Full multidimensional representation which does not use tensor product structure. The ordering of the...
Definition fe_base.hpp:154
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:161
Array< real_t > B
Basis functions evaluated at quadrature points.
Definition fe_base.hpp:189
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:174
const class FiniteElement * FE
The FiniteElement that created and owns this object.
Definition fe_base.hpp:141
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:178
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
Abstract class for all finite elements.
Definition fe_base.hpp:239
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition fe_base.cpp:365
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:316
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition mesh.hpp:2790
const Mesh * mesh
Definition mesh.hpp:2796
Vector J
Jacobians of the element transformations at all quadrature points.
Definition mesh.hpp:2829
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
Definition mesh.cpp:856
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition mesh.cpp:6961
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
bool use_tensor_products
Tensor product evaluation mode.
@ VALUES
Evaluate the values at quadrature points.
@ DERIVATIVES
Evaluate the derivatives at quadrature points.
@ PHYSICAL_DERIVATIVES
Evaluate the physical derivatives.
@ DETERMINANTS
Assuming the derivative at quadrature points form a matrix, this flag can be used to compute and stor...
QuadratureInterpolator(const FiniteElementSpace &fes, const IntegrationRule &ir)
void Mult(const Vector &e_vec, unsigned eval_flags, Vector &q_val, Vector &q_der, Vector &q_det) const
Interpolate the E-vector e_vec to quadrature points.
void Determinants(const Vector &e_vec, Vector &q_det) const
Compute the determinants of the derivatives (with respect to reference coordinates) of the E-vector e...
void MultTranspose(unsigned eval_flags, const Vector &q_val, const Vector &q_der, Vector &e_vec) const
Perform the transpose operation of Mult(). (TODO)
void Values(const Vector &e_vec, Vector &q_val) const
Interpolate the values of the E-vector e_vec at quadrature points.
void Derivatives(const Vector &e_vec, Vector &q_der) const
Interpolate the derivatives (with respect to reference coordinates) of the E-vector e_vec at quadratu...
QVectorLayout q_layout
Output Q-vector layout.
void PhysDerivatives(const Vector &e_vec, Vector &q_der) const
Interpolate the derivatives in physical space of the E-vector e_vec at quadrature points.
const IntegrationRule * IntRule
Not owned.
Vector d_buffer
Auxiliary device buffer.
const FiniteElementSpace * fespace
Not owned.
const QuadratureSpace * qspace
Not owned.
Class representing the storage layout of a QuadratureFunction.
Definition qspace.hpp:120
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition qspace.hpp:149
Class for finite elements with basis functions that return scalar values.
Definition fe_base.hpp:656
Vector data type.
Definition vector.hpp:80
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:474
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:136
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:482
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
MFEM_HOST_DEVICE void CalcInverse(const T *data, T *inv_data)
Return the inverse of a matrix with given size and data into the matrix with data inv_data.
Definition kernels.hpp:246
MFEM_HOST_DEVICE T Det(const T *data)
Compute the determinant of a square matrix of size dim with given data.
Definition kernels.hpp:237
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true,...
Definition device.hpp:337
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
void forall_2D(int N, int X, int Y, lambda &&body)
Definition forall.hpp:763
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition fespace.hpp:1345
QVectorLayout
Type describing possible layouts for Q-vectors.
Definition fespace.hpp:53
@ byNODES
NQPT x VDIM x NE (values) / NQPT x VDIM x DIM x NE (grads)
@ byVDIM
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
float real_t
Definition config.hpp:43
void forall(int N, lambda &&body)
Definition forall.hpp:754