MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_hcurl_kernels.cpp
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
13
14namespace mfem
15{
16
17namespace internal
18{
19
20void PAHcurlMassAssembleDiagonal2D(const int D1D,
21 const int Q1D,
22 const int NE,
23 const bool symmetric,
24 const Array<real_t> &bo,
25 const Array<real_t> &bc,
26 const Vector &pa_data,
27 Vector &diag)
28{
29 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
30 auto Bc = Reshape(bc.Read(), Q1D, D1D);
31 auto op = Reshape(pa_data.Read(), Q1D, Q1D, symmetric ? 3 : 4, NE);
32 auto D = Reshape(diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
33
34 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
35 {
36 constexpr static int VDIM = 2;
37 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
38
39 int osc = 0;
40
41 for (int c = 0; c < VDIM; ++c) // loop over x, y components
42 {
43 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
44 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
45
46 real_t mass[MAX_Q1D];
47
48 for (int dy = 0; dy < D1Dy; ++dy)
49 {
50 for (int qx = 0; qx < Q1D; ++qx)
51 {
52 mass[qx] = 0.0;
53 for (int qy = 0; qy < Q1D; ++qy)
54 {
55 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
56
57 mass[qx] += wy * wy * ((c == 0) ? op(qx,qy,0,e) :
58 op(qx,qy,symmetric ? 2 : 3, e));
59 }
60 }
61
62 for (int dx = 0; dx < D1Dx; ++dx)
63 {
64 for (int qx = 0; qx < Q1D; ++qx)
65 {
66 const real_t wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
67 D(dx + (dy * D1Dx) + osc, e) += mass[qx] * wx * wx;
68 }
69 }
70 }
71
72 osc += D1Dx * D1Dy;
73 } // loop c
74 }); // end of element loop
75}
76
77void PAHcurlMassAssembleDiagonal3D(const int D1D,
78 const int Q1D,
79 const int NE,
80 const bool symmetric,
81 const Array<real_t> &bo,
82 const Array<real_t> &bc,
83 const Vector &pa_data,
84 Vector &diag)
85{
86 MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D,
87 "Error: D1D > MAX_D1D");
88 MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().HCURL_MAX_Q1D,
89 "Error: Q1D > MAX_Q1D");
90 constexpr static int VDIM = 3;
91
92 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
93 auto Bc = Reshape(bc.Read(), Q1D, D1D);
94 auto op = Reshape(pa_data.Read(), Q1D, Q1D, Q1D, symmetric ? 6 : 9, NE);
95 auto D = Reshape(diag.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
96
97 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
98 {
99 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
100
101 int osc = 0;
102
103 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
104 {
105 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
106 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
107 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
108
109 const int opc = (c == 0) ? 0 : ((c == 1) ? (symmetric ? 3 : 4) :
110 (symmetric ? 5 : 8));
111
112 real_t mass[MAX_Q1D];
113
114 for (int dz = 0; dz < D1Dz; ++dz)
115 {
116 for (int dy = 0; dy < D1Dy; ++dy)
117 {
118 for (int qx = 0; qx < Q1D; ++qx)
119 {
120 mass[qx] = 0.0;
121 for (int qy = 0; qy < Q1D; ++qy)
122 {
123 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
124
125 for (int qz = 0; qz < Q1D; ++qz)
126 {
127 const real_t wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
128
129 mass[qx] += wy * wy * wz * wz * op(qx,qy,qz,opc,e);
130 }
131 }
132 }
133
134 for (int dx = 0; dx < D1Dx; ++dx)
135 {
136 for (int qx = 0; qx < Q1D; ++qx)
137 {
138 const real_t wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
139 D(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += mass[qx] * wx * wx;
140 }
141 }
142 }
143 }
144
145 osc += D1Dx * D1Dy * D1Dz;
146 } // loop c
147 }); // end of element loop
148}
149
150void PAHcurlMassApply2D(const int D1D,
151 const int Q1D,
152 const int NE,
153 const bool symmetric,
154 const Array<real_t> &bo,
155 const Array<real_t> &bc,
156 const Array<real_t> &bot,
157 const Array<real_t> &bct,
158 const Vector &pa_data,
159 const Vector &x,
160 Vector &y)
161{
162 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
163 auto Bc = Reshape(bc.Read(), Q1D, D1D);
164 auto Bot = Reshape(bot.Read(), D1D-1, Q1D);
165 auto Bct = Reshape(bct.Read(), D1D, Q1D);
166 auto op = Reshape(pa_data.Read(), Q1D, Q1D, symmetric ? 3 : 4, NE);
167 auto X = Reshape(x.Read(), 2*(D1D-1)*D1D, NE);
168 auto Y = Reshape(y.ReadWrite(), 2*(D1D-1)*D1D, NE);
169
170 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
171 {
172 constexpr static int VDIM = 2;
173 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
174 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
175
176 real_t mass[MAX_Q1D][MAX_Q1D][VDIM];
177
178 for (int qy = 0; qy < Q1D; ++qy)
179 {
180 for (int qx = 0; qx < Q1D; ++qx)
181 {
182 for (int c = 0; c < VDIM; ++c)
183 {
184 mass[qy][qx][c] = 0.0;
185 }
186 }
187 }
188
189 int osc = 0;
190
191 for (int c = 0; c < VDIM; ++c) // loop over x, y components
192 {
193 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
194 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
195
196 for (int dy = 0; dy < D1Dy; ++dy)
197 {
198 real_t massX[MAX_Q1D];
199 for (int qx = 0; qx < Q1D; ++qx)
200 {
201 massX[qx] = 0.0;
202 }
203
204 for (int dx = 0; dx < D1Dx; ++dx)
205 {
206 const real_t t = X(dx + (dy * D1Dx) + osc, e);
207 for (int qx = 0; qx < Q1D; ++qx)
208 {
209 massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
210 }
211 }
212
213 for (int qy = 0; qy < Q1D; ++qy)
214 {
215 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
216 for (int qx = 0; qx < Q1D; ++qx)
217 {
218 mass[qy][qx][c] += massX[qx] * wy;
219 }
220 }
221 }
222
223 osc += D1Dx * D1Dy;
224 } // loop (c) over components
225
226 // Apply D operator.
227 for (int qy = 0; qy < Q1D; ++qy)
228 {
229 for (int qx = 0; qx < Q1D; ++qx)
230 {
231 const real_t O11 = op(qx,qy,0,e);
232 const real_t O21 = op(qx,qy,1,e);
233 const real_t O12 = symmetric ? O21 : op(qx,qy,2,e);
234 const real_t O22 = symmetric ? op(qx,qy,2,e) : op(qx,qy,3,e);
235 const real_t massX = mass[qy][qx][0];
236 const real_t massY = mass[qy][qx][1];
237 mass[qy][qx][0] = (O11*massX)+(O12*massY);
238 mass[qy][qx][1] = (O21*massX)+(O22*massY);
239 }
240 }
241
242 for (int qy = 0; qy < Q1D; ++qy)
243 {
244 osc = 0;
245
246 for (int c = 0; c < VDIM; ++c) // loop over x, y components
247 {
248 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
249 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
250
251 real_t massX[MAX_D1D];
252 for (int dx = 0; dx < D1Dx; ++dx)
253 {
254 massX[dx] = 0.0;
255 }
256 for (int qx = 0; qx < Q1D; ++qx)
257 {
258 for (int dx = 0; dx < D1Dx; ++dx)
259 {
260 massX[dx] += mass[qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
261 }
262 }
263
264 for (int dy = 0; dy < D1Dy; ++dy)
265 {
266 const real_t wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
267
268 for (int dx = 0; dx < D1Dx; ++dx)
269 {
270 Y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
271 }
272 }
273
274 osc += D1Dx * D1Dy;
275 } // loop c
276 } // loop qy
277 }); // end of element loop
278}
279
280void PAHcurlMassApply3D(const int D1D,
281 const int Q1D,
282 const int NE,
283 const bool symmetric,
284 const Array<real_t> &bo,
285 const Array<real_t> &bc,
286 const Array<real_t> &bot,
287 const Array<real_t> &bct,
288 const Vector &pa_data,
289 const Vector &x,
290 Vector &y)
291{
292 MFEM_VERIFY(D1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D,
293 "Error: D1D > MAX_D1D");
294 MFEM_VERIFY(Q1D <= DeviceDofQuadLimits::Get().HCURL_MAX_Q1D,
295 "Error: Q1D > MAX_Q1D");
296 constexpr static int VDIM = 3;
297
298 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
299 auto Bc = Reshape(bc.Read(), Q1D, D1D);
300 auto Bot = Reshape(bot.Read(), D1D-1, Q1D);
301 auto Bct = Reshape(bct.Read(), D1D, Q1D);
302 auto op = Reshape(pa_data.Read(), Q1D, Q1D, Q1D, symmetric ? 6 : 9, NE);
303 auto X = Reshape(x.Read(), 3*(D1D-1)*D1D*D1D, NE);
304 auto Y = Reshape(y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
305
306 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
307 {
308 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
309 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
310
311 real_t mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
312
313 for (int qz = 0; qz < Q1D; ++qz)
314 {
315 for (int qy = 0; qy < Q1D; ++qy)
316 {
317 for (int qx = 0; qx < Q1D; ++qx)
318 {
319 for (int c = 0; c < VDIM; ++c)
320 {
321 mass[qz][qy][qx][c] = 0.0;
322 }
323 }
324 }
325 }
326
327 int osc = 0;
328
329 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
330 {
331 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
332 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
333 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
334
335 for (int dz = 0; dz < D1Dz; ++dz)
336 {
337 real_t massXY[MAX_Q1D][MAX_Q1D];
338 for (int qy = 0; qy < Q1D; ++qy)
339 {
340 for (int qx = 0; qx < Q1D; ++qx)
341 {
342 massXY[qy][qx] = 0.0;
343 }
344 }
345
346 for (int dy = 0; dy < D1Dy; ++dy)
347 {
348 real_t massX[MAX_Q1D];
349 for (int qx = 0; qx < Q1D; ++qx)
350 {
351 massX[qx] = 0.0;
352 }
353
354 for (int dx = 0; dx < D1Dx; ++dx)
355 {
356 const real_t t = X(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
357 for (int qx = 0; qx < Q1D; ++qx)
358 {
359 massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
360 }
361 }
362
363 for (int qy = 0; qy < Q1D; ++qy)
364 {
365 const real_t wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
366 for (int qx = 0; qx < Q1D; ++qx)
367 {
368 const real_t wx = massX[qx];
369 massXY[qy][qx] += wx * wy;
370 }
371 }
372 }
373
374 for (int qz = 0; qz < Q1D; ++qz)
375 {
376 const real_t wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
377 for (int qy = 0; qy < Q1D; ++qy)
378 {
379 for (int qx = 0; qx < Q1D; ++qx)
380 {
381 mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
382 }
383 }
384 }
385 }
386
387 osc += D1Dx * D1Dy * D1Dz;
388 } // loop (c) over components
389
390 // Apply D operator.
391 for (int qz = 0; qz < Q1D; ++qz)
392 {
393 for (int qy = 0; qy < Q1D; ++qy)
394 {
395 for (int qx = 0; qx < Q1D; ++qx)
396 {
397 const real_t O11 = op(qx,qy,qz,0,e);
398 const real_t O12 = op(qx,qy,qz,1,e);
399 const real_t O13 = op(qx,qy,qz,2,e);
400 const real_t O21 = symmetric ? O12 : op(qx,qy,qz,3,e);
401 const real_t O22 = symmetric ? op(qx,qy,qz,3,e) : op(qx,qy,qz,4,e);
402 const real_t O23 = symmetric ? op(qx,qy,qz,4,e) : op(qx,qy,qz,5,e);
403 const real_t O31 = symmetric ? O13 : op(qx,qy,qz,6,e);
404 const real_t O32 = symmetric ? O23 : op(qx,qy,qz,7,e);
405 const real_t O33 = symmetric ? op(qx,qy,qz,5,e) : op(qx,qy,qz,8,e);
406 const real_t massX = mass[qz][qy][qx][0];
407 const real_t massY = mass[qz][qy][qx][1];
408 const real_t massZ = mass[qz][qy][qx][2];
409 mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
410 mass[qz][qy][qx][1] = (O21*massX)+(O22*massY)+(O23*massZ);
411 mass[qz][qy][qx][2] = (O31*massX)+(O32*massY)+(O33*massZ);
412 }
413 }
414 }
415
416 for (int qz = 0; qz < Q1D; ++qz)
417 {
418 real_t massXY[MAX_D1D][MAX_D1D];
419
420 osc = 0;
421
422 for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
423 {
424 const int D1Dz = (c == 2) ? D1D - 1 : D1D;
425 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
426 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
427
428 for (int dy = 0; dy < D1Dy; ++dy)
429 {
430 for (int dx = 0; dx < D1Dx; ++dx)
431 {
432 massXY[dy][dx] = 0.0;
433 }
434 }
435 for (int qy = 0; qy < Q1D; ++qy)
436 {
437 real_t massX[MAX_D1D];
438 for (int dx = 0; dx < D1Dx; ++dx)
439 {
440 massX[dx] = 0;
441 }
442 for (int qx = 0; qx < Q1D; ++qx)
443 {
444 for (int dx = 0; dx < D1Dx; ++dx)
445 {
446 massX[dx] += mass[qz][qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
447 }
448 }
449 for (int dy = 0; dy < D1Dy; ++dy)
450 {
451 const real_t wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
452 for (int dx = 0; dx < D1Dx; ++dx)
453 {
454 massXY[dy][dx] += massX[dx] * wy;
455 }
456 }
457 }
458
459 for (int dz = 0; dz < D1Dz; ++dz)
460 {
461 const real_t wz = (c == 2) ? Bot(dz,qz) : Bct(dz,qz);
462 for (int dy = 0; dy < D1Dy; ++dy)
463 {
464 for (int dx = 0; dx < D1Dx; ++dx)
465 {
466 Y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += massXY[dy][dx] * wz;
467 }
468 }
469 }
470
471 osc += D1Dx * D1Dy * D1Dz;
472 } // loop c
473 } // loop qz
474 }); // end of element loop
475}
476
477void PACurlCurlSetup2D(const int Q1D,
478 const int NE,
479 const Array<real_t> &w,
480 const Vector &j,
481 Vector &coeff,
482 Vector &op)
483{
484 const int NQ = Q1D*Q1D;
485 auto W = w.Read();
486 auto J = Reshape(j.Read(), NQ, 2, 2, NE);
487 auto C = Reshape(coeff.Read(), NQ, NE);
488 auto y = Reshape(op.Write(), NQ, NE);
489 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
490 {
491 for (int q = 0; q < NQ; ++q)
492 {
493 const real_t J11 = J(q,0,0,e);
494 const real_t J21 = J(q,1,0,e);
495 const real_t J12 = J(q,0,1,e);
496 const real_t J22 = J(q,1,1,e);
497 const real_t detJ = (J11*J22)-(J21*J12);
498 y(q,e) = W[q] * C(q,e) / detJ;
499 }
500 });
501}
502
503void PACurlCurlSetup3D(const int Q1D,
504 const int coeffDim,
505 const int NE,
506 const Array<real_t> &w,
507 const Vector &j,
508 Vector &coeff,
509 Vector &op)
510{
511 const int NQ = Q1D*Q1D*Q1D;
512 const bool symmetric = (coeffDim != 9);
513 auto W = w.Read();
514 auto J = Reshape(j.Read(), NQ, 3, 3, NE);
515 auto C = Reshape(coeff.Read(), coeffDim, NQ, NE);
516 auto y = Reshape(op.Write(), NQ, symmetric ? 6 : 9, NE);
517
518 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
519 {
520 for (int q = 0; q < NQ; ++q)
521 {
522 const real_t J11 = J(q,0,0,e);
523 const real_t J21 = J(q,1,0,e);
524 const real_t J31 = J(q,2,0,e);
525 const real_t J12 = J(q,0,1,e);
526 const real_t J22 = J(q,1,1,e);
527 const real_t J32 = J(q,2,1,e);
528 const real_t J13 = J(q,0,2,e);
529 const real_t J23 = J(q,1,2,e);
530 const real_t J33 = J(q,2,2,e);
531 const real_t detJ = J11 * (J22 * J33 - J32 * J23) -
532 J21 * (J12 * J33 - J32 * J13) +
533 J31 * (J12 * J23 - J22 * J13);
534
535 const real_t c_detJ = W[q] / detJ;
536
537 if (coeffDim == 6 || coeffDim == 9) // Matrix coefficient version
538 {
539 // Set y to the 6 or 9 entries of J^T M J / det
540 const real_t M11 = C(0, q, e);
541 const real_t M12 = C(1, q, e);
542 const real_t M13 = C(2, q, e);
543 const real_t M21 = (!symmetric) ? C(3, q, e) : M12;
544 const real_t M22 = (!symmetric) ? C(4, q, e) : C(3, q, e);
545 const real_t M23 = (!symmetric) ? C(5, q, e) : C(4, q, e);
546 const real_t M31 = (!symmetric) ? C(6, q, e) : M13;
547 const real_t M32 = (!symmetric) ? C(7, q, e) : M23;
548 const real_t M33 = (!symmetric) ? C(8, q, e) : C(5, q, e);
549
550 // First compute R = MJ
551 const real_t R11 = M11*J11 + M12*J21 + M13*J31;
552 const real_t R12 = M11*J12 + M12*J22 + M13*J32;
553 const real_t R13 = M11*J13 + M12*J23 + M13*J33;
554 const real_t R21 = M21*J11 + M22*J21 + M23*J31;
555 const real_t R22 = M21*J12 + M22*J22 + M23*J32;
556 const real_t R23 = M21*J13 + M22*J23 + M23*J33;
557 const real_t R31 = M31*J11 + M32*J21 + M33*J31;
558 const real_t R32 = M31*J12 + M32*J22 + M33*J32;
559 const real_t R33 = M31*J13 + M32*J23 + M33*J33;
560
561 // Now set y to J^T R / det
562 y(q,0,e) = c_detJ * (J11*R11 + J21*R21 + J31*R31); // 1,1
563 const real_t Y12 = c_detJ * (J11*R12 + J21*R22 + J31*R32);
564 y(q,1,e) = Y12; // 1,2
565 y(q,2,e) = c_detJ * (J11*R13 + J21*R23 + J31*R33); // 1,3
566
567 const real_t Y21 = c_detJ * (J12*R11 + J22*R21 + J32*R31);
568 const real_t Y22 = c_detJ * (J12*R12 + J22*R22 + J32*R32);
569 const real_t Y23 = c_detJ * (J12*R13 + J22*R23 + J32*R33);
570
571 const real_t Y33 = c_detJ * (J13*R13 + J23*R23 + J33*R33);
572
573 y(q,3,e) = symmetric ? Y22 : Y21; // 2,2 or 2,1
574 y(q,4,e) = symmetric ? Y23 : Y22; // 2,3 or 2,2
575 y(q,5,e) = symmetric ? Y33 : Y23; // 3,3 or 2,3
576
577 if (!symmetric)
578 {
579 y(q,6,e) = c_detJ * (J13*R11 + J23*R21 + J33*R31); // 3,1
580 y(q,7,e) = c_detJ * (J13*R12 + J23*R22 + J33*R32); // 3,2
581 y(q,8,e) = Y33; // 3,3
582 }
583 }
584 else // Vector or scalar coefficient version
585 {
586 // Set y to the 6 entries of J^T D J / det^2
587 const real_t D1 = C(0, q, e);
588 const real_t D2 = coeffDim == 3 ? C(1, q, e) : D1;
589 const real_t D3 = coeffDim == 3 ? C(2, q, e) : D1;
590
591 y(q,0,e) = c_detJ * (D1*J11*J11 + D2*J21*J21 + D3*J31*J31); // 1,1
592 y(q,1,e) = c_detJ * (D1*J11*J12 + D2*J21*J22 + D3*J31*J32); // 1,2
593 y(q,2,e) = c_detJ * (D1*J11*J13 + D2*J21*J23 + D3*J31*J33); // 1,3
594 y(q,3,e) = c_detJ * (D1*J12*J12 + D2*J22*J22 + D3*J32*J32); // 2,2
595 y(q,4,e) = c_detJ * (D1*J12*J13 + D2*J22*J23 + D3*J32*J33); // 2,3
596 y(q,5,e) = c_detJ * (D1*J13*J13 + D2*J23*J23 + D3*J33*J33); // 3,3
597 }
598 }
599 });
600}
601
602void PACurlCurlAssembleDiagonal2D(const int D1D, const int Q1D, const bool,
603 const int NE, const Array<real_t> &bo,
604 const Array<real_t> &, const Array<real_t> &,
605 const Array<real_t> &gc,
606 const Vector &pa_data, Vector &diag)
607{
608 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
609 auto Gc = Reshape(gc.Read(), Q1D, D1D);
610 auto op = Reshape(pa_data.Read(), Q1D, Q1D, NE);
611 auto D = Reshape(diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
612
613 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
614 {
615 constexpr static int VDIM = 2;
616 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
617
618 int osc = 0;
619
620 for (int c = 0; c < VDIM; ++c) // loop over x, y components
621 {
622 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
623 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
624
625 real_t t[MAX_Q1D];
626
627 for (int dy = 0; dy < D1Dy; ++dy)
628 {
629 for (int qx = 0; qx < Q1D; ++qx)
630 {
631 t[qx] = 0.0;
632 for (int qy = 0; qy < Q1D; ++qy)
633 {
634 const real_t wy = (c == 1) ? Bo(qy,dy) : -Gc(qy,dy);
635 t[qx] += wy * wy * op(qx,qy,e);
636 }
637 }
638
639 for (int dx = 0; dx < D1Dx; ++dx)
640 {
641 for (int qx = 0; qx < Q1D; ++qx)
642 {
643 const real_t wx = ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
644 D(dx + (dy * D1Dx) + osc, e) += t[qx] * wx * wx;
645 }
646 }
647 }
648
649 osc += D1Dx * D1Dy;
650 } // loop c
651 }); // end of element loop
652}
653
654void PACurlCurlApply2D(const int D1D, const int Q1D, const bool, const int NE,
655 const Array<real_t> &bo, const Array<real_t> &,
656 const Array<real_t> &bot, const Array<real_t> &,
657 const Array<real_t> &gc, const Array<real_t> &gct,
658 const Vector &pa_data, const Vector &x, Vector &y,
659 const bool useAbs)
660{
661
662 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
663 auto Bot = Reshape(bot.Read(), D1D-1, Q1D);
664 auto Gc = Reshape(gc.Read(), Q1D, D1D);
665 auto Gct = Reshape(gct.Read(), D1D, Q1D);
666 auto op = Reshape(pa_data.Read(), Q1D, Q1D, NE);
667 auto X = Reshape(x.Read(), 2*(D1D-1)*D1D, NE);
668 auto Y = Reshape(y.ReadWrite(), 2*(D1D-1)*D1D, NE);
669
670 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
671 {
672 constexpr static int VDIM = 2;
673 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
674 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
675
676 real_t curl[MAX_Q1D][MAX_Q1D];
677
678 // curl[qy][qx] will be computed as du_y/dx - du_x/dy
679
680 for (int qy = 0; qy < Q1D; ++qy)
681 {
682 for (int qx = 0; qx < Q1D; ++qx)
683 {
684 curl[qy][qx] = 0.0;
685 }
686 }
687
688 int osc = 0;
689
690 for (int c = 0; c < VDIM; ++c) // loop over x, y components
691 {
692 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
693 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
694
695 for (int dy = 0; dy < D1Dy; ++dy)
696 {
697 real_t gradX[MAX_Q1D];
698 for (int qx = 0; qx < Q1D; ++qx)
699 {
700 gradX[qx] = 0;
701 }
702
703 for (int dx = 0; dx < D1Dx; ++dx)
704 {
705 const real_t t = X(dx + (dy * D1Dx) + osc, e);
706 for (int qx = 0; qx < Q1D; ++qx)
707 {
708 gradX[qx] += t * ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
709 }
710 }
711
712 for (int qy = 0; qy < Q1D; ++qy)
713 {
714 const int sign = useAbs ? 1 : -1;
715 const real_t wy = (c == 0) ? (sign*Gc(qy,dy)) : Bo(qy,dy);
716 for (int qx = 0; qx < Q1D; ++qx)
717 {
718 curl[qy][qx] += gradX[qx] * wy;
719 }
720 }
721 }
722
723 osc += D1Dx * D1Dy;
724 } // loop (c) over components
725
726 // Apply D operator.
727 for (int qy = 0; qy < Q1D; ++qy)
728 {
729 for (int qx = 0; qx < Q1D; ++qx)
730 {
731 curl[qy][qx] *= op(qx,qy,e);
732 }
733 }
734
735 for (int qy = 0; qy < Q1D; ++qy)
736 {
737 osc = 0;
738
739 for (int c = 0; c < VDIM; ++c) // loop over x, y components
740 {
741 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
742 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
743
744 real_t gradX[MAX_D1D];
745 for (int dx = 0; dx < D1Dx; ++dx)
746 {
747 gradX[dx] = 0.0;
748 }
749 for (int qx = 0; qx < Q1D; ++qx)
750 {
751 for (int dx = 0; dx < D1Dx; ++dx)
752 {
753 gradX[dx] += curl[qy][qx] * ((c == 0) ? Bot(dx,qx) : Gct(dx,qx));
754 }
755 }
756 for (int dy = 0; dy < D1Dy; ++dy)
757 {
758 const int sign = useAbs ? 1 : -1;
759 const real_t wy = (c == 0) ? (sign*Gct(dy,qy)) : Bot(dy,qy);
760
761 for (int dx = 0; dx < D1Dx; ++dx)
762 {
763 Y(dx + (dy * D1Dx) + osc, e) += gradX[dx] * wy;
764 }
765 }
766
767 osc += D1Dx * D1Dy;
768 } // loop c
769 } // loop qy
770 }); // end of element loop
771}
772
773void PAHcurlL2Setup2D(const int Q1D,
774 const int NE,
775 const Array<real_t> &w,
776 Vector &coeff,
777 Vector &op)
778{
779 const int NQ = Q1D*Q1D;
780 auto W = w.Read();
781 auto C = Reshape(coeff.Read(), NQ, NE);
782 auto y = Reshape(op.Write(), NQ, NE);
783 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
784 {
785 for (int q = 0; q < NQ; ++q)
786 {
787 y(q,e) = W[q] * C(q,e);
788 }
789 });
790}
791
792void PAHcurlL2Setup3D(const int NQ,
793 const int coeffDim,
794 const int NE,
795 const Array<real_t> &w,
796 Vector &coeff,
797 Vector &op)
798{
799 auto W = w.Read();
800 auto C = Reshape(coeff.Read(), coeffDim, NQ, NE);
801 auto y = Reshape(op.Write(), coeffDim, NQ, NE);
802
803 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
804 {
805 for (int q = 0; q < NQ; ++q)
806 {
807 for (int c=0; c<coeffDim; ++c)
808 {
809 y(c,q,e) = W[q] * C(c,q,e);
810 }
811 }
812 });
813}
814
815void PAHcurlL2Apply2D(const int D1D,
816 const int D1Dtest,
817 const int Q1D,
818 const int NE,
819 const Array<real_t> &bo,
820 const Array<real_t> &bot,
821 const Array<real_t> &bt,
822 const Array<real_t> &gc,
823 const Vector &pa_data,
824 const Vector &x, // trial = H(curl)
825 Vector &y) // test = L2 or H1
826{
827 const int H1 = (D1Dtest == D1D);
828
829 MFEM_VERIFY(y.Size() == NE*D1Dtest*D1Dtest, "Test vector of wrong dimension");
830
831 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
832 auto Bot = Reshape(bot.Read(), D1D-1, Q1D);
833 auto Bt = Reshape(bt.Read(), D1D, Q1D);
834 auto Gc = Reshape(gc.Read(), Q1D, D1D);
835 auto op = Reshape(pa_data.Read(), Q1D, Q1D, NE);
836 auto X = Reshape(x.Read(), 2*(D1D-1)*D1D, NE);
837 auto Y = Reshape(y.ReadWrite(), D1Dtest, D1Dtest, NE);
838
839 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
840 {
841 constexpr static int VDIM = 2;
842 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
843 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
844
845 real_t curl[MAX_Q1D][MAX_Q1D];
846
847 // curl[qy][qx] will be computed as du_y/dx - du_x/dy
848
849 for (int qy = 0; qy < Q1D; ++qy)
850 {
851 for (int qx = 0; qx < Q1D; ++qx)
852 {
853 curl[qy][qx] = 0.0;
854 }
855 }
856
857 int osc = 0;
858
859 for (int c = 0; c < VDIM; ++c) // loop over x, y components
860 {
861 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
862 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
863
864 for (int dy = 0; dy < D1Dy; ++dy)
865 {
866 real_t gradX[MAX_Q1D];
867 for (int qx = 0; qx < Q1D; ++qx)
868 {
869 gradX[qx] = 0;
870 }
871
872 for (int dx = 0; dx < D1Dx; ++dx)
873 {
874 const real_t t = X(dx + (dy * D1Dx) + osc, e);
875 for (int qx = 0; qx < Q1D; ++qx)
876 {
877 gradX[qx] += t * ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
878 }
879 }
880
881 for (int qy = 0; qy < Q1D; ++qy)
882 {
883 const real_t wy = (c == 0) ? -Gc(qy,dy) : Bo(qy,dy);
884 for (int qx = 0; qx < Q1D; ++qx)
885 {
886 curl[qy][qx] += gradX[qx] * wy;
887 }
888 }
889 }
890
891 osc += D1Dx * D1Dy;
892 } // loop (c) over components
893
894 // Apply D operator.
895 for (int qy = 0; qy < Q1D; ++qy)
896 {
897 for (int qx = 0; qx < Q1D; ++qx)
898 {
899 curl[qy][qx] *= op(qx,qy,e);
900 }
901 }
902
903 for (int qy = 0; qy < Q1D; ++qy)
904 {
905 real_t sol_x[MAX_D1D];
906 for (int dx = 0; dx < D1Dtest; ++dx)
907 {
908 sol_x[dx] = 0.0;
909 }
910 for (int qx = 0; qx < Q1D; ++qx)
911 {
912 const real_t s = curl[qy][qx];
913 for (int dx = 0; dx < D1Dtest; ++dx)
914 {
915 sol_x[dx] += s * ((H1 == 1) ? Bt(dx,qx) : Bot(dx,qx));
916 }
917 }
918 for (int dy = 0; dy < D1Dtest; ++dy)
919 {
920 const real_t wy = (H1 == 1) ? Bt(dy,qy) : Bot(dy,qy);
921
922 for (int dx = 0; dx < D1Dtest; ++dx)
923 {
924 Y(dx,dy,e) += sol_x[dx] * wy;
925 }
926 }
927 } // loop qy
928 }); // end of element loop
929}
930
931void PAHcurlL2ApplyTranspose2D(const int D1D,
932 const int D1Dtest,
933 const int Q1D,
934 const int NE,
935 const Array<real_t> &bo,
936 const Array<real_t> &bot,
937 const Array<real_t> &b,
938 const Array<real_t> &gct,
939 const Vector &pa_data,
940 const Vector &x, // trial = H(curl)
941 Vector &y) // test = L2 or H1
942{
943 const int H1 = (D1Dtest == D1D);
944
945 MFEM_VERIFY(x.Size() == NE*D1Dtest*D1Dtest, "Test vector of wrong dimension");
946
947 auto Bo = Reshape(bo.Read(), Q1D, D1D-1);
948 auto B = Reshape(b.Read(), Q1D, D1D);
949 auto Bot = Reshape(bot.Read(), D1D-1, Q1D);
950 auto Gct = Reshape(gct.Read(), D1D, Q1D);
951 auto op = Reshape(pa_data.Read(), Q1D, Q1D, NE);
952 auto X = Reshape(x.Read(), D1Dtest, D1Dtest, NE);
953 auto Y = Reshape(y.ReadWrite(), 2*(D1D-1)*D1D, NE);
954
955 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
956 {
957 constexpr static int VDIM = 2;
958 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
959 constexpr static int MAX_Q1D = DofQuadLimits::HCURL_MAX_Q1D;
960
961 real_t mass[MAX_Q1D][MAX_Q1D];
962
963 // Zero-order term in L2 or H1 test space
964
965 for (int qy = 0; qy < Q1D; ++qy)
966 {
967 for (int qx = 0; qx < Q1D; ++qx)
968 {
969 mass[qy][qx] = 0.0;
970 }
971 }
972
973 for (int dy = 0; dy < D1Dtest; ++dy)
974 {
975 real_t sol_x[MAX_Q1D];
976 for (int qy = 0; qy < Q1D; ++qy)
977 {
978 sol_x[qy] = 0.0;
979 }
980 for (int dx = 0; dx < D1Dtest; ++dx)
981 {
982 const real_t s = X(dx,dy,e);
983 for (int qx = 0; qx < Q1D; ++qx)
984 {
985 sol_x[qx] += s * ((H1 == 1) ? B(qx,dx) : Bo(qx,dx));
986 }
987 }
988 for (int qy = 0; qy < Q1D; ++qy)
989 {
990 const real_t d2q = (H1 == 1) ? B(qy,dy) : Bo(qy,dy);
991 for (int qx = 0; qx < Q1D; ++qx)
992 {
993 mass[qy][qx] += d2q * sol_x[qx];
994 }
995 }
996 }
997
998 // Apply D operator.
999 for (int qy = 0; qy < Q1D; ++qy)
1000 {
1001 for (int qx = 0; qx < Q1D; ++qx)
1002 {
1003 mass[qy][qx] *= op(qx,qy,e);
1004 }
1005 }
1006
1007 for (int qy = 0; qy < Q1D; ++qy)
1008 {
1009 int osc = 0;
1010
1011 for (int c = 0; c < VDIM; ++c) // loop over x, y components
1012 {
1013 const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1014 const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1015
1016 real_t gradX[MAX_D1D];
1017 for (int dx = 0; dx < D1Dx; ++dx)
1018 {
1019 gradX[dx] = 0.0;
1020 }
1021 for (int qx = 0; qx < Q1D; ++qx)
1022 {
1023 for (int dx = 0; dx < D1Dx; ++dx)
1024 {
1025 gradX[dx] += mass[qy][qx] * ((c == 0) ? Bot(dx,qx) : Gct(dx,qx));
1026 }
1027 }
1028 for (int dy = 0; dy < D1Dy; ++dy)
1029 {
1030 const real_t wy = (c == 0) ? -Gct(dy,qy) : Bot(dy,qy);
1031
1032 for (int dx = 0; dx < D1Dx; ++dx)
1033 {
1034 Y(dx + (dy * D1Dx) + osc, e) += gradX[dx] * wy;
1035 }
1036 }
1037
1038 osc += D1Dx * D1Dy;
1039 } // loop c
1040 } // loop qy
1041 }); // end of element loop
1042}
1043
1044} // namespace internal
1045
1046} // namespace mfem
real_t b
Definition lissajous.cpp:42
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