MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_hcurlhdiv_kernels.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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#ifndef MFEM_BILININTEG_HCURLHDIV_KERNELS_HPP
13#define MFEM_BILININTEG_HCURLHDIV_KERNELS_HPP
14
20#include "../bilininteg.hpp"
21
22/// \cond DO_NOT_DOCUMENT
23namespace mfem
24{
25
26namespace internal
27{
28
29// PA H(curl)-H(div) Mass Apply 2D kernel
30void PAHcurlHdivMassSetup2D(const int Q1D,
31 const int coeffDim,
32 const int NE,
33 const bool transpose,
34 const Array<real_t> &w_,
35 const Vector &j,
36 Vector &coeff_,
37 Vector &op);
38
39// PA H(curl)-H(div) Mass Assemble 3D kernel
40void PAHcurlHdivMassSetup3D(const int Q1D,
41 const int coeffDim,
42 const int NE,
43 const bool transpose,
44 const Array<real_t> &w_,
45 const Vector &j,
46 Vector &coeff_,
47 Vector &op);
48
49// PA H(curl)-H(div) Mass Apply 2D kernel
50void PAHcurlHdivMassApply2D(const int D1D,
51 const int D1Dtest,
52 const int Q1D,
53 const int NE,
54 const bool scalarCoeff,
55 const bool trialHcurl,
56 const bool transpose,
57 const Array<real_t> &Bo_,
58 const Array<real_t> &Bc_,
59 const Array<real_t> &Bot_,
60 const Array<real_t> &Bct_,
61 const Vector &op_,
62 const Vector &x_,
63 Vector &y_);
64
65// PA H(curl)-H(div) Mass Apply 3D kernel
66void PAHcurlHdivMassApply3D(const int D1D,
67 const int D1Dtest,
68 const int Q1D,
69 const int NE,
70 const bool scalarCoeff,
71 const bool trialHcurl,
72 const bool transpose,
73 const Array<real_t> &Bo_,
74 const Array<real_t> &Bc_,
75 const Array<real_t> &Bot_,
76 const Array<real_t> &Bct_,
77 const Vector &op_,
78 const Vector &x_,
79 Vector &y_);
80
81// PA H(curl)-H(div) Curl Apply 3D kernel
82template<int T_D1D = 0, int T_D1D_TEST = 0, int T_Q1D = 0>
83inline void PAHcurlHdivApply3D(const int d1d,
84 const int d1dtest,
85 const int q1d,
86 const int NE,
87 const Array<real_t> &bo,
88 const Array<real_t> &bc,
89 const Array<real_t> &bot,
90 const Array<real_t> &bct,
91 const Array<real_t> &gc,
92 const Vector &pa_data,
93 const Vector &x,
94 Vector &y)
95{
96 MFEM_VERIFY(T_D1D || d1d <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D,
97 "Error: d1d > HCURL_MAX_D1D");
98 MFEM_VERIFY(T_D1D_TEST || d1dtest <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D,
99 "Error: d1dtest > HCURL_MAX_D1D");
100 MFEM_VERIFY(T_Q1D || q1d <= DeviceDofQuadLimits::Get().HCURL_MAX_Q1D,
101 "Error: q1d > HCURL_MAX_Q1D");
102 const int D1D = T_D1D ? T_D1D : d1d;
103 const int D1Dtest = T_D1D_TEST ? T_D1D_TEST : d1dtest;
104 const int Q1D = T_Q1D ? T_Q1D : q1d;
105
106 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
107 auto Bc = Reshape(bc.Read(), Q1D, D1D);
108 auto Bot = Reshape(bot.Read(), D1Dtest-1, Q1D);
109 auto Bct = Reshape(bct.Read(), D1Dtest, Q1D);
110 auto Gc = Reshape(gc.Read(), Q1D, D1D);
111 auto op = Reshape(pa_data.Read(), Q1D, Q1D, Q1D, 6, NE);
112 auto X = Reshape(x.Read(), 3*(D1D-1)*D1D*D1D, NE);
113 auto Y = Reshape(y.ReadWrite(), 3*(D1Dtest-1)*(D1Dtest-1)*D1Dtest, NE);
114
115 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
116 {
117 // Using Piola transformations (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u}
118 // for u in H(curl) and w = (1 / det (dF)) dF \hat{w} for w in H(div), we get
119 // (\nabla\times u) \cdot w = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{w}
120 // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
121 // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
122 // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
123
124 constexpr int VDIM = 3;
125 constexpr int MD1D = T_D1D ? T_D1D :
126 DofQuadLimits::HCURL_MAX_D1D; // Assuming HDIV_MAX_D1D <= HCURL_MAX_D1D
127 constexpr int MQ1D = T_Q1D ? T_Q1D : DofQuadLimits::HCURL_MAX_Q1D;
128 const int D1D = T_D1D ? T_D1D : d1d;
129 const int D1Dtest = T_D1D_TEST ? T_D1D_TEST : d1dtest;
130 const int Q1D = T_Q1D ? T_Q1D : q1d;
131
132 real_t curl[MQ1D][MQ1D][MQ1D][VDIM];
133 // curl[qz][qy][qx] will be computed as the vector curl at each quadrature point.
134
135 for (int qz = 0; qz < Q1D; ++qz)
136 {
137 for (int qy = 0; qy < Q1D; ++qy)
138 {
139 for (int qx = 0; qx < Q1D; ++qx)
140 {
141 for (int c = 0; c < VDIM; ++c)
142 {
143 curl[qz][qy][qx][c] = 0.0;
144 }
145 }
146 }
147 }
148
149 // We treat x, y, z components separately for optimization specific to each.
150
151 int osc = 0;
152
153 {
154 // x component
155 const int D1Dz = D1D;
156 const int D1Dy = D1D;
157 const int D1Dx = D1D - 1;
158
159 for (int dz = 0; dz < D1Dz; ++dz)
160 {
161 real_t gradXY[MQ1D][MQ1D][2];
162 for (int qy = 0; qy < Q1D; ++qy)
163 {
164 for (int qx = 0; qx < Q1D; ++qx)
165 {
166 for (int d = 0; d < 2; ++d)
167 {
168 gradXY[qy][qx][d] = 0.0;
169 }
170 }
171 }
172
173 for (int dy = 0; dy < D1Dy; ++dy)
174 {
175 real_t massX[MQ1D];
176 for (int qx = 0; qx < Q1D; ++qx)
177 {
178 massX[qx] = 0.0;
179 }
180
181 for (int dx = 0; dx < D1Dx; ++dx)
182 {
183 const real_t t = X(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
184 for (int qx = 0; qx < Q1D; ++qx)
185 {
186 massX[qx] += t * Bo(qx,dx);
187 }
188 }
189
190 for (int qy = 0; qy < Q1D; ++qy)
191 {
192 const real_t wy = Bc(qy,dy);
193 const real_t wDy = Gc(qy,dy);
194 for (int qx = 0; qx < Q1D; ++qx)
195 {
196 const real_t wx = massX[qx];
197 gradXY[qy][qx][0] += wx * wDy;
198 gradXY[qy][qx][1] += wx * wy;
199 }
200 }
201 }
202
203 for (int qz = 0; qz < Q1D; ++qz)
204 {
205 const real_t wz = Bc(qz,dz);
206 const real_t wDz = Gc(qz,dz);
207 for (int qy = 0; qy < Q1D; ++qy)
208 {
209 for (int qx = 0; qx < Q1D; ++qx)
210 {
211 // \hat{\nabla}\times\hat{u} is [0, (u_0)_{x_2}, -(u_0)_{x_1}]
212 curl[qz][qy][qx][1] += gradXY[qy][qx][1] * wDz; // (u_0)_{x_2}
213 curl[qz][qy][qx][2] -= gradXY[qy][qx][0] * wz; // -(u_0)_{x_1}
214 }
215 }
216 }
217 }
218
219 osc += D1Dx * D1Dy * D1Dz;
220 }
221
222 {
223 // y component
224 const int D1Dz = D1D;
225 const int D1Dy = D1D - 1;
226 const int D1Dx = D1D;
227
228 for (int dz = 0; dz < D1Dz; ++dz)
229 {
230 real_t gradXY[MQ1D][MQ1D][2];
231 for (int qy = 0; qy < Q1D; ++qy)
232 {
233 for (int qx = 0; qx < Q1D; ++qx)
234 {
235 for (int d = 0; d < 2; ++d)
236 {
237 gradXY[qy][qx][d] = 0.0;
238 }
239 }
240 }
241
242 for (int dx = 0; dx < D1Dx; ++dx)
243 {
244 real_t massY[MQ1D];
245 for (int qy = 0; qy < Q1D; ++qy)
246 {
247 massY[qy] = 0.0;
248 }
249
250 for (int dy = 0; dy < D1Dy; ++dy)
251 {
252 const real_t t = X(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
253 for (int qy = 0; qy < Q1D; ++qy)
254 {
255 massY[qy] += t * Bo(qy,dy);
256 }
257 }
258
259 for (int qx = 0; qx < Q1D; ++qx)
260 {
261 const real_t wx = Bc(qx,dx);
262 const real_t wDx = Gc(qx,dx);
263 for (int qy = 0; qy < Q1D; ++qy)
264 {
265 const real_t wy = massY[qy];
266 gradXY[qy][qx][0] += wDx * wy;
267 gradXY[qy][qx][1] += wx * wy;
268 }
269 }
270 }
271
272 for (int qz = 0; qz < Q1D; ++qz)
273 {
274 const real_t wz = Bc(qz,dz);
275 const real_t wDz = Gc(qz,dz);
276 for (int qy = 0; qy < Q1D; ++qy)
277 {
278 for (int qx = 0; qx < Q1D; ++qx)
279 {
280 // \hat{\nabla}\times\hat{u} is [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
281 curl[qz][qy][qx][0] -= gradXY[qy][qx][1] * wDz; // -(u_1)_{x_2}
282 curl[qz][qy][qx][2] += gradXY[qy][qx][0] * wz; // (u_1)_{x_0}
283 }
284 }
285 }
286 }
287
288 osc += D1Dx * D1Dy * D1Dz;
289 }
290
291 {
292 // z component
293 const int D1Dz = D1D - 1;
294 const int D1Dy = D1D;
295 const int D1Dx = D1D;
296
297 for (int dx = 0; dx < D1Dx; ++dx)
298 {
299 real_t gradYZ[MQ1D][MQ1D][2];
300 for (int qz = 0; qz < Q1D; ++qz)
301 {
302 for (int qy = 0; qy < Q1D; ++qy)
303 {
304 for (int d = 0; d < 2; ++d)
305 {
306 gradYZ[qz][qy][d] = 0.0;
307 }
308 }
309 }
310
311 for (int dy = 0; dy < D1Dy; ++dy)
312 {
313 real_t massZ[MQ1D];
314 for (int qz = 0; qz < Q1D; ++qz)
315 {
316 massZ[qz] = 0.0;
317 }
318
319 for (int dz = 0; dz < D1Dz; ++dz)
320 {
321 const real_t t = X(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
322 for (int qz = 0; qz < Q1D; ++qz)
323 {
324 massZ[qz] += t * Bo(qz,dz);
325 }
326 }
327
328 for (int qy = 0; qy < Q1D; ++qy)
329 {
330 const real_t wy = Bc(qy,dy);
331 const real_t wDy = Gc(qy,dy);
332 for (int qz = 0; qz < Q1D; ++qz)
333 {
334 const real_t wz = massZ[qz];
335 gradYZ[qz][qy][0] += wz * wy;
336 gradYZ[qz][qy][1] += wz * wDy;
337 }
338 }
339 }
340
341 for (int qx = 0; qx < Q1D; ++qx)
342 {
343 const real_t wx = Bc(qx,dx);
344 const real_t wDx = Gc(qx,dx);
345
346 for (int qy = 0; qy < Q1D; ++qy)
347 {
348 for (int qz = 0; qz < Q1D; ++qz)
349 {
350 // \hat{\nabla}\times\hat{u} is [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
351 curl[qz][qy][qx][0] += gradYZ[qz][qy][1] * wx; // (u_2)_{x_1}
352 curl[qz][qy][qx][1] -= gradYZ[qz][qy][0] * wDx; // -(u_2)_{x_0}
353 }
354 }
355 }
356 }
357 }
358
359 // Apply D operator.
360 for (int qz = 0; qz < Q1D; ++qz)
361 {
362 for (int qy = 0; qy < Q1D; ++qy)
363 {
364 for (int qx = 0; qx < Q1D; ++qx)
365 {
366 const real_t O11 = op(qx,qy,qz,0,e);
367 const real_t O12 = op(qx,qy,qz,1,e);
368 const real_t O13 = op(qx,qy,qz,2,e);
369 const real_t O22 = op(qx,qy,qz,3,e);
370 const real_t O23 = op(qx,qy,qz,4,e);
371 const real_t O33 = op(qx,qy,qz,5,e);
372
373 const real_t c1 = (O11 * curl[qz][qy][qx][0]) + (O12 * curl[qz][qy][qx][1]) +
374 (O13 * curl[qz][qy][qx][2]);
375 const real_t c2 = (O12 * curl[qz][qy][qx][0]) + (O22 * curl[qz][qy][qx][1]) +
376 (O23 * curl[qz][qy][qx][2]);
377 const real_t c3 = (O13 * curl[qz][qy][qx][0]) + (O23 * curl[qz][qy][qx][1]) +
378 (O33 * curl[qz][qy][qx][2]);
379
380 curl[qz][qy][qx][0] = c1;
381 curl[qz][qy][qx][1] = c2;
382 curl[qz][qy][qx][2] = c3;
383 }
384 }
385 }
386
387 for (int qz = 0; qz < Q1D; ++qz)
388 {
389 real_t massXY[MD1D][MD1D];
390
391 osc = 0;
392
393 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
394 {
395 const int D1Dz = (c == 2) ? D1Dtest : D1Dtest - 1;
396 const int D1Dy = (c == 1) ? D1Dtest : D1Dtest - 1;
397 const int D1Dx = (c == 0) ? D1Dtest : D1Dtest - 1;
398
399 for (int dy = 0; dy < D1Dy; ++dy)
400 {
401 for (int dx = 0; dx < D1Dx; ++dx)
402 {
403 massXY[dy][dx] = 0;
404 }
405 }
406 for (int qy = 0; qy < Q1D; ++qy)
407 {
408 real_t massX[MD1D];
409 for (int dx = 0; dx < D1Dx; ++dx)
410 {
411 massX[dx] = 0;
412 }
413 for (int qx = 0; qx < Q1D; ++qx)
414 {
415 for (int dx = 0; dx < D1Dx; ++dx)
416 {
417 massX[dx] += curl[qz][qy][qx][c] *
418 ((c == 0) ? Bct(dx,qx) : Bot(dx,qx));
419 }
420 }
421 for (int dy = 0; dy < D1Dy; ++dy)
422 {
423 const real_t wy = (c == 1) ? Bct(dy,qy) : Bot(dy,qy);
424 for (int dx = 0; dx < D1Dx; ++dx)
425 {
426 massXY[dy][dx] += massX[dx] * wy;
427 }
428 }
429 }
430
431 for (int dz = 0; dz < D1Dz; ++dz)
432 {
433 const real_t wz = (c == 2) ? Bct(dz,qz) : Bot(dz,qz);
434 for (int dy = 0; dy < D1Dy; ++dy)
435 {
436 for (int dx = 0; dx < D1Dx; ++dx)
437 {
438 Y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) +=
439 massXY[dy][dx] * wz;
440 }
441 }
442 }
443
444 osc += D1Dx * D1Dy * D1Dz;
445 } // loop c
446 } // loop qz
447 }); // end of element loop
448}
449
450// PA H(curl)-H(div) Curl Apply Transpose 3D kernel
451template<int T_D1D = 0, int T_D1D_TEST = 0, int T_Q1D = 0>
452inline void PAHcurlHdivApplyTranspose3D(const int d1d,
453 const int d1dtest,
454 const int q1d,
455 const int NE,
456 const Array<real_t> &bo,
457 const Array<real_t> &bc,
458 const Array<real_t> &bot,
459 const Array<real_t> &bct,
460 const Array<real_t> &gct,
461 const Vector &pa_data,
462 const Vector &x,
463 Vector &y)
464{
465 MFEM_VERIFY(T_D1D || d1d <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D,
466 "Error: d1d > HCURL_MAX_D1D");
467 MFEM_VERIFY(T_D1D_TEST || d1dtest <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D,
468 "Error: d1dtest > HCURL_MAX_D1D");
469 MFEM_VERIFY(T_Q1D || q1d <= DeviceDofQuadLimits::Get().HCURL_MAX_Q1D,
470 "Error: q1d > HCURL_MAX_Q1D");
471 const int D1D = T_D1D ? T_D1D : d1d;
472 const int D1Dtest = T_D1D_TEST ? T_D1D_TEST : d1dtest;
473 const int Q1D = T_Q1D ? T_Q1D : q1d;
474
475 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
476 auto Bc = Reshape(bc.Read(), Q1D, D1D);
477 auto Bot = Reshape(bot.Read(), D1Dtest-1, Q1D);
478 auto Bct = Reshape(bct.Read(), D1Dtest, Q1D);
479 auto Gct = Reshape(gct.Read(), D1D, Q1D);
480 auto op = Reshape(pa_data.Read(), Q1D, Q1D, Q1D, 6, NE);
481 auto X = Reshape(x.Read(), 3*(D1Dtest-1)*(D1Dtest-1)*D1Dtest, NE);
482 auto Y = Reshape(y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
483
484 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
485 {
486 // Using Piola transformations (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u}
487 // for u in H(curl) and w = (1 / det (dF)) dF \hat{w} for w in H(div), we get
488 // (\nabla\times u) \cdot w = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{w}
489 // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
490 // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
491 // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
492
493 constexpr int VDIM = 3;
494 constexpr int MD1D = T_D1D ? T_D1D :
495 DofQuadLimits::HCURL_MAX_D1D; // Assuming HDIV_MAX_D1D <= HCURL_MAX_D1D
496 constexpr int MQ1D = T_Q1D ? T_Q1D : DofQuadLimits::HCURL_MAX_Q1D;
497 const int D1D = T_D1D ? T_D1D : d1d;
498 const int D1Dtest = T_D1D_TEST ? T_D1D_TEST : d1dtest;
499 const int Q1D = T_Q1D ? T_Q1D : q1d;
500
501 real_t mass[MQ1D][MQ1D][MQ1D][VDIM];
502
503 for (int qz = 0; qz < Q1D; ++qz)
504 {
505 for (int qy = 0; qy < Q1D; ++qy)
506 {
507 for (int qx = 0; qx < Q1D; ++qx)
508 {
509 for (int c = 0; c < VDIM; ++c)
510 {
511 mass[qz][qy][qx][c] = 0.0;
512 }
513 }
514 }
515 }
516
517 int osc = 0;
518
519 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
520 {
521 const int D1Dz = (c == 2) ? D1Dtest : D1Dtest - 1;
522 const int D1Dy = (c == 1) ? D1Dtest : D1Dtest - 1;
523 const int D1Dx = (c == 0) ? D1Dtest : D1Dtest - 1;
524
525 for (int dz = 0; dz < D1Dz; ++dz)
526 {
527 real_t massXY[MQ1D][MQ1D];
528 for (int qy = 0; qy < Q1D; ++qy)
529 {
530 for (int qx = 0; qx < Q1D; ++qx)
531 {
532 massXY[qy][qx] = 0.0;
533 }
534 }
535
536 for (int dy = 0; dy < D1Dy; ++dy)
537 {
538 real_t massX[MQ1D];
539 for (int qx = 0; qx < Q1D; ++qx)
540 {
541 massX[qx] = 0.0;
542 }
543
544 for (int dx = 0; dx < D1Dx; ++dx)
545 {
546 const real_t t = X(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
547 for (int qx = 0; qx < Q1D; ++qx)
548 {
549 massX[qx] += t * ((c == 0) ? Bc(qx,dx) : Bo(qx,dx));
550 }
551 }
552
553 for (int qy = 0; qy < Q1D; ++qy)
554 {
555 const real_t wy = (c == 1) ? Bc(qy,dy) : Bo(qy,dy);
556 for (int qx = 0; qx < Q1D; ++qx)
557 {
558 const real_t wx = massX[qx];
559 massXY[qy][qx] += wx * wy;
560 }
561 }
562 }
563
564 for (int qz = 0; qz < Q1D; ++qz)
565 {
566 const real_t wz = (c == 2) ? Bc(qz,dz) : Bo(qz,dz);
567 for (int qy = 0; qy < Q1D; ++qy)
568 {
569 for (int qx = 0; qx < Q1D; ++qx)
570 {
571 mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
572 }
573 }
574 }
575 }
576
577 osc += D1Dx * D1Dy * D1Dz;
578 } // loop (c) over components
579
580 // Apply D operator.
581 for (int qz = 0; qz < Q1D; ++qz)
582 {
583 for (int qy = 0; qy < Q1D; ++qy)
584 {
585 for (int qx = 0; qx < Q1D; ++qx)
586 {
587 const real_t O11 = op(qx,qy,qz,0,e);
588 const real_t O12 = op(qx,qy,qz,1,e);
589 const real_t O13 = op(qx,qy,qz,2,e);
590 const real_t O22 = op(qx,qy,qz,3,e);
591 const real_t O23 = op(qx,qy,qz,4,e);
592 const real_t O33 = op(qx,qy,qz,5,e);
593 const real_t massX = mass[qz][qy][qx][0];
594 const real_t massY = mass[qz][qy][qx][1];
595 const real_t massZ = mass[qz][qy][qx][2];
596 mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
597 mass[qz][qy][qx][1] = (O12*massX)+(O22*massY)+(O23*massZ);
598 mass[qz][qy][qx][2] = (O13*massX)+(O23*massY)+(O33*massZ);
599 }
600 }
601 }
602
603 // x component
604 osc = 0;
605 {
606 const int D1Dz = D1D;
607 const int D1Dy = D1D;
608 const int D1Dx = D1D - 1;
609
610 for (int qz = 0; qz < Q1D; ++qz)
611 {
612 real_t gradXY12[MD1D][MD1D];
613 real_t gradXY21[MD1D][MD1D];
614
615 for (int dy = 0; dy < D1Dy; ++dy)
616 {
617 for (int dx = 0; dx < D1Dx; ++dx)
618 {
619 gradXY12[dy][dx] = 0.0;
620 gradXY21[dy][dx] = 0.0;
621 }
622 }
623 for (int qy = 0; qy < Q1D; ++qy)
624 {
625 real_t massX[MD1D][2];
626 for (int dx = 0; dx < D1Dx; ++dx)
627 {
628 for (int n = 0; n < 2; ++n)
629 {
630 massX[dx][n] = 0.0;
631 }
632 }
633 for (int qx = 0; qx < Q1D; ++qx)
634 {
635 for (int dx = 0; dx < D1Dx; ++dx)
636 {
637 const real_t wx = Bot(dx,qx);
638
639 massX[dx][0] += wx * mass[qz][qy][qx][1];
640 massX[dx][1] += wx * mass[qz][qy][qx][2];
641 }
642 }
643 for (int dy = 0; dy < D1Dy; ++dy)
644 {
645 const real_t wy = Bct(dy,qy);
646 const real_t wDy = Gct(dy,qy);
647
648 for (int dx = 0; dx < D1Dx; ++dx)
649 {
650 gradXY21[dy][dx] += massX[dx][0] * wy;
651 gradXY12[dy][dx] += massX[dx][1] * wDy;
652 }
653 }
654 }
655
656 for (int dz = 0; dz < D1Dz; ++dz)
657 {
658 const real_t wz = Bct(dz,qz);
659 const real_t wDz = Gct(dz,qz);
660 for (int dy = 0; dy < D1Dy; ++dy)
661 {
662 for (int dx = 0; dx < D1Dx; ++dx)
663 {
664 // \hat{\nabla}\times\hat{u} is [0, (u_0)_{x_2}, -(u_0)_{x_1}]
665 // (u_0)_{x_2} * (op * curl)_1 - (u_0)_{x_1} * (op * curl)_2
666 Y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
667 e) += (gradXY21[dy][dx] * wDz) - (gradXY12[dy][dx] * wz);
668 }
669 }
670 }
671 } // loop qz
672
673 osc += D1Dx * D1Dy * D1Dz;
674 }
675
676 // y component
677 {
678 const int D1Dz = D1D;
679 const int D1Dy = D1D - 1;
680 const int D1Dx = D1D;
681
682 for (int qz = 0; qz < Q1D; ++qz)
683 {
684 real_t gradXY02[MD1D][MD1D];
685 real_t gradXY20[MD1D][MD1D];
686
687 for (int dy = 0; dy < D1Dy; ++dy)
688 {
689 for (int dx = 0; dx < D1Dx; ++dx)
690 {
691 gradXY02[dy][dx] = 0.0;
692 gradXY20[dy][dx] = 0.0;
693 }
694 }
695 for (int qx = 0; qx < Q1D; ++qx)
696 {
697 real_t massY[MD1D][2];
698 for (int dy = 0; dy < D1Dy; ++dy)
699 {
700 massY[dy][0] = 0.0;
701 massY[dy][1] = 0.0;
702 }
703 for (int qy = 0; qy < Q1D; ++qy)
704 {
705 for (int dy = 0; dy < D1Dy; ++dy)
706 {
707 const real_t wy = Bot(dy,qy);
708
709 massY[dy][0] += wy * mass[qz][qy][qx][2];
710 massY[dy][1] += wy * mass[qz][qy][qx][0];
711 }
712 }
713 for (int dx = 0; dx < D1Dx; ++dx)
714 {
715 const real_t wx = Bct(dx,qx);
716 const real_t wDx = Gct(dx,qx);
717
718 for (int dy = 0; dy < D1Dy; ++dy)
719 {
720 gradXY02[dy][dx] += massY[dy][0] * wDx;
721 gradXY20[dy][dx] += massY[dy][1] * wx;
722 }
723 }
724 }
725
726 for (int dz = 0; dz < D1Dz; ++dz)
727 {
728 const real_t wz = Bct(dz,qz);
729 const real_t wDz = Gct(dz,qz);
730 for (int dy = 0; dy < D1Dy; ++dy)
731 {
732 for (int dx = 0; dx < D1Dx; ++dx)
733 {
734 // \hat{\nabla}\times\hat{u} is [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
735 // -(u_1)_{x_2} * (op * curl)_0 + (u_1)_{x_0} * (op * curl)_2
736 Y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
737 e) += (-gradXY20[dy][dx] * wDz) + (gradXY02[dy][dx] * wz);
738 }
739 }
740 }
741 } // loop qz
742
743 osc += D1Dx * D1Dy * D1Dz;
744 }
745
746 // z component
747 {
748 const int D1Dz = D1D - 1;
749 const int D1Dy = D1D;
750 const int D1Dx = D1D;
751
752 for (int qx = 0; qx < Q1D; ++qx)
753 {
754 real_t gradYZ01[MD1D][MD1D];
755 real_t gradYZ10[MD1D][MD1D];
756
757 for (int dy = 0; dy < D1Dy; ++dy)
758 {
759 for (int dz = 0; dz < D1Dz; ++dz)
760 {
761 gradYZ01[dz][dy] = 0.0;
762 gradYZ10[dz][dy] = 0.0;
763 }
764 }
765 for (int qy = 0; qy < Q1D; ++qy)
766 {
767 real_t massZ[MD1D][2];
768 for (int dz = 0; dz < D1Dz; ++dz)
769 {
770 for (int n = 0; n < 2; ++n)
771 {
772 massZ[dz][n] = 0.0;
773 }
774 }
775 for (int qz = 0; qz < Q1D; ++qz)
776 {
777 for (int dz = 0; dz < D1Dz; ++dz)
778 {
779 const real_t wz = Bot(dz,qz);
780
781 massZ[dz][0] += wz * mass[qz][qy][qx][0];
782 massZ[dz][1] += wz * mass[qz][qy][qx][1];
783 }
784 }
785 for (int dy = 0; dy < D1Dy; ++dy)
786 {
787 const real_t wy = Bct(dy,qy);
788 const real_t wDy = Gct(dy,qy);
789
790 for (int dz = 0; dz < D1Dz; ++dz)
791 {
792 gradYZ01[dz][dy] += wy * massZ[dz][1];
793 gradYZ10[dz][dy] += wDy * massZ[dz][0];
794 }
795 }
796 }
797
798 for (int dx = 0; dx < D1Dx; ++dx)
799 {
800 const real_t wx = Bct(dx,qx);
801 const real_t wDx = Gct(dx,qx);
802
803 for (int dy = 0; dy < D1Dy; ++dy)
804 {
805 for (int dz = 0; dz < D1Dz; ++dz)
806 {
807 // \hat{\nabla}\times\hat{u} is [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
808 // (u_2)_{x_1} * (op * curl)_0 - (u_2)_{x_0} * (op * curl)_1
809 Y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
810 e) += (gradYZ10[dz][dy] * wx) - (gradYZ01[dz][dy] * wDx);
811 }
812 }
813 }
814 } // loop qx
815 }
816 }); // end of element loop
817}
818
819} // namespace internal
820
821} // namespace mfem
822
823/// \endcond DO_NOT_DOCUMENT
824
825#endif
mfem::real_t real_t
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:138
void forall(int N, lambda &&body)
Definition forall.hpp:839
static const DeviceDofQuadLimits & Get()
Return a const reference to the DeviceDofQuadLimits singleton.
Definition forall.hpp:128