MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
bilininteg_interp_pa.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#include "../bilininteg.hpp"
14#include "../gridfunc.hpp"
15#include "../qfunction.hpp"
16
17namespace mfem
18{
19
20// Apply to x corresponding to DOFs in H^1 (domain) the (topological) gradient
21// to get a dof in H(curl) (range). You can think of the range as the "test" space
22// and the domain as the "trial" space, but there's no integration.
23static void PAHcurlApplyGradient2D(const int c_dofs1D,
24 const int o_dofs1D,
25 const int NE,
26 const Array<real_t> &B_,
27 const Array<real_t> &G_,
28 const Vector &x_,
29 Vector &y_)
30{
31 auto B = Reshape(B_.Read(), c_dofs1D, c_dofs1D);
32 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
33
34 auto x = Reshape(x_.Read(), c_dofs1D, c_dofs1D, NE);
35 auto y = Reshape(y_.ReadWrite(), 2 * c_dofs1D * o_dofs1D, NE);
36
37 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().MAX_D1D &&
38 o_dofs1D <= c_dofs1D, "");
39
40 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
41 {
42 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
43 real_t w[MAX_D1D][MAX_D1D];
44
45 // horizontal part
46 for (int dx = 0; dx < c_dofs1D; ++dx)
47 {
48 for (int ey = 0; ey < c_dofs1D; ++ey)
49 {
50 w[dx][ey] = 0.0;
51 for (int dy = 0; dy < c_dofs1D; ++dy)
52 {
53 w[dx][ey] += B(ey, dy) * x(dx, dy, e);
54 }
55 }
56 }
57
58 for (int ey = 0; ey < c_dofs1D; ++ey)
59 {
60 for (int ex = 0; ex < o_dofs1D; ++ex)
61 {
62 real_t s = 0.0;
63 for (int dx = 0; dx < c_dofs1D; ++dx)
64 {
65 s += G(ex, dx) * w[dx][ey];
66 }
67 const int local_index = ey*o_dofs1D + ex;
68 y(local_index, e) += s;
69 }
70 }
71
72 // vertical part
73 for (int dx = 0; dx < c_dofs1D; ++dx)
74 {
75 for (int ey = 0; ey < o_dofs1D; ++ey)
76 {
77 w[dx][ey] = 0.0;
78 for (int dy = 0; dy < c_dofs1D; ++dy)
79 {
80 w[dx][ey] += G(ey, dy) * x(dx, dy, e);
81 }
82 }
83 }
84
85 for (int ey = 0; ey < o_dofs1D; ++ey)
86 {
87 for (int ex = 0; ex < c_dofs1D; ++ex)
88 {
89 real_t s = 0.0;
90 for (int dx = 0; dx < c_dofs1D; ++dx)
91 {
92 s += B(ex, dx) * w[dx][ey];
93 }
94 const int local_index = c_dofs1D * o_dofs1D + ey*c_dofs1D + ex;
95 y(local_index, e) += s;
96 }
97 }
98 });
99}
100
101// Specialization of PAHcurlApplyGradient2D to the case where B is identity
102static void PAHcurlApplyGradient2DBId(const int c_dofs1D,
103 const int o_dofs1D,
104 const int NE,
105 const Array<real_t> &G_,
106 const Vector &x_,
107 Vector &y_)
108{
109 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
110
111 auto x = Reshape(x_.Read(), c_dofs1D, c_dofs1D, NE);
112 auto y = Reshape(y_.ReadWrite(), 2 * c_dofs1D * o_dofs1D, NE);
113
114 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().MAX_D1D &&
115 o_dofs1D <= c_dofs1D, "");
116
117 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
118 {
119 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
120 real_t w[MAX_D1D][MAX_D1D];
121
122 // horizontal part
123 for (int dx = 0; dx < c_dofs1D; ++dx)
124 {
125 for (int ey = 0; ey < c_dofs1D; ++ey)
126 {
127 const int dy = ey;
128 w[dx][ey] = x(dx, dy, e);
129 }
130 }
131
132 for (int ey = 0; ey < c_dofs1D; ++ey)
133 {
134 for (int ex = 0; ex < o_dofs1D; ++ex)
135 {
136 real_t s = 0.0;
137 for (int dx = 0; dx < c_dofs1D; ++dx)
138 {
139 s += G(ex, dx) * w[dx][ey];
140 }
141 const int local_index = ey*o_dofs1D + ex;
142 y(local_index, e) += s;
143 }
144 }
145
146 // vertical part
147 for (int dx = 0; dx < c_dofs1D; ++dx)
148 {
149 for (int ey = 0; ey < o_dofs1D; ++ey)
150 {
151 w[dx][ey] = 0.0;
152 for (int dy = 0; dy < c_dofs1D; ++dy)
153 {
154 w[dx][ey] += G(ey, dy) * x(dx, dy, e);
155 }
156 }
157 }
158
159 for (int ey = 0; ey < o_dofs1D; ++ey)
160 {
161 for (int ex = 0; ex < c_dofs1D; ++ex)
162 {
163 const int dx = ex;
164 const real_t s = w[dx][ey];
165 const int local_index = c_dofs1D * o_dofs1D + ey*c_dofs1D + ex;
166 y(local_index, e) += s;
167 }
168 }
169 });
170}
171
172static void PAHcurlApplyGradientTranspose2D(
173 const int c_dofs1D, const int o_dofs1D, const int NE,
174 const Array<real_t> &B_, const Array<real_t> &G_,
175 const Vector &x_, Vector &y_)
176{
177 auto B = Reshape(B_.Read(), c_dofs1D, c_dofs1D);
178 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
179
180 auto x = Reshape(x_.Read(), 2 * c_dofs1D * o_dofs1D, NE);
181 auto y = Reshape(y_.ReadWrite(), c_dofs1D, c_dofs1D, NE);
182
183 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
184 o_dofs1D <= c_dofs1D, "");
185
186 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
187 {
188 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
189 real_t w[MAX_D1D][MAX_D1D];
190
191 // horizontal part (open x, closed y)
192 for (int dy = 0; dy < c_dofs1D; ++dy)
193 {
194 for (int ex = 0; ex < o_dofs1D; ++ex)
195 {
196 w[dy][ex] = 0.0;
197 for (int ey = 0; ey < c_dofs1D; ++ey)
198 {
199 const int local_index = ey*o_dofs1D + ex;
200 w[dy][ex] += B(ey, dy) * x(local_index, e);
201 }
202 }
203 }
204
205 for (int dy = 0; dy < c_dofs1D; ++dy)
206 {
207 for (int dx = 0; dx < c_dofs1D; ++dx)
208 {
209 real_t s = 0.0;
210 for (int ex = 0; ex < o_dofs1D; ++ex)
211 {
212 s += G(ex, dx) * w[dy][ex];
213 }
214 y(dx, dy, e) += s;
215 }
216 }
217
218 // vertical part (open y, closed x)
219 for (int dy = 0; dy < c_dofs1D; ++dy)
220 {
221 for (int ex = 0; ex < c_dofs1D; ++ex)
222 {
223 w[dy][ex] = 0.0;
224 for (int ey = 0; ey < o_dofs1D; ++ey)
225 {
226 const int local_index = c_dofs1D * o_dofs1D + ey*c_dofs1D + ex;
227 w[dy][ex] += G(ey, dy) * x(local_index, e);
228 }
229 }
230 }
231
232 for (int dy = 0; dy < c_dofs1D; ++dy)
233 {
234 for (int dx = 0; dx < c_dofs1D; ++dx)
235 {
236 real_t s = 0.0;
237 for (int ex = 0; ex < c_dofs1D; ++ex)
238 {
239 s += B(ex, dx) * w[dy][ex];
240 }
241 y(dx, dy, e) += s;
242 }
243 }
244 });
245}
246
247// Specialization of PAHcurlApplyGradientTranspose2D to the case where
248// B is identity
249static void PAHcurlApplyGradientTranspose2DBId(
250 const int c_dofs1D, const int o_dofs1D, const int NE,
251 const Array<real_t> &G_,
252 const Vector &x_, Vector &y_)
253{
254 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
255
256 auto x = Reshape(x_.Read(), 2 * c_dofs1D * o_dofs1D, NE);
257 auto y = Reshape(y_.ReadWrite(), c_dofs1D, c_dofs1D, NE);
258
259 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
260 o_dofs1D <= c_dofs1D, "");
261
262 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
263 {
264 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
265 real_t w[MAX_D1D][MAX_D1D];
266
267 // horizontal part (open x, closed y)
268 for (int dy = 0; dy < c_dofs1D; ++dy)
269 {
270 for (int ex = 0; ex < o_dofs1D; ++ex)
271 {
272 const int ey = dy;
273 const int local_index = ey*o_dofs1D + ex;
274 w[dy][ex] = x(local_index, e);
275 }
276 }
277
278 for (int dy = 0; dy < c_dofs1D; ++dy)
279 {
280 for (int dx = 0; dx < c_dofs1D; ++dx)
281 {
282 real_t s = 0.0;
283 for (int ex = 0; ex < o_dofs1D; ++ex)
284 {
285 s += G(ex, dx) * w[dy][ex];
286 }
287 y(dx, dy, e) += s;
288 }
289 }
290
291 // vertical part (open y, closed x)
292 for (int dy = 0; dy < c_dofs1D; ++dy)
293 {
294 for (int ex = 0; ex < c_dofs1D; ++ex)
295 {
296 w[dy][ex] = 0.0;
297 for (int ey = 0; ey < o_dofs1D; ++ey)
298 {
299 const int local_index = c_dofs1D * o_dofs1D + ey*c_dofs1D + ex;
300 w[dy][ex] += G(ey, dy) * x(local_index, e);
301 }
302 }
303 }
304
305 for (int dy = 0; dy < c_dofs1D; ++dy)
306 {
307 for (int dx = 0; dx < c_dofs1D; ++dx)
308 {
309 const int ex = dx;
310 const real_t s = w[dy][ex];
311 y(dx, dy, e) += s;
312 }
313 }
314 });
315}
316
317static void PAHcurlApplyGradient3D(const int c_dofs1D,
318 const int o_dofs1D,
319 const int NE,
320 const Array<real_t> &B_,
321 const Array<real_t> &G_,
322 const Vector &x_,
323 Vector &y_)
324{
325 auto B = Reshape(B_.Read(), c_dofs1D, c_dofs1D);
326 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
327
328 auto x = Reshape(x_.Read(), c_dofs1D, c_dofs1D, c_dofs1D, NE);
329 auto y = Reshape(y_.ReadWrite(), (3 * c_dofs1D * c_dofs1D * o_dofs1D), NE);
330
331 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
332 o_dofs1D <= c_dofs1D, "");
333
334 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
335 {
336 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
337 real_t w1[MAX_D1D][MAX_D1D][MAX_D1D];
338 real_t w2[MAX_D1D][MAX_D1D][MAX_D1D];
339
340 // ---
341 // dofs that point parallel to x-axis (open in x, closed in y, z)
342 // ---
343
344 // contract in z
345 for (int ez = 0; ez < c_dofs1D; ++ez)
346 {
347 for (int dx = 0; dx < c_dofs1D; ++dx)
348 {
349 for (int dy = 0; dy < c_dofs1D; ++dy)
350 {
351 w1[dx][dy][ez] = 0.0;
352 for (int dz = 0; dz < c_dofs1D; ++dz)
353 {
354 w1[dx][dy][ez] += B(ez, dz) * x(dx, dy, dz, e);
355 }
356 }
357 }
358 }
359
360 // contract in y
361 for (int ez = 0; ez < c_dofs1D; ++ez)
362 {
363 for (int ey = 0; ey < c_dofs1D; ++ey)
364 {
365 for (int dx = 0; dx < c_dofs1D; ++dx)
366 {
367 w2[dx][ey][ez] = 0.0;
368 for (int dy = 0; dy < c_dofs1D; ++dy)
369 {
370 w2[dx][ey][ez] += B(ey, dy) * w1[dx][dy][ez];
371 }
372 }
373 }
374 }
375
376 // contract in x
377 for (int ez = 0; ez < c_dofs1D; ++ez)
378 {
379 for (int ey = 0; ey < c_dofs1D; ++ey)
380 {
381 for (int ex = 0; ex < o_dofs1D; ++ex)
382 {
383 real_t s = 0.0;
384 for (int dx = 0; dx < c_dofs1D; ++dx)
385 {
386 s += G(ex, dx) * w2[dx][ey][ez];
387 }
388 const int local_index = ez*c_dofs1D*o_dofs1D + ey*o_dofs1D + ex;
389 y(local_index, e) += s;
390 }
391 }
392 }
393
394 // ---
395 // dofs that point parallel to y-axis (open in y, closed in x, z)
396 // ---
397
398 // contract in z
399 for (int ez = 0; ez < c_dofs1D; ++ez)
400 {
401 for (int dx = 0; dx < c_dofs1D; ++dx)
402 {
403 for (int dy = 0; dy < c_dofs1D; ++dy)
404 {
405 w1[dx][dy][ez] = 0.0;
406 for (int dz = 0; dz < c_dofs1D; ++dz)
407 {
408 w1[dx][dy][ez] += B(ez, dz) * x(dx, dy, dz, e);
409 }
410 }
411 }
412 }
413
414 // contract in y
415 for (int ez = 0; ez < c_dofs1D; ++ez)
416 {
417 for (int ey = 0; ey < o_dofs1D; ++ey)
418 {
419 for (int dx = 0; dx < c_dofs1D; ++dx)
420 {
421 w2[dx][ey][ez] = 0.0;
422 for (int dy = 0; dy < c_dofs1D; ++dy)
423 {
424 w2[dx][ey][ez] += G(ey, dy) * w1[dx][dy][ez];
425 }
426 }
427 }
428 }
429
430 // contract in x
431 for (int ez = 0; ez < c_dofs1D; ++ez)
432 {
433 for (int ey = 0; ey < o_dofs1D; ++ey)
434 {
435 for (int ex = 0; ex < c_dofs1D; ++ex)
436 {
437 real_t s = 0.0;
438 for (int dx = 0; dx < c_dofs1D; ++dx)
439 {
440 s += B(ex, dx) * w2[dx][ey][ez];
441 }
442 const int local_index = c_dofs1D*c_dofs1D*o_dofs1D +
443 ez*c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
444 y(local_index, e) += s;
445 }
446 }
447 }
448
449 // ---
450 // dofs that point parallel to z-axis (open in z, closed in x, y)
451 // ---
452
453 // contract in z
454 for (int ez = 0; ez < o_dofs1D; ++ez)
455 {
456 for (int dx = 0; dx < c_dofs1D; ++dx)
457 {
458 for (int dy = 0; dy < c_dofs1D; ++dy)
459 {
460 w1[dx][dy][ez] = 0.0;
461 for (int dz = 0; dz < c_dofs1D; ++dz)
462 {
463 w1[dx][dy][ez] += G(ez, dz) * x(dx, dy, dz, e);
464 }
465 }
466 }
467 }
468
469 // contract in y
470 for (int ez = 0; ez < o_dofs1D; ++ez)
471 {
472 for (int ey = 0; ey < c_dofs1D; ++ey)
473 {
474 for (int dx = 0; dx < c_dofs1D; ++dx)
475 {
476 w2[dx][ey][ez] = 0.0;
477 for (int dy = 0; dy < c_dofs1D; ++dy)
478 {
479 w2[dx][ey][ez] += B(ey, dy) * w1[dx][dy][ez];
480 }
481 }
482 }
483 }
484
485 // contract in x
486 for (int ez = 0; ez < o_dofs1D; ++ez)
487 {
488 for (int ey = 0; ey < c_dofs1D; ++ey)
489 {
490 for (int ex = 0; ex < c_dofs1D; ++ex)
491 {
492 real_t s = 0.0;
493 for (int dx = 0; dx < c_dofs1D; ++dx)
494 {
495 s += B(ex, dx) * w2[dx][ey][ez];
496 }
497 const int local_index = 2*c_dofs1D*c_dofs1D*o_dofs1D +
498 ez*c_dofs1D*c_dofs1D + ey*c_dofs1D + ex;
499 y(local_index, e) += s;
500 }
501 }
502 }
503 });
504}
505
506// Specialization of PAHcurlApplyGradient3D to the case where B is identity
507static void PAHcurlApplyGradient3DBId(const int c_dofs1D,
508 const int o_dofs1D,
509 const int NE,
510 const Array<real_t> &G_,
511 const Vector &x_,
512 Vector &y_)
513{
514 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
515
516 auto x = Reshape(x_.Read(), c_dofs1D, c_dofs1D, c_dofs1D, NE);
517 auto y = Reshape(y_.ReadWrite(), (3 * c_dofs1D * c_dofs1D * o_dofs1D), NE);
518
519 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
520 o_dofs1D <= c_dofs1D, "");
521
522 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
523 {
524 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
525
526 real_t w1[MAX_D1D][MAX_D1D][MAX_D1D];
527 real_t w2[MAX_D1D][MAX_D1D][MAX_D1D];
528
529 // ---
530 // dofs that point parallel to x-axis (open in x, closed in y, z)
531 // ---
532
533 // contract in z
534 for (int ez = 0; ez < c_dofs1D; ++ez)
535 {
536 for (int dx = 0; dx < c_dofs1D; ++dx)
537 {
538 for (int dy = 0; dy < c_dofs1D; ++dy)
539 {
540 const int dz = ez;
541 w1[dx][dy][ez] = x(dx, dy, dz, e);
542 }
543 }
544 }
545
546 // contract in y
547 for (int ez = 0; ez < c_dofs1D; ++ez)
548 {
549 for (int ey = 0; ey < c_dofs1D; ++ey)
550 {
551 for (int dx = 0; dx < c_dofs1D; ++dx)
552 {
553 const int dy = ey;
554 w2[dx][ey][ez] = w1[dx][dy][ez];
555 }
556 }
557 }
558
559 // contract in x
560 for (int ez = 0; ez < c_dofs1D; ++ez)
561 {
562 for (int ey = 0; ey < c_dofs1D; ++ey)
563 {
564 for (int ex = 0; ex < o_dofs1D; ++ex)
565 {
566 real_t s = 0.0;
567 for (int dx = 0; dx < c_dofs1D; ++dx)
568 {
569 s += G(ex, dx) * w2[dx][ey][ez];
570 }
571 const int local_index = ez*c_dofs1D*o_dofs1D + ey*o_dofs1D + ex;
572 y(local_index, e) += s;
573 }
574 }
575 }
576
577 // ---
578 // dofs that point parallel to y-axis (open in y, closed in x, z)
579 // ---
580
581 // contract in z
582 for (int ez = 0; ez < c_dofs1D; ++ez)
583 {
584 for (int dx = 0; dx < c_dofs1D; ++dx)
585 {
586 for (int dy = 0; dy < c_dofs1D; ++dy)
587 {
588 const int dz = ez;
589 w1[dx][dy][ez] = x(dx, dy, dz, e);
590 }
591 }
592 }
593
594 // contract in y
595 for (int ez = 0; ez < c_dofs1D; ++ez)
596 {
597 for (int ey = 0; ey < o_dofs1D; ++ey)
598 {
599 for (int dx = 0; dx < c_dofs1D; ++dx)
600 {
601 w2[dx][ey][ez] = 0.0;
602 for (int dy = 0; dy < c_dofs1D; ++dy)
603 {
604 w2[dx][ey][ez] += G(ey, dy) * w1[dx][dy][ez];
605 }
606 }
607 }
608 }
609
610 // contract in x
611 for (int ez = 0; ez < c_dofs1D; ++ez)
612 {
613 for (int ey = 0; ey < o_dofs1D; ++ey)
614 {
615 for (int ex = 0; ex < c_dofs1D; ++ex)
616 {
617 const int dx = ex;
618 const real_t s = w2[dx][ey][ez];
619 const int local_index = c_dofs1D*c_dofs1D*o_dofs1D +
620 ez*c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
621 y(local_index, e) += s;
622 }
623 }
624 }
625
626 // ---
627 // dofs that point parallel to z-axis (open in z, closed in x, y)
628 // ---
629
630 // contract in z
631 for (int ez = 0; ez < o_dofs1D; ++ez)
632 {
633 for (int dx = 0; dx < c_dofs1D; ++dx)
634 {
635 for (int dy = 0; dy < c_dofs1D; ++dy)
636 {
637 w1[dx][dy][ez] = 0.0;
638 for (int dz = 0; dz < c_dofs1D; ++dz)
639 {
640 w1[dx][dy][ez] += G(ez, dz) * x(dx, dy, dz, e);
641 }
642 }
643 }
644 }
645
646 // contract in y
647 for (int ez = 0; ez < o_dofs1D; ++ez)
648 {
649 for (int ey = 0; ey < c_dofs1D; ++ey)
650 {
651 for (int dx = 0; dx < c_dofs1D; ++dx)
652 {
653 const int dy = ey;
654 w2[dx][ey][ez] = w1[dx][dy][ez];
655 }
656 }
657 }
658
659 // contract in x
660 for (int ez = 0; ez < o_dofs1D; ++ez)
661 {
662 for (int ey = 0; ey < c_dofs1D; ++ey)
663 {
664 for (int ex = 0; ex < c_dofs1D; ++ex)
665 {
666 const int dx = ex;
667 const real_t s = w2[dx][ey][ez];
668 const int local_index = 2*c_dofs1D*c_dofs1D*o_dofs1D +
669 ez*c_dofs1D*c_dofs1D + ey*c_dofs1D + ex;
670 y(local_index, e) += s;
671 }
672 }
673 }
674 });
675}
676
677static void PAHcurlApplyGradientTranspose3D(
678 const int c_dofs1D, const int o_dofs1D, const int NE,
679 const Array<real_t> &B_, const Array<real_t> &G_,
680 const Vector &x_, Vector &y_)
681{
682 auto B = Reshape(B_.Read(), c_dofs1D, c_dofs1D);
683 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
684
685 auto x = Reshape(x_.Read(), (3 * c_dofs1D * c_dofs1D * o_dofs1D), NE);
686 auto y = Reshape(y_.ReadWrite(), c_dofs1D, c_dofs1D, c_dofs1D, NE);
687
688 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
689 o_dofs1D <= c_dofs1D, "");
690
691 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
692 {
693 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
694 real_t w1[MAX_D1D][MAX_D1D][MAX_D1D];
695 real_t w2[MAX_D1D][MAX_D1D][MAX_D1D];
696 // ---
697 // dofs that point parallel to x-axis (open in x, closed in y, z)
698 // ---
699
700 // contract in z
701 for (int dz = 0; dz < c_dofs1D; ++dz)
702 {
703 for (int ex = 0; ex < o_dofs1D; ++ex)
704 {
705 for (int ey = 0; ey < c_dofs1D; ++ey)
706 {
707 w1[ex][ey][dz] = 0.0;
708 for (int ez = 0; ez < c_dofs1D; ++ez)
709 {
710 const int local_index = ez*c_dofs1D*o_dofs1D + ey*o_dofs1D + ex;
711 w1[ex][ey][dz] += B(ez, dz) * x(local_index, e);
712 }
713 }
714 }
715 }
716
717 // contract in y
718 for (int dz = 0; dz < c_dofs1D; ++dz)
719 {
720 for (int dy = 0; dy < c_dofs1D; ++dy)
721 {
722 for (int ex = 0; ex < o_dofs1D; ++ex)
723 {
724 w2[ex][dy][dz] = 0.0;
725 for (int ey = 0; ey < c_dofs1D; ++ey)
726 {
727 w2[ex][dy][dz] += B(ey, dy) * w1[ex][ey][dz];
728 }
729 }
730 }
731 }
732
733 // contract in x
734 for (int dz = 0; dz < c_dofs1D; ++dz)
735 {
736 for (int dy = 0; dy < c_dofs1D; ++dy)
737 {
738 for (int dx = 0; dx < c_dofs1D; ++dx)
739 {
740 real_t s = 0.0;
741 for (int ex = 0; ex < o_dofs1D; ++ex)
742 {
743 s += G(ex, dx) * w2[ex][dy][dz];
744 }
745 y(dx, dy, dz, e) += s;
746 }
747 }
748 }
749
750 // ---
751 // dofs that point parallel to y-axis (open in y, closed in x, z)
752 // ---
753
754 // contract in z
755 for (int dz = 0; dz < c_dofs1D; ++dz)
756 {
757 for (int ex = 0; ex < c_dofs1D; ++ex)
758 {
759 for (int ey = 0; ey < o_dofs1D; ++ey)
760 {
761 w1[ex][ey][dz] = 0.0;
762 for (int ez = 0; ez < c_dofs1D; ++ez)
763 {
764 const int local_index = c_dofs1D*c_dofs1D*o_dofs1D +
765 ez*c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
766 w1[ex][ey][dz] += B(ez, dz) * x(local_index, e);
767 }
768 }
769 }
770 }
771
772 // contract in y
773 for (int dz = 0; dz < c_dofs1D; ++dz)
774 {
775 for (int dy = 0; dy < c_dofs1D; ++dy)
776 {
777 for (int ex = 0; ex < c_dofs1D; ++ex)
778 {
779 w2[ex][dy][dz] = 0.0;
780 for (int ey = 0; ey < o_dofs1D; ++ey)
781 {
782 w2[ex][dy][dz] += G(ey, dy) * w1[ex][ey][dz];
783 }
784 }
785 }
786 }
787
788 // contract in x
789 for (int dz = 0; dz < c_dofs1D; ++dz)
790 {
791 for (int dy = 0; dy < c_dofs1D; ++dy)
792 {
793 for (int dx = 0; dx < c_dofs1D; ++dx)
794 {
795 real_t s = 0.0;
796 for (int ex = 0; ex < c_dofs1D; ++ex)
797 {
798 s += B(ex, dx) * w2[ex][dy][dz];
799 }
800 y(dx, dy, dz, e) += s;
801 }
802 }
803 }
804
805 // ---
806 // dofs that point parallel to z-axis (open in z, closed in x, y)
807 // ---
808
809 // contract in z
810 for (int dz = 0; dz < c_dofs1D; ++dz)
811 {
812 for (int ex = 0; ex < c_dofs1D; ++ex)
813 {
814 for (int ey = 0; ey < c_dofs1D; ++ey)
815 {
816 w1[ex][ey][dz] = 0.0;
817 for (int ez = 0; ez < o_dofs1D; ++ez)
818 {
819 const int local_index = 2*c_dofs1D*c_dofs1D*o_dofs1D +
820 ez*c_dofs1D*c_dofs1D + ey*c_dofs1D + ex;
821 w1[ex][ey][dz] += G(ez, dz) * x(local_index, e);
822 }
823 }
824 }
825 }
826
827 // contract in y
828 for (int dz = 0; dz < c_dofs1D; ++dz)
829 {
830 for (int dy = 0; dy < c_dofs1D; ++dy)
831 {
832 for (int ex = 0; ex < c_dofs1D; ++ex)
833 {
834 w2[ex][dy][dz] = 0.0;
835 for (int ey = 0; ey < c_dofs1D; ++ey)
836 {
837 w2[ex][dy][dz] += B(ey, dy) * w1[ex][ey][dz];
838 }
839 }
840 }
841 }
842
843 // contract in x
844 for (int dz = 0; dz < c_dofs1D; ++dz)
845 {
846 for (int dy = 0; dy < c_dofs1D; ++dy)
847 {
848 for (int dx = 0; dx < c_dofs1D; ++dx)
849 {
850 real_t s = 0.0;
851 for (int ex = 0; ex < c_dofs1D; ++ex)
852 {
853 s += B(ex, dx) * w2[ex][dy][dz];
854 }
855 y(dx, dy, dz, e) += s;
856 }
857 }
858 }
859 });
860}
861
862// Specialization of PAHcurlApplyGradientTranspose3D to the case where
863// B is identity
864static void PAHcurlApplyGradientTranspose3DBId(
865 const int c_dofs1D, const int o_dofs1D, const int NE,
866 const Array<real_t> &G_,
867 const Vector &x_, Vector &y_)
868{
869 auto G = Reshape(G_.Read(), o_dofs1D, c_dofs1D);
870
871 auto x = Reshape(x_.Read(), (3 * c_dofs1D * c_dofs1D * o_dofs1D), NE);
872 auto y = Reshape(y_.ReadWrite(), c_dofs1D, c_dofs1D, c_dofs1D, NE);
873
874 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
875 o_dofs1D <= c_dofs1D, "");
876
877 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
878 {
879 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
880
881 real_t w1[MAX_D1D][MAX_D1D][MAX_D1D];
882 real_t w2[MAX_D1D][MAX_D1D][MAX_D1D];
883 // ---
884 // dofs that point parallel to x-axis (open in x, closed in y, z)
885 // ---
886
887 // contract in z
888 for (int dz = 0; dz < c_dofs1D; ++dz)
889 {
890 for (int ex = 0; ex < o_dofs1D; ++ex)
891 {
892 for (int ey = 0; ey < c_dofs1D; ++ey)
893 {
894 const int ez = dz;
895 const int local_index = ez*c_dofs1D*o_dofs1D + ey*o_dofs1D + ex;
896 w1[ex][ey][dz] = x(local_index, e);
897 }
898 }
899 }
900
901 // contract in y
902 for (int dz = 0; dz < c_dofs1D; ++dz)
903 {
904 for (int dy = 0; dy < c_dofs1D; ++dy)
905 {
906 for (int ex = 0; ex < o_dofs1D; ++ex)
907 {
908 const int ey = dy;
909 w2[ex][dy][dz] = w1[ex][ey][dz];
910 }
911 }
912 }
913
914 // contract in x
915 for (int dz = 0; dz < c_dofs1D; ++dz)
916 {
917 for (int dy = 0; dy < c_dofs1D; ++dy)
918 {
919 for (int dx = 0; dx < c_dofs1D; ++dx)
920 {
921 real_t s = 0.0;
922 for (int ex = 0; ex < o_dofs1D; ++ex)
923 {
924 s += G(ex, dx) * w2[ex][dy][dz];
925 }
926 y(dx, dy, dz, e) += s;
927 }
928 }
929 }
930
931 // ---
932 // dofs that point parallel to y-axis (open in y, closed in x, z)
933 // ---
934
935 // contract in z
936 for (int dz = 0; dz < c_dofs1D; ++dz)
937 {
938 for (int ex = 0; ex < c_dofs1D; ++ex)
939 {
940 for (int ey = 0; ey < o_dofs1D; ++ey)
941 {
942 const int ez = dz;
943 const int local_index = c_dofs1D*c_dofs1D*o_dofs1D +
944 ez*c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
945 w1[ex][ey][dz] = x(local_index, e);
946 }
947 }
948 }
949
950 // contract in y
951 for (int dz = 0; dz < c_dofs1D; ++dz)
952 {
953 for (int dy = 0; dy < c_dofs1D; ++dy)
954 {
955 for (int ex = 0; ex < c_dofs1D; ++ex)
956 {
957 w2[ex][dy][dz] = 0.0;
958 for (int ey = 0; ey < o_dofs1D; ++ey)
959 {
960 w2[ex][dy][dz] += G(ey, dy) * w1[ex][ey][dz];
961 }
962 }
963 }
964 }
965
966 // contract in x
967 for (int dz = 0; dz < c_dofs1D; ++dz)
968 {
969 for (int dy = 0; dy < c_dofs1D; ++dy)
970 {
971 for (int dx = 0; dx < c_dofs1D; ++dx)
972 {
973 const int ex = dx;
974 real_t s = w2[ex][dy][dz];
975 y(dx, dy, dz, e) += s;
976 }
977 }
978 }
979
980 // ---
981 // dofs that point parallel to z-axis (open in z, closed in x, y)
982 // ---
983
984 // contract in z
985 for (int dz = 0; dz < c_dofs1D; ++dz)
986 {
987 for (int ex = 0; ex < c_dofs1D; ++ex)
988 {
989 for (int ey = 0; ey < c_dofs1D; ++ey)
990 {
991 w1[ex][ey][dz] = 0.0;
992 for (int ez = 0; ez < o_dofs1D; ++ez)
993 {
994 const int local_index = 2*c_dofs1D*c_dofs1D*o_dofs1D +
995 ez*c_dofs1D*c_dofs1D + ey*c_dofs1D + ex;
996 w1[ex][ey][dz] += G(ez, dz) * x(local_index, e);
997 }
998 }
999 }
1000 }
1001
1002 // contract in y
1003 for (int dz = 0; dz < c_dofs1D; ++dz)
1004 {
1005 for (int dy = 0; dy < c_dofs1D; ++dy)
1006 {
1007 for (int ex = 0; ex < c_dofs1D; ++ex)
1008 {
1009 const int ey = dy;
1010 w2[ex][dy][dz] = w1[ex][ey][dz];
1011 }
1012 }
1013 }
1014
1015 // contract in x
1016 for (int dz = 0; dz < c_dofs1D; ++dz)
1017 {
1018 for (int dy = 0; dy < c_dofs1D; ++dy)
1019 {
1020 for (int dx = 0; dx < c_dofs1D; ++dx)
1021 {
1022 const int ex = dx;
1023 real_t s = w2[ex][dy][dz];
1024 y(dx, dy, dz, e) += s;
1025 }
1026 }
1027 }
1028 });
1029}
1030
1032 const FiniteElementSpace &test_fes)
1033{
1034 // Assumes tensor-product elements, with a vector test space and H^1 trial space.
1035 Mesh *mesh = trial_fes.GetMesh();
1036 const FiniteElement *trial_fel = trial_fes.GetTypicalFE();
1037 const FiniteElement *test_fel = test_fes.GetTypicalFE();
1038
1039 const NodalTensorFiniteElement *trial_el =
1040 dynamic_cast<const NodalTensorFiniteElement*>(trial_fel);
1041 MFEM_VERIFY(trial_el != NULL, "Only NodalTensorFiniteElement is supported!");
1042
1043 const VectorTensorFiniteElement *test_el =
1044 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
1045 MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!");
1046
1047 const int dims = trial_el->GetDim();
1048 MFEM_VERIFY(dims == 2 || dims == 3, "Bad dimension!");
1049 dim = mesh->Dimension();
1050 MFEM_VERIFY(dim == 2 || dim == 3, "Bad dimension!");
1051 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(),
1052 "Orders do not match!");
1053 ne = trial_fes.GetNE();
1054
1055 const int order = trial_el->GetOrder();
1056 dofquad_fe = new H1_SegmentElement(order, trial_el->GetBasisType());
1058 mfem::IntegrationRule closed_ir;
1059 closed_ir.SetSize(order + 1);
1060 qf1d.GaussLobatto(order + 1, &closed_ir);
1061 mfem::IntegrationRule open_ir;
1062 open_ir.SetSize(order);
1063 qf1d.GaussLegendre(order, &open_ir);
1064
1065 maps_O_C = &dofquad_fe->GetDofToQuad(open_ir, DofToQuad::TENSOR);
1066 o_dofs1D = maps_O_C->nqpt;
1067 if (trial_el->GetBasisType() == BasisType::GaussLobatto)
1068 {
1069 B_id = true;
1070 c_dofs1D = maps_O_C->ndof;
1071 }
1072 else
1073 {
1074 B_id = false;
1075 maps_C_C = &dofquad_fe->GetDofToQuad(closed_ir, DofToQuad::TENSOR);
1076 c_dofs1D = maps_C_C->nqpt;
1077 }
1078}
1079
1081{
1082 if (dim == 3)
1083 {
1084 if (B_id)
1085 {
1086 PAHcurlApplyGradient3DBId(c_dofs1D, o_dofs1D, ne,
1087 maps_O_C->G, x, y);
1088 }
1089 else
1090 {
1091 PAHcurlApplyGradient3D(c_dofs1D, o_dofs1D, ne, maps_C_C->B,
1092 maps_O_C->G, x, y);
1093 }
1094 }
1095 else if (dim == 2)
1096 {
1097 if (B_id)
1098 {
1099 PAHcurlApplyGradient2DBId(c_dofs1D, o_dofs1D, ne,
1100 maps_O_C->G, x, y);
1101 }
1102 else
1103 {
1104 PAHcurlApplyGradient2D(c_dofs1D, o_dofs1D, ne, maps_C_C->B, maps_O_C->G,
1105 x, y);
1106 }
1107 }
1108 else
1109 {
1110 mfem_error("Bad dimension!");
1111 }
1112}
1113
1115{
1116 if (dim == 3)
1117 {
1118 if (B_id)
1119 {
1120 PAHcurlApplyGradientTranspose3DBId(c_dofs1D, o_dofs1D, ne,
1121 maps_O_C->G, x, y);
1122 }
1123 else
1124 {
1125 PAHcurlApplyGradientTranspose3D(c_dofs1D, o_dofs1D, ne, maps_C_C->B,
1126 maps_O_C->G, x, y);
1127 }
1128 }
1129 else if (dim == 2)
1130 {
1131 if (B_id)
1132 {
1133 PAHcurlApplyGradientTranspose2DBId(c_dofs1D, o_dofs1D, ne,
1134 maps_O_C->G, x, y);
1135 }
1136 else
1137 {
1138 PAHcurlApplyGradientTranspose2D(c_dofs1D, o_dofs1D, ne, maps_C_C->B,
1139 maps_O_C->G, x, y);
1140 }
1141 }
1142 else
1143 {
1144 mfem_error("Bad dimension!");
1145 }
1146}
1147
1148static void PAHcurlVecH1IdentityApply2D(const int c_dofs1D,
1149 const int o_dofs1D,
1150 const int NE,
1151 const Array<real_t> &Bclosed,
1152 const Array<real_t> &Bopen,
1153 const Vector &pa_data,
1154 const Vector &x_,
1155 Vector &y_)
1156{
1157 auto Bc = Reshape(Bclosed.Read(), c_dofs1D, c_dofs1D);
1158 auto Bo = Reshape(Bopen.Read(), o_dofs1D, c_dofs1D);
1159
1160 auto x = Reshape(x_.Read(), c_dofs1D, c_dofs1D, 2, NE);
1161 auto y = Reshape(y_.ReadWrite(), (2 * c_dofs1D * o_dofs1D), NE);
1162
1163 auto vk = Reshape(pa_data.Read(), 2, (2 * c_dofs1D * o_dofs1D), NE);
1164
1165 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
1166 o_dofs1D <= c_dofs1D, "");
1167
1168 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1169 {
1170 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
1171
1172 real_t w[2][MAX_D1D][MAX_D1D];
1173
1174 // dofs that point parallel to x-axis (open in x, closed in y)
1175
1176 // contract in y
1177 for (int ey = 0; ey < c_dofs1D; ++ey)
1178 {
1179 for (int dx = 0; dx < c_dofs1D; ++dx)
1180 {
1181 for (int j=0; j<2; ++j)
1182 {
1183 w[j][dx][ey] = 0.0;
1184 for (int dy = 0; dy < c_dofs1D; ++dy)
1185 {
1186 w[j][dx][ey] += Bc(ey, dy) * x(dx, dy, j, e);
1187 }
1188 }
1189 }
1190 }
1191
1192 // contract in x
1193 for (int ey = 0; ey < c_dofs1D; ++ey)
1194 {
1195 for (int ex = 0; ex < o_dofs1D; ++ex)
1196 {
1197 for (int j=0; j<2; ++j)
1198 {
1199 real_t s = 0.0;
1200 for (int dx = 0; dx < c_dofs1D; ++dx)
1201 {
1202 s += Bo(ex, dx) * w[j][dx][ey];
1203 }
1204 const int local_index = ey*o_dofs1D + ex;
1205 y(local_index, e) += s * vk(j, local_index, e);
1206 }
1207 }
1208 }
1209
1210 // dofs that point parallel to y-axis (open in y, closed in x)
1211
1212 // contract in y
1213 for (int ey = 0; ey < o_dofs1D; ++ey)
1214 {
1215 for (int dx = 0; dx < c_dofs1D; ++dx)
1216 {
1217 for (int j=0; j<2; ++j)
1218 {
1219 w[j][dx][ey] = 0.0;
1220 for (int dy = 0; dy < c_dofs1D; ++dy)
1221 {
1222 w[j][dx][ey] += Bo(ey, dy) * x(dx, dy, j, e);
1223 }
1224 }
1225 }
1226 }
1227
1228 // contract in x
1229 for (int ey = 0; ey < o_dofs1D; ++ey)
1230 {
1231 for (int ex = 0; ex < c_dofs1D; ++ex)
1232 {
1233 for (int j=0; j<2; ++j)
1234 {
1235 real_t s = 0.0;
1236 for (int dx = 0; dx < c_dofs1D; ++dx)
1237 {
1238 s += Bc(ex, dx) * w[j][dx][ey];
1239 }
1240 const int local_index = c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
1241 y(local_index, e) += s * vk(j, local_index, e);
1242 }
1243 }
1244 }
1245 });
1246}
1247
1248static void PAHcurlVecH1IdentityApplyTranspose2D(const int c_dofs1D,
1249 const int o_dofs1D,
1250 const int NE,
1251 const Array<real_t> &Bclosed,
1252 const Array<real_t> &Bopen,
1253 const Vector &pa_data,
1254 const Vector &x_,
1255 Vector &y_)
1256{
1257 auto Bc = Reshape(Bclosed.Read(), c_dofs1D, c_dofs1D);
1258 auto Bo = Reshape(Bopen.Read(), o_dofs1D, c_dofs1D);
1259
1260 auto x = Reshape(x_.Read(), (2 * c_dofs1D * o_dofs1D), NE);
1261 auto y = Reshape(y_.ReadWrite(), c_dofs1D, c_dofs1D, 2, NE);
1262
1263 auto vk = Reshape(pa_data.Read(), 2, (2 * c_dofs1D * o_dofs1D), NE);
1264
1265 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().HCURL_MAX_D1D &&
1266 o_dofs1D <= c_dofs1D, "");
1267
1268 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1269 {
1270 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
1271
1272 real_t w[2][MAX_D1D][MAX_D1D];
1273
1274 // dofs that point parallel to x-axis (open in x, closed in y)
1275
1276 // contract in x
1277 for (int ey = 0; ey < c_dofs1D; ++ey)
1278 {
1279 for (int dx = 0; dx < c_dofs1D; ++dx)
1280 {
1281 for (int j=0; j<2; ++j) { w[j][dx][ey] = 0.0; }
1282 }
1283 for (int ex = 0; ex < o_dofs1D; ++ex)
1284 {
1285 const int local_index = ey*o_dofs1D + ex;
1286 const real_t xd = x(local_index, e);
1287
1288 for (int dx = 0; dx < c_dofs1D; ++dx)
1289 {
1290 for (int j=0; j<2; ++j)
1291 {
1292 w[j][dx][ey] += Bo(ex, dx) * xd * vk(j, local_index, e);
1293 }
1294 }
1295 }
1296 }
1297
1298 // contract in y
1299 for (int dx = 0; dx < c_dofs1D; ++dx)
1300 {
1301 for (int dy = 0; dy < c_dofs1D; ++dy)
1302 {
1303 for (int j=0; j<2; ++j)
1304 {
1305 real_t s = 0.0;
1306 for (int ey = 0; ey < c_dofs1D; ++ey)
1307 {
1308 s += w[j][dx][ey] * Bc(ey, dy);
1309 }
1310 y(dx, dy, j, e) += s;
1311 }
1312 }
1313 }
1314
1315 // dofs that point parallel to y-axis (open in y, closed in x)
1316
1317 // contract in x
1318 for (int ey = 0; ey < o_dofs1D; ++ey)
1319 {
1320 for (int dx = 0; dx < c_dofs1D; ++dx)
1321 {
1322 for (int j=0; j<2; ++j) { w[j][dx][ey] = 0.0; }
1323 }
1324 for (int ex = 0; ex < c_dofs1D; ++ex)
1325 {
1326 const int local_index = c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
1327 const real_t xd = x(local_index, e);
1328 for (int dx = 0; dx < c_dofs1D; ++dx)
1329 {
1330 for (int j=0; j<2; ++j)
1331 {
1332 w[j][dx][ey] += Bc(ex, dx) * xd * vk(j, local_index, e);
1333 }
1334 }
1335 }
1336 }
1337
1338 // contract in y
1339 for (int dx = 0; dx < c_dofs1D; ++dx)
1340 {
1341 for (int dy = 0; dy < c_dofs1D; ++dy)
1342 {
1343 for (int j=0; j<2; ++j)
1344 {
1345 real_t s = 0.0;
1346 for (int ey = 0; ey < o_dofs1D; ++ey)
1347 {
1348 s += w[j][dx][ey] * Bo(ey, dy);
1349 }
1350 y(dx, dy, j, e) += s;
1351 }
1352 }
1353 }
1354 });
1355}
1356
1357static void PAHcurlVecH1IdentityApply3D(const int c_dofs1D,
1358 const int o_dofs1D,
1359 const int NE,
1360 const Array<real_t> &Bclosed,
1361 const Array<real_t> &Bopen,
1362 const Vector &pa_data,
1363 const Vector &x_,
1364 Vector &y_)
1365{
1366 auto Bc = Reshape(Bclosed.Read(), c_dofs1D, c_dofs1D);
1367 auto Bo = Reshape(Bopen.Read(), o_dofs1D, c_dofs1D);
1368
1369 auto x = Reshape(x_.Read(), c_dofs1D, c_dofs1D, c_dofs1D, 3, NE);
1370 auto y = Reshape(y_.ReadWrite(), (3 * c_dofs1D * c_dofs1D * o_dofs1D), NE);
1371
1372 auto vk = Reshape(pa_data.Read(), 3, (3 * c_dofs1D * c_dofs1D * o_dofs1D),
1373 NE);
1374 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().MAX_D1D &&
1375 o_dofs1D <= c_dofs1D, "");
1376
1377 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1378 {
1379 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
1380
1381 real_t w1[3][MAX_D1D][MAX_D1D][MAX_D1D];
1382 real_t w2[3][MAX_D1D][MAX_D1D][MAX_D1D];
1383
1384 // dofs that point parallel to x-axis (open in x, closed in y, z)
1385
1386 // contract in z
1387 for (int ez = 0; ez < c_dofs1D; ++ez)
1388 {
1389 for (int dx = 0; dx < c_dofs1D; ++dx)
1390 {
1391 for (int dy = 0; dy < c_dofs1D; ++dy)
1392 {
1393 for (int j=0; j<3; ++j)
1394 {
1395 w1[j][dx][dy][ez] = 0.0;
1396 for (int dz = 0; dz < c_dofs1D; ++dz)
1397 {
1398 w1[j][dx][dy][ez] += Bc(ez, dz) * x(dx, dy, dz, j, e);
1399 }
1400 }
1401 }
1402 }
1403 }
1404
1405 // contract in y
1406 for (int ez = 0; ez < c_dofs1D; ++ez)
1407 {
1408 for (int ey = 0; ey < c_dofs1D; ++ey)
1409 {
1410 for (int dx = 0; dx < c_dofs1D; ++dx)
1411 {
1412 for (int j=0; j<3; ++j)
1413 {
1414 w2[j][dx][ey][ez] = 0.0;
1415 for (int dy = 0; dy < c_dofs1D; ++dy)
1416 {
1417 w2[j][dx][ey][ez] += Bc(ey, dy) * w1[j][dx][dy][ez];
1418 }
1419 }
1420 }
1421 }
1422 }
1423
1424 // contract in x
1425 for (int ez = 0; ez < c_dofs1D; ++ez)
1426 {
1427 for (int ey = 0; ey < c_dofs1D; ++ey)
1428 {
1429 for (int ex = 0; ex < o_dofs1D; ++ex)
1430 {
1431 for (int j=0; j<3; ++j)
1432 {
1433 real_t s = 0.0;
1434 for (int dx = 0; dx < c_dofs1D; ++dx)
1435 {
1436 s += Bo(ex, dx) * w2[j][dx][ey][ez];
1437 }
1438 const int local_index = ez*c_dofs1D*o_dofs1D + ey*o_dofs1D + ex;
1439 y(local_index, e) += s * vk(j, local_index, e);
1440 }
1441 }
1442 }
1443 }
1444
1445 // dofs that point parallel to y-axis (open in y, closed in x, z)
1446
1447 // contract in z
1448 for (int ez = 0; ez < c_dofs1D; ++ez)
1449 {
1450 for (int dx = 0; dx < c_dofs1D; ++dx)
1451 {
1452 for (int dy = 0; dy < c_dofs1D; ++dy)
1453 {
1454 for (int j=0; j<3; ++j)
1455 {
1456 w1[j][dx][dy][ez] = 0.0;
1457 for (int dz = 0; dz < c_dofs1D; ++dz)
1458 {
1459 w1[j][dx][dy][ez] += Bc(ez, dz) * x(dx, dy, dz, j, e);
1460 }
1461 }
1462 }
1463 }
1464 }
1465
1466 // contract in y
1467 for (int ez = 0; ez < c_dofs1D; ++ez)
1468 {
1469 for (int ey = 0; ey < o_dofs1D; ++ey)
1470 {
1471 for (int dx = 0; dx < c_dofs1D; ++dx)
1472 {
1473 for (int j=0; j<3; ++j)
1474 {
1475 w2[j][dx][ey][ez] = 0.0;
1476 for (int dy = 0; dy < c_dofs1D; ++dy)
1477 {
1478 w2[j][dx][ey][ez] += Bo(ey, dy) * w1[j][dx][dy][ez];
1479 }
1480 }
1481 }
1482 }
1483 }
1484
1485 // contract in x
1486 for (int ez = 0; ez < c_dofs1D; ++ez)
1487 {
1488 for (int ey = 0; ey < o_dofs1D; ++ey)
1489 {
1490 for (int ex = 0; ex < c_dofs1D; ++ex)
1491 {
1492 for (int j=0; j<3; ++j)
1493 {
1494 real_t s = 0.0;
1495 for (int dx = 0; dx < c_dofs1D; ++dx)
1496 {
1497 s += Bc(ex, dx) * w2[j][dx][ey][ez];
1498 }
1499 const int local_index = c_dofs1D*c_dofs1D*o_dofs1D +
1500 ez*c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
1501 y(local_index, e) += s * vk(j, local_index, e);
1502 }
1503 }
1504 }
1505 }
1506
1507 // dofs that point parallel to z-axis (open in z, closed in x, y)
1508
1509 // contract in z
1510 for (int ez = 0; ez < o_dofs1D; ++ez)
1511 {
1512 for (int dx = 0; dx < c_dofs1D; ++dx)
1513 {
1514 for (int dy = 0; dy < c_dofs1D; ++dy)
1515 {
1516 for (int j=0; j<3; ++j)
1517 {
1518 w1[j][dx][dy][ez] = 0.0;
1519 for (int dz = 0; dz < c_dofs1D; ++dz)
1520 {
1521 w1[j][dx][dy][ez] += Bo(ez, dz) * x(dx, dy, dz, j, e);
1522 }
1523 }
1524 }
1525 }
1526 }
1527
1528 // contract in y
1529 for (int ez = 0; ez < o_dofs1D; ++ez)
1530 {
1531 for (int ey = 0; ey < c_dofs1D; ++ey)
1532 {
1533 for (int dx = 0; dx < c_dofs1D; ++dx)
1534 {
1535 for (int j=0; j<3; ++j)
1536 {
1537 w2[j][dx][ey][ez] = 0.0;
1538 for (int dy = 0; dy < c_dofs1D; ++dy)
1539 {
1540 w2[j][dx][ey][ez] += Bc(ey, dy) * w1[j][dx][dy][ez];
1541 }
1542 }
1543 }
1544 }
1545 }
1546
1547 // contract in x
1548 for (int ez = 0; ez < o_dofs1D; ++ez)
1549 {
1550 for (int ey = 0; ey < c_dofs1D; ++ey)
1551 {
1552 for (int ex = 0; ex < c_dofs1D; ++ex)
1553 {
1554 for (int j=0; j<3; ++j)
1555 {
1556 real_t s = 0.0;
1557 for (int dx = 0; dx < c_dofs1D; ++dx)
1558 {
1559 s += Bc(ex, dx) * w2[j][dx][ey][ez];
1560 }
1561 const int local_index = 2*c_dofs1D*c_dofs1D*o_dofs1D +
1562 ez*c_dofs1D*c_dofs1D + ey*c_dofs1D + ex;
1563 y(local_index, e) += s * vk(j, local_index, e);
1564 }
1565 }
1566 }
1567 }
1568 });
1569}
1570
1571static void PAHcurlVecH1IdentityApplyTranspose3D(const int c_dofs1D,
1572 const int o_dofs1D,
1573 const int NE,
1574 const Array<real_t> &Bclosed,
1575 const Array<real_t> &Bopen,
1576 const Vector &pa_data,
1577 const Vector &x_,
1578 Vector &y_)
1579{
1580 auto Bc = Reshape(Bclosed.Read(), c_dofs1D, c_dofs1D);
1581 auto Bo = Reshape(Bopen.Read(), o_dofs1D, c_dofs1D);
1582
1583 auto x = Reshape(x_.Read(), (3 * c_dofs1D * c_dofs1D * o_dofs1D), NE);
1584 auto y = Reshape(y_.ReadWrite(), c_dofs1D, c_dofs1D, c_dofs1D, 3, NE);
1585
1586 auto vk = Reshape(pa_data.Read(), 3, (3 * c_dofs1D * c_dofs1D * o_dofs1D),
1587 NE);
1588
1589 MFEM_VERIFY(c_dofs1D <= DeviceDofQuadLimits::Get().MAX_D1D &&
1590 o_dofs1D <= c_dofs1D, "");
1591
1592 mfem::forall(NE, [=] MFEM_HOST_DEVICE (int e)
1593 {
1594 constexpr static int MAX_D1D = DofQuadLimits::HCURL_MAX_D1D;
1595
1596 real_t w1[3][MAX_D1D][MAX_D1D][MAX_D1D];
1597 real_t w2[3][MAX_D1D][MAX_D1D][MAX_D1D];
1598
1599 // dofs that point parallel to x-axis (open in x, closed in y, z)
1600
1601 // contract in x
1602 for (int ez = 0; ez < c_dofs1D; ++ez)
1603 {
1604 for (int ey = 0; ey < c_dofs1D; ++ey)
1605 {
1606 for (int j=0; j<3; ++j)
1607 {
1608 for (int dx = 0; dx < c_dofs1D; ++dx)
1609 {
1610 w2[j][dx][ey][ez] = 0.0;
1611 }
1612 for (int ex = 0; ex < o_dofs1D; ++ex)
1613 {
1614 const int local_index = ez*c_dofs1D*o_dofs1D + ey*o_dofs1D + ex;
1615 const real_t xv = x(local_index, e) * vk(j, local_index, e);
1616 for (int dx = 0; dx < c_dofs1D; ++dx)
1617 {
1618 w2[j][dx][ey][ez] += xv * Bo(ex, dx);
1619 }
1620 }
1621 }
1622 }
1623 }
1624
1625 // contract in y
1626 for (int ez = 0; ez < c_dofs1D; ++ez)
1627 {
1628 for (int dx = 0; dx < c_dofs1D; ++dx)
1629 {
1630 for (int dy = 0; dy < c_dofs1D; ++dy)
1631 {
1632 for (int j=0; j<3; ++j)
1633 {
1634 w1[j][dx][dy][ez] = 0.0;
1635 for (int ey = 0; ey < c_dofs1D; ++ey)
1636 {
1637 w1[j][dx][dy][ez] += w2[j][dx][ey][ez] * Bc(ey, dy);
1638 }
1639 }
1640 }
1641 }
1642 }
1643
1644 // contract in z
1645 for (int dx = 0; dx < c_dofs1D; ++dx)
1646 {
1647 for (int dy = 0; dy < c_dofs1D; ++dy)
1648 {
1649 for (int dz = 0; dz < c_dofs1D; ++dz)
1650 {
1651 for (int j=0; j<3; ++j)
1652 {
1653 real_t s = 0.0;
1654 for (int ez = 0; ez < c_dofs1D; ++ez)
1655 {
1656 s += w1[j][dx][dy][ez] * Bc(ez, dz);
1657 }
1658 y(dx, dy, dz, j, e) += s;
1659 }
1660 }
1661 }
1662 }
1663
1664 // dofs that point parallel to y-axis (open in y, closed in x, z)
1665
1666 // contract in x
1667 for (int ez = 0; ez < c_dofs1D; ++ez)
1668 {
1669 for (int ey = 0; ey < o_dofs1D; ++ey)
1670 {
1671 for (int j=0; j<3; ++j)
1672 {
1673 for (int dx = 0; dx < c_dofs1D; ++dx)
1674 {
1675 w2[j][dx][ey][ez] = 0.0;
1676 }
1677 for (int ex = 0; ex < c_dofs1D; ++ex)
1678 {
1679 const int local_index = c_dofs1D*c_dofs1D*o_dofs1D +
1680 ez*c_dofs1D*o_dofs1D + ey*c_dofs1D + ex;
1681 const real_t xv = x(local_index, e) * vk(j, local_index, e);
1682 for (int dx = 0; dx < c_dofs1D; ++dx)
1683 {
1684 w2[j][dx][ey][ez] += xv * Bc(ex, dx);
1685 }
1686 }
1687 }
1688 }
1689 }
1690
1691 // contract in y
1692 for (int ez = 0; ez < c_dofs1D; ++ez)
1693 {
1694 for (int dx = 0; dx < c_dofs1D; ++dx)
1695 {
1696 for (int dy = 0; dy < c_dofs1D; ++dy)
1697 {
1698 for (int j=0; j<3; ++j)
1699 {
1700 w1[j][dx][dy][ez] = 0.0;
1701 for (int ey = 0; ey < o_dofs1D; ++ey)
1702 {
1703 w1[j][dx][dy][ez] += w2[j][dx][ey][ez] * Bo(ey, dy);
1704 }
1705 }
1706 }
1707 }
1708 }
1709
1710 // contract in z
1711 for (int dx = 0; dx < c_dofs1D; ++dx)
1712 {
1713 for (int dy = 0; dy < c_dofs1D; ++dy)
1714 {
1715 for (int dz = 0; dz < c_dofs1D; ++dz)
1716 {
1717 for (int j=0; j<3; ++j)
1718 {
1719 real_t s = 0.0;
1720 for (int ez = 0; ez < c_dofs1D; ++ez)
1721 {
1722 s += w1[j][dx][dy][ez] * Bc(ez, dz);
1723 }
1724 y(dx, dy, dz, j, e) += s;
1725 }
1726 }
1727 }
1728 }
1729
1730 // dofs that point parallel to z-axis (open in z, closed in x, y)
1731
1732 // contract in x
1733 for (int ez = 0; ez < o_dofs1D; ++ez)
1734 {
1735 for (int ey = 0; ey < c_dofs1D; ++ey)
1736 {
1737 for (int j=0; j<3; ++j)
1738 {
1739 for (int dx = 0; dx < c_dofs1D; ++dx)
1740 {
1741 w2[j][dx][ey][ez] = 0.0;
1742 }
1743 for (int ex = 0; ex < c_dofs1D; ++ex)
1744 {
1745 const int local_index = 2*c_dofs1D*c_dofs1D*o_dofs1D +
1746 ez*c_dofs1D*c_dofs1D + ey*c_dofs1D + ex;
1747 const real_t xv = x(local_index, e) * vk(j, local_index, e);
1748 for (int dx = 0; dx < c_dofs1D; ++dx)
1749 {
1750 w2[j][dx][ey][ez] += xv * Bc(ex, dx);
1751 }
1752 }
1753 }
1754 }
1755 }
1756
1757 // contract in y
1758 for (int ez = 0; ez < o_dofs1D; ++ez)
1759 {
1760 for (int dx = 0; dx < c_dofs1D; ++dx)
1761 {
1762 for (int dy = 0; dy < c_dofs1D; ++dy)
1763 {
1764 for (int j=0; j<3; ++j)
1765 {
1766 w1[j][dx][dy][ez] = 0.0;
1767 for (int ey = 0; ey < c_dofs1D; ++ey)
1768 {
1769 w1[j][dx][dy][ez] += w2[j][dx][ey][ez] * Bc(ey, dy);
1770 }
1771 }
1772 }
1773 }
1774 }
1775
1776 // contract in z
1777 for (int dx = 0; dx < c_dofs1D; ++dx)
1778 {
1779 for (int dy = 0; dy < c_dofs1D; ++dy)
1780 {
1781 for (int dz = 0; dz < c_dofs1D; ++dz)
1782 {
1783 for (int j=0; j<3; ++j)
1784 {
1785 real_t s = 0.0;
1786 for (int ez = 0; ez < o_dofs1D; ++ez)
1787 {
1788 s += w1[j][dx][dy][ez] * Bo(ez, dz);
1789 }
1790 y(dx, dy, dz, j, e) += s;
1791 }
1792 }
1793 }
1794 }
1795 });
1796}
1797
1799 const FiniteElementSpace &test_fes)
1800{
1801 // Assumes tensor-product elements, with a vector test space and H^1 trial space.
1802 Mesh *mesh = trial_fes.GetMesh();
1803 const FiniteElement *trial_fel = trial_fes.GetTypicalFE();
1804 const FiniteElement *test_fel = test_fes.GetTypicalFE();
1805
1806 const NodalTensorFiniteElement *trial_el =
1807 dynamic_cast<const NodalTensorFiniteElement*>(trial_fel);
1808 MFEM_VERIFY(trial_el != NULL, "Only NodalTensorFiniteElement is supported!");
1809
1810 const VectorTensorFiniteElement *test_el =
1811 dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
1812 MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!");
1813
1814 const int dims = trial_el->GetDim();
1815 MFEM_VERIFY(dims == 2 || dims == 3, "");
1816
1817 dim = mesh->Dimension();
1818 MFEM_VERIFY(dim == 2 || dim == 3, "");
1819
1820 MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
1821
1822 MFEM_VERIFY(vdim == 1, "vdim != 1 with PA is not supported yet!");
1823
1824 ne = trial_fes.GetNE();
1825
1826 const int order = trial_el->GetOrder();
1827 dofquad_fe.reset(new H1_SegmentElement(order));
1829 mfem::IntegrationRule closed_ir;
1830 closed_ir.SetSize(order + 1);
1831 qf1d.GaussLobatto(order + 1, &closed_ir);
1832 mfem::IntegrationRule open_ir;
1833 open_ir.SetSize(order);
1834 qf1d.GaussLegendre(order, &open_ir);
1835
1836 maps_C_C = &dofquad_fe->GetDofToQuad(closed_ir, DofToQuad::TENSOR);
1837 maps_O_C = &dofquad_fe->GetDofToQuad(open_ir, DofToQuad::TENSOR);
1838
1839 o_dofs1D = maps_O_C->nqpt;
1840 c_dofs1D = maps_C_C->nqpt;
1841 MFEM_VERIFY(maps_O_C->ndof == c_dofs1D &&
1842 maps_C_C->ndof == c_dofs1D, "Discrepancy in the number of DOFs");
1843
1844 const int ndof_test = (dim == 3) ? 3 * c_dofs1D * c_dofs1D * o_dofs1D
1845 : 2 * c_dofs1D * o_dofs1D;
1846
1847 const IntegrationRule & Nodes = test_el->GetNodes();
1848
1849 pa_data.SetSize(dim * ndof_test * ne, Device::GetMemoryType());
1850 auto op = Reshape(pa_data.HostWrite(), dim, ndof_test, ne);
1851
1852 const Array<int> &dofmap = test_el->GetDofMap();
1853
1854 if (dim == 3)
1855 {
1856 // Note that ND_HexahedronElement uses 6 vectors in tk rather than 3, with
1857 // the last 3 having negative signs. Here the signs are all positive, as
1858 // signs are applied in ElementRestriction.
1859
1860 const real_t tk[9] = { 1.,0.,0., 0.,1.,0., 0.,0.,1. };
1861
1862 for (int c=0; c<3; ++c)
1863 {
1864 for (int i=0; i<ndof_test/3; ++i)
1865 {
1866 const int d = (c*ndof_test/3) + i;
1867 // ND_HexahedronElement sets dof2tk = (dofmap < 0) ? 3+c : c, but here
1868 // no signs should be applied due to ElementRestriction.
1869 const int dof2tk = c;
1870 const int id = (dofmap[d] >= 0) ? dofmap[d] : -1 - dofmap[d];
1871
1872 for (int e=0; e<ne; ++e)
1873 {
1874 real_t v[3];
1876 tr->SetIntPoint(&Nodes.IntPoint(id));
1877 tr->Jacobian().Mult(tk + dof2tk*dim, v);
1878
1879 for (int j=0; j<3; ++j)
1880 {
1881 op(j,d,e) = v[j];
1882 }
1883 }
1884 }
1885 }
1886 }
1887 else // 2D case
1888 {
1889 const real_t tk[4] = { 1.,0., 0.,1. };
1890 for (int c=0; c<2; ++c)
1891 {
1892 for (int i=0; i<ndof_test/2; ++i)
1893 {
1894 const int d = (c*ndof_test/2) + i;
1895 // ND_QuadrilateralElement sets dof2tk = (dofmap < 0) ? 2+c : c, but here
1896 // no signs should be applied due to ElementRestriction.
1897 const int dof2tk = c;
1898 const int id = (dofmap[d] >= 0) ? dofmap[d] : -1 - dofmap[d];
1899
1900 for (int e=0; e<ne; ++e)
1901 {
1902 real_t v[2];
1904 tr->SetIntPoint(&Nodes.IntPoint(id));
1905 tr->Jacobian().Mult(tk + dof2tk*dim, v);
1906
1907 for (int j=0; j<2; ++j)
1908 {
1909 op(j,d,e) = v[j];
1910 }
1911 }
1912 }
1913 }
1914 }
1915}
1916
1918{
1919 if (dim == 3)
1920 {
1921 PAHcurlVecH1IdentityApply3D(c_dofs1D, o_dofs1D, ne, maps_C_C->B, maps_O_C->B,
1922 pa_data, x, y);
1923 }
1924 else if (dim == 2)
1925 {
1926 PAHcurlVecH1IdentityApply2D(c_dofs1D, o_dofs1D, ne, maps_C_C->B, maps_O_C->B,
1927 pa_data, x, y);
1928 }
1929 else
1930 {
1931 mfem_error("Bad dimension!");
1932 }
1933}
1934
1936{
1937 if (dim == 3)
1938 {
1939 PAHcurlVecH1IdentityApplyTranspose3D(c_dofs1D, o_dofs1D, ne, maps_C_C->B,
1940 maps_O_C->B, pa_data, x, y);
1941 }
1942 else if (dim == 2)
1943 {
1944 PAHcurlVecH1IdentityApplyTranspose2D(c_dofs1D, o_dofs1D, ne, maps_C_C->B,
1945 maps_O_C->B, pa_data, x, y);
1946 }
1947 else
1948 {
1949 mfem_error("Bad dimension!");
1950 }
1951}
1952
1953} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:337
@ GaussLobatto
Closed type.
Definition fe_base.hpp:36
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition device.hpp:278
Array< real_t > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition fe_base.hpp:214
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:165
Array< real_t > B
Basis functions evaluated at quadrature points.
Definition fe_base.hpp:193
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:178
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:182
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:891
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
const FiniteElement * GetTypicalFE() const
Return GetFE(0) if the local mesh is not empty; otherwise return a typical FE based on the Geometry t...
Definition fespace.cpp:3871
Abstract class for all finite elements.
Definition fe_base.hpp:244
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition fe_base.cpp:365
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:338
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:321
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:400
void AddMultPA(const Vector &x, Vector &y) const override
Method for partially assembled action.
void AddMultTransposePA(const Vector &x, Vector &y) const override
Method for partially assembled transposed action.
void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
Setup method for PA data.
Arbitrary order H1 elements in 1D.
Definition fe_h1.hpp:23
void AddMultTransposePA(const Vector &x, Vector &y) const override
Method for partially assembled transposed action.
void AddMultPA(const Vector &x, Vector &y) const override
Method for partially assembled action.
void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
Mesh data type.
Definition mesh.hpp:64
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:357
A Class that defines 1-D numerical quadrature rules on [0,1].
Definition intrules.hpp:376
static void GaussLegendre(const int np, IntegrationRule *ir)
Definition intrules.cpp:424
static void GaussLobatto(const int np, IntegrationRule *ir)
Definition intrules.cpp:512
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition fe_base.hpp:1273
Vector data type.
Definition vector.hpp:82
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:494
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:510
void mfem_error(const char *msg)
Definition error.cpp:154
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