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