MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
normal_deriv_restriction.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
13#include "fespace.hpp"
14#include "pgridfunc.hpp"
15#include "fe/face_map_utils.hpp"
16#include "../general/forall.hpp"
17
18namespace mfem
19{
20
21/// Compute the face index to volume index map "face_to_vol" in 2D
22static void NormalDerivativeSetupFaceIndexMap2D(
23 int nf, int d, const Array<int>& face_to_elem, Array<int>& face_to_vol)
24{
25 const auto f2e = Reshape(face_to_elem.HostRead(), 2, 2, nf);
26 auto f2v = Reshape(face_to_vol.HostWrite(), d, 2, nf);
27
28 for (int f = 0; f < nf; ++f)
29 {
30 const int fid0 = f2e(0, 1, f);
31 const int fid1 = f2e(1, 1, f);
32 for (int side = 0; side < 2; ++side)
33 {
34 const int el = f2e(side, 0, f);
35
36 if (el < 0)
37 {
38 for (int p = 0; p < d; ++p)
39 {
40 f2v(p, side, f) = -1;
41 }
42 }
43 else
44 {
45 for (int p = 0; p < d; ++p)
46 {
47 int i, j;
48 internal::FaceIdxToVolIdx2D(p, d, fid0, fid1, side, i, j);
49
50 f2v(p, side, f) = i + d * j;
51 }
52 }
53 }
54 }
55}
56
57/// Compute the face index to volume index map "face_to_vol" in 3D
58static void NormalDerivativeSetupFaceIndexMap3D(
59 int nf, int d, const Array<int>& face_to_elem, Array<int>& face_to_vol)
60{
61 const auto f2e = Reshape(face_to_elem.HostRead(), 2, 3, nf);
62 auto f2v = Reshape(face_to_vol.HostWrite(), d*d, 2, nf);
63
64 for (int f = 0; f < nf; ++f)
65 {
66 const int fid0 = f2e(0, 1, f);
67 const int fid1 = f2e(1, 1, f);
68 for (int side = 0; side < 2; ++side)
69 {
70 const int el = f2e(side, 0, f);
71 const int orientation = f2e(side, 2, f);
72
73 if (el < 0)
74 {
75 for (int p = 0; p < d*d; ++p)
76 {
77 f2v(p, side, f) = -1;
78 }
79 }
80 else
81 {
82 for (int p = 0; p < d*d; ++p)
83 {
84 int i, j, k; // 3D lexicographic index of quad point p
85 internal::FaceIdxToVolIdx3D(p, d, fid0, fid1, side, orientation, i, j, k);
86
87 f2v(p, side, f) = i + d * (j + d * k);
88 }
89 }
90 }
91 }
92}
93
95 const FiniteElementSpace &fes_,
96 const ElementDofOrdering f_ordering,
97 const FaceType face_type_)
98 : fes(fes_),
99 face_type(face_type_),
100 dim(fes.GetMesh()->Dimension()),
101 nf(fes.GetNFbyType(face_type)),
102 ne(fes.GetNE())
103{
104 MFEM_VERIFY(f_ordering == ElementDofOrdering::LEXICOGRAPHIC,
105 "Non-lexicographic ordering not currently supported in "
106 "L2NormalDerivativeFaceRestriction.");
107
108 Mesh &mesh = *fes.GetMesh();
109
110 const FiniteElement &fe = *fes.GetFE(0);
111 const int d = fe.GetDofToQuad(fe.GetNodes(), DofToQuad::TENSOR).ndof;
112
113 if (dim == 2)
114 {
115 // (el0, el1, fid0, fid1)
117 face_to_vol.SetSize(2 * nf * d);
118 }
119 else if (dim == 3)
120 {
121 // (el0, el1, fid0, fid1, or0, or1)
123 face_to_vol.SetSize(2 * nf * d * d);
124 }
125 else
126 {
127 MFEM_ABORT("Unsupported dimension.");
128 }
129 auto f2e = Reshape(face_to_elem.HostWrite(), 2, (dim == 2) ? 2 : 3, nf);
130
131 // Populate the face_to_elem array. The elem_indicator will be used to count
132 // the number of elements that are adjacent to faces of the given type.
133 Array<int> elem_indicator(ne);
134 elem_indicator = 0;
135
136 int f_ind = 0;
137 for (int f = 0; f < fes.GetNF(); ++f)
138 {
140
141 if (face.IsOfFaceType(face_type))
142 {
143 f2e(0, 0, f_ind) = face.element[0].index;
144 f2e(0, 1, f_ind) = face.element[0].local_face_id;
145 if (dim == 3)
146 {
147 f2e(0, 2, f_ind) = face.element[0].orientation;
148 }
149
150 elem_indicator[face.element[0].index] = 1;
151
153 {
154 const int el_idx_1 = face.element[1].index;
155 if (face.IsShared())
156 {
157 // Indicate shared face by index >= ne
158 f2e(1, 0, f_ind) = ne + el_idx_1;
159 }
160 else
161 {
162 // Face is not shared
163 f2e(1, 0, f_ind) = el_idx_1;
164 elem_indicator[el_idx_1] = 1;
165 }
166 f2e(1, 1, f_ind) = face.element[1].local_face_id;
167
168 if (dim == 3)
169 {
170 f2e(1, 2, f_ind) = face.element[1].orientation;
171 }
172 }
173 else
174 {
175 f2e(1, 0, f_ind) = -1;
176 f2e(1, 1, f_ind) = -1;
177
178 if (dim == 3)
179 {
180 f2e(1, 2, f_ind) = -1;
181 }
182 }
183
184 f_ind++;
185 }
186 }
187
188 // evaluate face to vol map
189 if (dim == 2)
190 {
191 NormalDerivativeSetupFaceIndexMap2D(nf, d, face_to_elem, face_to_vol);
192 }
193 else if (dim == 3)
194 {
195 NormalDerivativeSetupFaceIndexMap3D(nf, d, face_to_elem, face_to_vol);
196 }
197
198 // Number of elements adjacent to faces of face_type
199 ne_type = elem_indicator.Sum();
200
201 // In 2D: (el, f0,f1,f2,f3, s0,s1,s2,s3)
202 // In 3D: (el, f0,f1,f2,f3,f4,f5, s0,s1,s2,s3,s4,s5)
203 const int elem_data_sz = (dim == 2) ? 9 : 13;
204
205 elem_to_face.SetSize(elem_data_sz * ne_type);
206 elem_to_face = -1;
207
208 auto e2f = Reshape(elem_to_face.HostWrite(), elem_data_sz, ne_type);
209 elem_indicator.PartialSum();
210
211 const int nsides = (face_type == FaceType::Interior) ? 2 : 1;
212 const int side_begin = (dim == 2) ? 5 : 7;
213 for (int f = 0; f < nf; ++f)
214 {
215 for (int side = 0; side < nsides; ++side)
216 {
217 const int el = f2e(side, 0, f);
218 // Skip shared faces
219 if (el < ne)
220 {
221 const int face_id = f2e(side, 1, f);
222
223 const int e = elem_indicator[el] - 1;
224 e2f(0, e) = el;
225 e2f(1 + face_id, e) = f;
226 e2f(side_begin + face_id, e) = side;
227 }
228 }
229 }
230}
231
233{
234 if (nf == 0) { return; }
235 switch (dim)
236 {
237 case 2:
238 {
239 const int d1d = fes.GetElementOrder(0) + 1;
240 switch (d1d)
241 {
242 case 1: Mult2D<1>(x, y); break;
243 case 2: Mult2D<2>(x, y); break;
244 case 3: Mult2D<3>(x, y); break;
245 case 4: Mult2D<4>(x, y); break;
246 case 5: Mult2D<5>(x, y); break;
247 case 6: Mult2D<6>(x, y); break;
248 case 7: Mult2D<7>(x, y); break;
249 case 8: Mult2D<8>(x, y); break;
250 default: Mult2D(x, y); break;
251 }
252 }
253 break;
254 case 3:
255 {
256 const int d1d = fes.GetElementOrder(0) + 1;
257 switch (d1d)
258 {
259 case 1: Mult3D<1>(x, y); break;
260 case 2: Mult3D<2>(x, y); break;
261 case 3: Mult3D<3>(x, y); break;
262 case 4: Mult3D<4>(x, y); break;
263 case 5: Mult3D<5>(x, y); break;
264 case 6: Mult3D<6>(x, y); break;
265 case 7: Mult3D<7>(x, y); break;
266 case 8: Mult3D<8>(x, y); break;
267 default: Mult3D(x, y); break; // fallback
268 }
269 break;
270 }
271 default: MFEM_ABORT("Dimension not supported."); break;
272 }
273}
274
276 const Vector &x, Vector &y, const real_t a) const
277{
278 if (nf == 0) { return; }
279 switch (dim)
280 {
281 case 2:
282 {
283 const int d1d = fes.GetElementOrder(0) + 1;
284 switch (d1d)
285 {
286 case 1: AddMultTranspose2D<1>(x, y, a); break;
287 case 2: AddMultTranspose2D<2>(x, y, a); break;
288 case 3: AddMultTranspose2D<3>(x, y, a); break;
289 case 4: AddMultTranspose2D<4>(x, y, a); break;
290 case 5: AddMultTranspose2D<5>(x, y, a); break;
291 case 6: AddMultTranspose2D<6>(x, y, a); break;
292 case 7: AddMultTranspose2D<7>(x, y, a); break;
293 case 8: AddMultTranspose2D<8>(x, y, a); break;
294 default: AddMultTranspose2D(x, y, a); break;
295 }
296 }
297 break;
298 case 3:
299 {
300 const int d1d = fes.GetElementOrder(0) + 1;
301 switch (d1d)
302 {
303 case 1: AddMultTranspose3D<1>(x, y, a); break;
304 case 2: AddMultTranspose3D<2>(x, y, a); break;
305 case 3: AddMultTranspose3D<3>(x, y, a); break;
306 case 4: AddMultTranspose3D<4>(x, y, a); break;
307 case 5: AddMultTranspose3D<5>(x, y, a); break;
308 case 6: AddMultTranspose3D<6>(x, y, a); break;
309 case 7: AddMultTranspose3D<7>(x, y, a); break;
310 case 8: AddMultTranspose3D<8>(x, y, a); break;
311 default: AddMultTranspose3D(x, y, a); break; // fallback
312 }
313 break;
314 }
315 default: MFEM_ABORT("Not yet implemented"); break;
316 }
317}
318
319template <int T_D1D>
321{
322 const int vd = fes.GetVDim();
323 const bool t = fes.GetOrdering() == Ordering::byVDIM;
324 const int num_elem = ne;
325
326 const FiniteElement &fe = *fes.GetFE(0);
327 const DofToQuad &maps = fe.GetDofToQuad(fe.GetNodes(), DofToQuad::TENSOR);
328
329 const int q = maps.nqpt;
330 const int d = maps.ndof;
331
332 Vector face_nbr_data = GetLVectorFaceNbrData(fes, x, face_type);
333 const int ne_shared = face_nbr_data.Size() / d / d / vd;
334
335 MFEM_VERIFY(q == d, "");
336 MFEM_VERIFY(T_D1D == d || T_D1D == 0, "");
337
338 // derivative of 1D basis function
339 const auto G_ = Reshape(maps.G.Read(), q, d);
340 // (el0, el1, fid0, fid1)
341 const auto f2e = Reshape(face_to_elem.Read(), 2, 2, nf);
342
343 const auto f2v = Reshape(face_to_vol.Read(), q, 2, nf);
344
345 // if byvdim, d_x has shape (vdim, nddof, nddof, ne)
346 // otherwise, d_x has shape (nddof, nddof, ne, vdim)
347 const auto d_x = Reshape(x.Read(), t?vd:d, d, t?d:ne, t?ne:vd);
348 const auto d_x_shared = Reshape(face_nbr_data.Read(),
349 t?vd:d, d, t?d:ne_shared, t?ne_shared:vd);
350 auto d_y = Reshape(y.Write(), q, vd, 2, nf);
351
352 mfem::forall_2D(nf, 2, q, [=] MFEM_HOST_DEVICE (int f) -> void
353 {
354 constexpr int MD = (T_D1D) ? T_D1D : DofQuadLimits::MAX_D1D;
355
356 MFEM_SHARED real_t G_s[MD*MD];
357 DeviceMatrix G(G_s, q, d);
358
359 MFEM_SHARED int E[2];
360 MFEM_SHARED int FID[2];
361 MFEM_SHARED int F2V[2][MD];
362
363 if (MFEM_THREAD_ID(x) == 0)
364 {
365 MFEM_FOREACH_THREAD(j, y, d)
366 {
367 for (int i = 0; i < q; ++i)
368 {
369 G(i, j) = G_(i, j);
370 }
371 }
372 }
373
374 MFEM_FOREACH_THREAD(side, x, 2)
375 {
376 if (MFEM_THREAD_ID(y) == 0)
377 {
378 E[side] = f2e(side, 0, f);
379 FID[side] = f2e(side, 1, f);
380 }
381
382 MFEM_FOREACH_THREAD(j, y, d)
383 {
384 F2V[side][j] = f2v(j, side, f);
385 }
386 }
387 MFEM_SYNC_THREAD;
388
389 MFEM_FOREACH_THREAD(side, x, 2)
390 {
391 const int el = E[side];
392 const bool shared = (el >= num_elem);
393 const auto &d_x_e = shared ? d_x_shared : d_x;
394 const int el_idx = shared ? el - num_elem : el;
395
396 const int face_id = FID[side];
397
398 MFEM_FOREACH_THREAD(p, y, q)
399 {
400 if (el < 0)
401 {
402 for (int c = 0; c < vd; ++c)
403 {
404 d_y(p, c, side, f) = 0.0;
405 }
406 }
407 else
408 {
409 const int ij = F2V[side][p];
410 const int i = ij % q;
411 const int j = ij / q;
412
413 for (int c=0; c < vd; ++c)
414 {
415 real_t grad_n = 0;
416 for (int kk=0; kk < d; ++kk)
417 {
418 const int k = (face_id == 0 || face_id == 2) ? i : kk;
419 const int l = (face_id == 0 || face_id == 2) ? kk : j;
420 const real_t g = (face_id == 0 || face_id == 2) ? G(j,l) : G(i,k);
421 grad_n += g * d_x_e(t?c:k, t?k:l, t?l:el_idx, t?el_idx:c);
422 }
423 d_y(p, c, side, f) = grad_n;
424 }
425 }
426 }
427 }
428 });
429}
430
431template <int T_D1D>
433{
434 const int vd = fes.GetVDim();
435 const bool t = fes.GetOrdering() == Ordering::byVDIM;
436 const int num_elem = ne;
437
438 const FiniteElement &fe = *fes.GetFE(0);
439 const DofToQuad &maps = fe.GetDofToQuad(fe.GetNodes(), DofToQuad::TENSOR);
440
441 const int q = maps.nqpt;
442 const int d = maps.ndof;
443 const int q2d = q * q;
444
445 Vector face_nbr_data = GetLVectorFaceNbrData(fes, x, face_type);
446 const int ne_shared = face_nbr_data.Size() / d / d / d / vd;
447
448 MFEM_VERIFY(q == d, "");
449 MFEM_VERIFY(T_D1D == d || T_D1D == 0, "");
450
451 const auto G_ = Reshape(maps.G.Read(), q, d);
452 // (el0, el1, fid0, fid1, or0, or1)
453 const auto f2e = Reshape(face_to_elem.Read(), 2, 3, nf);
454 const auto f2v = Reshape(face_to_vol.Read(), q2d, 2, nf);
455
456 // t ? (vdim, d, d, d, ne) : (d, d, d, ne, vdim)
457 const auto d_x = Reshape(x.Read(), t?vd:d, d, d, t?d:ne, t?ne:vd);
458 const auto d_x_shared = Reshape(face_nbr_data.Read(),
459 t?vd:d, d, d, t?d:ne_shared, t?ne_shared:vd);
460 auto d_y = Reshape(y.Write(), q2d, vd, 2, nf);
461
462 mfem::forall_2D(nf, q2d, 2, [=] MFEM_HOST_DEVICE (int f) -> void
463 {
464 static constexpr int MD = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
465
466 MFEM_SHARED real_t G_s[MD*MD];
467 DeviceMatrix G(G_s, d, q);
468
469 MFEM_SHARED int E[2];
470 MFEM_SHARED int FID[2];
471 MFEM_SHARED int F2V[2][MD*MD];
472
473 // Load G matrix into shared memory
474 if (MFEM_THREAD_ID(y) == 0)
475 {
476 MFEM_FOREACH_THREAD(j, x, d*q)
477 {
478 const int p = j % q;
479 const int k = j / q;
480 G(k, p) = G_(p, k);
481 }
482 }
483
484 MFEM_FOREACH_THREAD(side, y, 2)
485 {
486 if (MFEM_THREAD_ID(x) == 0)
487 {
488 E[side] = f2e(side, 0, f);
489 FID[side] = f2e(side, 1, f);
490 }
491 MFEM_FOREACH_THREAD(j, x, q2d)
492 {
493 F2V[side][j] = f2v(j, side, f);
494 }
495 }
496 MFEM_SYNC_THREAD;
497
498 MFEM_FOREACH_THREAD(side, y, 2)
499 {
500 const int el = E[side];
501 const bool shared = (el >= num_elem);
502 const auto &d_x_e = shared ? d_x_shared : d_x;
503 const int el_idx = shared ? el - num_elem : el;
504
505 const int face_id = FID[side];
506
507 // Is this face parallel to the x-y plane in reference coordinates?
508 const bool xy_plane = (face_id == 0 || face_id == 5);
509 const bool xz_plane = (face_id == 1 || face_id == 3);
510 const bool yz_plane = (face_id == 2 || face_id == 4);
511
512 MFEM_FOREACH_THREAD(p, x, q2d)
513 {
514 if (el_idx < 0)
515 {
516 for (int c = 0; c < vd; ++c)
517 {
518 d_y(p, c, side, f) = 0.0;
519 }
520 }
521 else
522 {
523 const int ijk = F2V[side][p];
524 const int k = ijk / q2d;
525 const int i = ijk % q;
526 const int j = (ijk - q2d*k) / q;
527
528 // the fixed 1D index of the normal component of the face
529 // quadrature point
530 const int g_row = yz_plane ? i : xz_plane ? j : k;
531
532 for (int c = 0; c < vd; ++c)
533 {
534 real_t grad_n = 0.0;
535
536 for (int kk = 0; kk < d; ++kk)
537 {
538 // (l, m, n) 3D lexicographic index of interior points used
539 // in evaluating normal derivatives
540 const int l = yz_plane ? kk : i;
541 const int m = xz_plane ? kk : j;
542 const int n = xy_plane ? kk : k;
543
544 const real_t g = G(kk, g_row);
545
546 grad_n += g * d_x_e(t?c:l, t?l:m, t?m:n, t?n:el_idx, t?el_idx:c);
547 }
548 d_y(p, c, side, f) = grad_n;
549 }
550 }
551 }
552 }
553 });
554}
555
556template <int T_D1D>
558 const Vector &y, Vector &x, const real_t a) const
559{
560 const int vd = fes.GetVDim();
561 const bool t = fes.GetOrdering() == Ordering::byVDIM;
562
563 const FiniteElement &fe = *fes.GetFE(0);
564 const DofToQuad &maps = fe.GetDofToQuad(fe.GetNodes(), DofToQuad::TENSOR);
565
566 const int q = maps.nqpt;
567 const int d = maps.ndof;
568
569 // derivative of 1D basis function
570 auto G_ = Reshape(maps.G.Read(), q, d);
571
572 // entries of e2f: (el,f0,f1,f2,f3,s0,s1,s2,s3)
573 auto e2f = Reshape(elem_to_face.Read(), 9, ne_type);
574
575 auto f2v = Reshape(face_to_vol.Read(), d, 2, nf);
576
577 // if byvdim, d_x has shape (vdim, nddof, nddof, ne)
578 // otherwise, d_x has shape (nddof, nddof, ne, vdim)
579 auto d_x = Reshape(x.ReadWrite(), t?vd:d, d, t?d:ne, t?ne:vd);
580 auto d_y = Reshape(y.Read(), q, vd, 2, nf);
581
582 mfem::forall_2D(ne_type, d, d, [=] MFEM_HOST_DEVICE (int e)
583 {
584 constexpr int MD = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
585
586 MFEM_SHARED real_t y_s[MD];
587 MFEM_SHARED int pp[MD];
588 MFEM_SHARED int jj;
589 if (MFEM_THREAD_ID(x) == 0 && MFEM_THREAD_ID(y) == 0) { jj = 0; }
590
591 MFEM_SHARED real_t BG[MD*MD];
592 DeviceMatrix G(BG, q, d);
593
594 MFEM_SHARED real_t x_s[MD*MD];
595 DeviceMatrix xx(x_s, d, d);
596
597 MFEM_SHARED int el; // global element index
598 MFEM_SHARED int faces[4];
599 MFEM_SHARED int sides[4];
600
601 MFEM_FOREACH_THREAD(i,x,d)
602 {
603 MFEM_FOREACH_THREAD(p,y,q)
604 {
605 G(p,i) = a * G_(p,i);
606 xx(p,i) = 0.0;
607 }
608 }
609
610 if (MFEM_THREAD_ID(y) == 0)
611 {
612 if (MFEM_THREAD_ID(x) == 0)
613 {
614 el = e2f(0, e);
615 }
616
617 MFEM_FOREACH_THREAD(i, x, 4)
618 {
619 faces[i] = e2f(1 + i, e);
620 sides[i] = e2f(5 + i, e);
621 }
622 }
623 MFEM_SYNC_THREAD;
624
625 for (int face_id=0; face_id < 4; ++face_id)
626 {
627 const int f = faces[face_id];
628
629 if (f < 0) { continue; }
630
631 const int side = sides[face_id];
632
633 if (MFEM_THREAD_ID(y) == 0)
634 {
635 MFEM_FOREACH_THREAD(p,x,d)
636 {
637 y_s[p] = d_y(p, 0, side, f);
638
639 const int ij = f2v(p, side, f);
640 const int i = ij % q;
641 const int j = ij / q;
642
643 pp[(face_id == 0 || face_id == 2) ? i : j] = p;
644 if (MFEM_THREAD_ID(x) == 0)
645 {
646 jj = (face_id == 0 || face_id == 2) ? j : i;
647 }
648 }
649 }
650 MFEM_SYNC_THREAD;
651
652 MFEM_FOREACH_THREAD(k,x,d)
653 {
654 MFEM_FOREACH_THREAD(l,y,d)
655 {
656 const int p = (face_id == 0 || face_id == 2) ? pp[k] : pp[l];
657 const int kk = (face_id == 0 || face_id == 2) ? l : k;
658 const real_t g = G(jj, kk);
659 xx(k,l) += g * y_s[p];
660 }
661 }
662 }
663 MFEM_SYNC_THREAD;
664
665 MFEM_FOREACH_THREAD(k,x,d)
666 {
667 MFEM_FOREACH_THREAD(l,y,d)
668 {
669 const int c = 0;
670 d_x(t?c:k, t?k:l, t?l:el, t?el:c) += xx(k,l);
671 }
672 }
673 });
674}
675
676template <int T_D1D>
678 const Vector &y, Vector &x, const real_t a) const
679{
680 const int vd = fes.GetVDim();
681 const bool t = fes.GetOrdering() == Ordering::byVDIM;
682
683 MFEM_VERIFY(vd == 1, "vdim > 1 not supported.");
684
685 const FiniteElement &fe = *fes.GetFE(0);
686 const DofToQuad &maps = fe.GetDofToQuad(fe.GetNodes(), DofToQuad::TENSOR);
687
688 const int q = maps.nqpt;
689 const int d = maps.ndof;
690 const int q2d = q * q;
691
692 MFEM_VERIFY(q == d, "");
693 MFEM_VERIFY(T_D1D == d || T_D1D == 0, "");
694
695 auto G_ = Reshape(maps.G.Read(), q, d);
696
697 // (el, f0,f1,f2,f3,f4,f5, s0,s1,s2,s3,s4,s5)
698 auto e2f = Reshape(elem_to_face.Read(), 13, ne_type);
699
700 auto f2v = Reshape(face_to_vol.Read(), q2d, 2, nf);
701
702 auto d_x = Reshape(x.ReadWrite(), t?vd:d, d, d, t?d:ne, t?ne:vd);
703 const auto d_y = Reshape(y.Read(), q2d, vd, 2, nf);
704
705 mfem::forall_2D(ne_type, q, q, [=] MFEM_HOST_DEVICE (int e) -> void
706 {
707 static constexpr int MD = T_D1D ? T_D1D : DofQuadLimits::MAX_D1D;
708
709 MFEM_SHARED int pp[MD][MD];
710 MFEM_SHARED real_t y_s[MD*MD];
711 MFEM_SHARED int jj;
712 if (MFEM_THREAD_ID(x) == 0 && MFEM_THREAD_ID(y) == 0) { jj = 0; }
713
714 MFEM_SHARED real_t xx_s[MD*MD*MD];
715 auto xx = Reshape(xx_s, d, d, d);
716
717 MFEM_SHARED real_t G_s[MD*MD];
718 DeviceMatrix G(G_s, q, d);
719
720 MFEM_SHARED int el;
721 MFEM_SHARED int faces[6];
722 MFEM_SHARED int sides[6];
723
724 // Load G into shared memory
725 MFEM_FOREACH_THREAD(j, x, d)
726 {
727 MFEM_FOREACH_THREAD(i, y, q)
728 {
729 G(i, j) = a * G_(i, j);
730 G(i, j) = a * G_(i, j);
731 G(i, j) = a * G_(i, j);
732 }
733 }
734
735 if (MFEM_THREAD_ID(y) == 0)
736 {
737 if (MFEM_THREAD_ID(x) == 0)
738 {
739 el = e2f(0, e); // global element index
740 }
741
742 MFEM_FOREACH_THREAD(i, x, 6)
743 {
744 faces[i] = e2f(1 + i, e);
745 sides[i] = e2f(7 + i, e);
746 }
747 }
748
749 MFEM_FOREACH_THREAD(k, x, d)
750 {
751 MFEM_FOREACH_THREAD(j, y, d)
752 {
753 for (int i = 0; i < d; ++i)
754 {
755 xx(i, j, k) = 0.0;
756 }
757 }
758 }
759 MFEM_SYNC_THREAD;
760
761 for (int face_id = 0; face_id < 6; ++face_id)
762 {
763 const int f = faces[face_id];
764
765 if (f < 0)
766 {
767 continue;
768 }
769
770 const int side = sides[face_id];
771
772 // is this face parallel to the x-y plane in reference coordinates?
773 const bool xy_plane = (face_id == 0 || face_id == 5);
774 const bool xz_plane = (face_id == 1 || face_id == 3);
775
776 MFEM_FOREACH_THREAD(p1, x, q)
777 {
778 MFEM_FOREACH_THREAD(p2, y, q)
779 {
780 const int p = p1 + q * p2;
781 y_s[p] = d_y(p, 0, side, f);
782
783 const int ijk = f2v(p, side, f);
784 const int k = ijk / q2d;
785 const int i = ijk % q;
786 const int j = (ijk - q2d*k) / q;
787
788 pp[(xy_plane || xz_plane) ? i : j][(xy_plane) ? j : k] = p;
789 if (MFEM_THREAD_ID(x) == 0 && MFEM_THREAD_ID(y) == 0)
790 {
791 jj = (xy_plane) ? k : (xz_plane) ? j : i;
792 }
793 }
794 }
795 MFEM_SYNC_THREAD;
796
797 MFEM_FOREACH_THREAD(n, x, d)
798 {
799 MFEM_FOREACH_THREAD(m, y, d)
800 {
801 for (int l = 0; l < d; ++l)
802 {
803 const int p = (xy_plane) ? pp[l][m] : (xz_plane) ? pp[l][n] : pp[m][n];
804 const int kk = (xy_plane) ? n : (xz_plane) ? m : l;
805 const real_t g = G(jj, kk);
806 xx(l, m, n) += g * y_s[p];
807 }
808 }
809 }
810 }
811 MFEM_SYNC_THREAD;
812
813 // map back to global array
814 MFEM_FOREACH_THREAD(n, x, d)
815 {
816 MFEM_FOREACH_THREAD(m, y, d)
817 {
818 for (int l = 0; l < d; ++l)
819 {
820 const int c = 0;
821 d_x(t?c:l, t?l:m, t?m:n, t?n:el, t?el:c) += xx(l, m, n);
822 }
823 }
824 }
825 });
826}
827
828} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
T * HostWrite()
Shortcut for mfem::Write(a.GetMemory(), a.Size(), false).
Definition array.hpp:329
A basic generic Tensor class, appropriate for use on the GPU.
Definition dtensor.hpp:82
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Definition fe_base.hpp:137
Array< real_t > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition fe_base.hpp:210
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:161
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
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
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition fespace.hpp:746
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
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:725
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
Definition fespace.cpp:176
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
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:395
void Mult(const Vector &x, Vector &y) const
Computes the normal derivatives on the face_type faces of the mesh.
Array< int > elem_to_face
Element-wise information array.
L2NormalDerivativeFaceRestriction(const FiniteElementSpace &fes_, const ElementDofOrdering f_ordering, const FaceType face_type_)
Constructor.
const int nf
Number of faces of the given face_type.
Array< int > face_to_elem
Face-wise information array.
const FiniteElementSpace & fes
The L2 finite element space.
int ne_type
Number of elements with faces of type face type.
const FaceType face_type
Face type: either boundary or interior.
void Mult2D(const Vector &x, Vector &y) const
void Mult3D(const Vector &x, Vector &y) const
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const
Computes the transpose of the action of Mult(), accumulating into y with coefficient a.
Array< int > face_to_vol
maps face index to volume index
void AddMultTranspose2D(const Vector &x, Vector &y, const real_t a) const
void AddMultTranspose3D(const Vector &x, Vector &y, const real_t a) const
Mesh data type.
Definition mesh.hpp:56
FaceInformation GetFaceInformation(int f) const
Definition mesh.cpp:1169
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 real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:490
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
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
Mesh * GetMesh(int type)
Definition ex29.cpp:218
real_t a
Definition lissajous.cpp:41
Vector GetLVectorFaceNbrData(const FiniteElementSpace &fes, const Vector &x, FaceType ftype)
Return the face-neighbor data given the L-vector x.
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
float real_t
Definition config.hpp:43
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:75
@ LEXICOGRAPHIC
Lexicographic ordering for tensor-product FiniteElements.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
FaceType
Definition mesh.hpp:47
real_t p(const Vector &x, real_t t)
RefCoord t[3]
This structure is used as a human readable output format that deciphers the information contained in ...
Definition mesh.hpp:1894
bool IsOfFaceType(FaceType type) const
Return true if the face is of the same type as type.
Definition mesh.hpp:1941
struct mfem::Mesh::FaceInformation::@13 element[2]
bool IsShared() const
Return true if the face is a shared interior face which is NOT a master nonconforming face.
Definition mesh.hpp:1919