MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_rt.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
12// Raviart-Thomas Finite Element classes
13
14#include "fe_rt.hpp"
15#include "face_map_utils.hpp"
16#include "../coefficient.hpp"
17
18namespace mfem
19{
20
21using namespace std;
22
23const real_t RT_QuadrilateralElement::nk[8] =
24{ 0., -1., 1., 0., 0., 1., -1., 0. };
25
27 const int cb_type,
28 const int ob_type)
29 : VectorTensorFiniteElement(2, 2*(p + 1)*(p + 2), p + 1, cb_type, ob_type,
30 H_DIV, DofMapType::L2_DOF_MAP),
31 dof2nk(dof),
32 cp(poly1d.ClosedPoints(p + 1, cb_type))
33{
34 if (obasis1d.IsIntegratedType()) { is_nodal = false; }
35
37
38 const real_t *op = poly1d.OpenPoints(p, ob_type);
39 const int dof2 = dof/2;
40
41#ifndef MFEM_THREAD_SAFE
42 shape_cx.SetSize(p + 2);
43 shape_ox.SetSize(p + 1);
44 shape_cy.SetSize(p + 2);
45 shape_oy.SetSize(p + 1);
46 dshape_cx.SetSize(p + 2);
47 dshape_cy.SetSize(p + 2);
48#endif
49
50 // edges
51 int o = 0;
52 for (int i = 0; i <= p; i++) // (0,1)
53 {
54 dof_map[1*dof2 + i + 0*(p + 1)] = o++;
55 }
56 for (int i = 0; i <= p; i++) // (1,2)
57 {
58 dof_map[0*dof2 + (p + 1) + i*(p + 2)] = o++;
59 }
60 for (int i = 0; i <= p; i++) // (2,3)
61 {
62 dof_map[1*dof2 + (p - i) + (p + 1)*(p + 1)] = o++;
63 }
64 for (int i = 0; i <= p; i++) // (3,0)
65 {
66 dof_map[0*dof2 + 0 + (p - i)*(p + 2)] = o++;
67 }
68
69 // interior
70 for (int j = 0; j <= p; j++) // x-components
71 for (int i = 1; i <= p; i++)
72 {
73 dof_map[0*dof2 + i + j*(p + 2)] = o++;
74 }
75 for (int j = 1; j <= p; j++) // y-components
76 for (int i = 0; i <= p; i++)
77 {
78 dof_map[1*dof2 + i + j*(p + 1)] = o++;
79 }
80
81 // dof orientations
82 // x-components
83 for (int j = 0; j <= p; j++)
84 for (int i = 0; i <= p/2; i++)
85 {
86 int idx = 0*dof2 + i + j*(p + 2);
87 dof_map[idx] = -1 - dof_map[idx];
88 }
89 if (p%2 == 1)
90 for (int j = p/2 + 1; j <= p; j++)
91 {
92 int idx = 0*dof2 + (p/2 + 1) + j*(p + 2);
93 dof_map[idx] = -1 - dof_map[idx];
94 }
95 // y-components
96 for (int j = 0; j <= p/2; j++)
97 for (int i = 0; i <= p; i++)
98 {
99 int idx = 1*dof2 + i + j*(p + 1);
100 dof_map[idx] = -1 - dof_map[idx];
101 }
102 if (p%2 == 1)
103 for (int i = 0; i <= p/2; i++)
104 {
105 int idx = 1*dof2 + i + (p/2 + 1)*(p + 1);
106 dof_map[idx] = -1 - dof_map[idx];
107 }
108
109 o = 0;
110 for (int j = 0; j <= p; j++)
111 for (int i = 0; i <= p + 1; i++)
112 {
113 int idx;
114 if ((idx = dof_map[o++]) < 0)
115 {
116 idx = -1 - idx;
117 dof2nk[idx] = 3;
118 }
119 else
120 {
121 dof2nk[idx] = 1;
122 }
123 Nodes.IntPoint(idx).Set2(cp[i], op[j]);
124 }
125 for (int j = 0; j <= p + 1; j++)
126 for (int i = 0; i <= p; i++)
127 {
128 int idx;
129 if ((idx = dof_map[o++]) < 0)
130 {
131 idx = -1 - idx;
132 dof2nk[idx] = 0;
133 }
134 else
135 {
136 dof2nk[idx] = 2;
137 }
138 Nodes.IntPoint(idx).Set2(op[i], cp[j]);
139 }
140}
141
143 DenseMatrix &shape) const
144{
145 const int pp1 = order;
146
147#ifdef MFEM_THREAD_SAFE
148 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
149#endif
150
152 {
153#ifdef MFEM_THREAD_SAFE
154 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
155#endif
156 basis1d.Eval(ip.x, shape_cx, dshape_cx);
157 basis1d.Eval(ip.y, shape_cy, dshape_cy);
159 obasis1d.EvalIntegrated(dshape_cx, shape_ox);
160 obasis1d.EvalIntegrated(dshape_cy, shape_oy);
161 }
162 else
163 {
164 basis1d.Eval(ip.x, shape_cx);
165 basis1d.Eval(ip.y, shape_cy);
166 obasis1d.Eval(ip.x, shape_ox);
167 obasis1d.Eval(ip.y, shape_oy);
168 }
169
170 int o = 0;
171 for (int j = 0; j < pp1; j++)
172 for (int i = 0; i <= pp1; i++)
173 {
174 int idx, s;
175 if ((idx = dof_map[o++]) < 0)
176 {
177 idx = -1 - idx, s = -1;
178 }
179 else
180 {
181 s = +1;
182 }
183 shape(idx,0) = s*shape_cx(i)*shape_oy(j);
184 shape(idx,1) = 0.;
185 }
186 for (int j = 0; j <= pp1; j++)
187 for (int i = 0; i < pp1; i++)
188 {
189 int idx, s;
190 if ((idx = dof_map[o++]) < 0)
191 {
192 idx = -1 - idx, s = -1;
193 }
194 else
195 {
196 s = +1;
197 }
198 shape(idx,0) = 0.;
199 shape(idx,1) = s*shape_ox(i)*shape_cy(j);
200 }
201}
202
204 Vector &divshape) const
205{
206 const int pp1 = order;
207
208#ifdef MFEM_THREAD_SAFE
209 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
210 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
211#endif
212
213 basis1d.Eval(ip.x, shape_cx, dshape_cx);
214 basis1d.Eval(ip.y, shape_cy, dshape_cy);
216 {
218 obasis1d.EvalIntegrated(dshape_cx, shape_ox);
219 obasis1d.EvalIntegrated(dshape_cy, shape_oy);
220 }
221 else
222 {
223 obasis1d.Eval(ip.x, shape_ox);
224 obasis1d.Eval(ip.y, shape_oy);
225 }
226
227 int o = 0;
228 for (int j = 0; j < pp1; j++)
229 for (int i = 0; i <= pp1; i++)
230 {
231 int idx, s;
232 if ((idx = dof_map[o++]) < 0)
233 {
234 idx = -1 - idx, s = -1;
235 }
236 else
237 {
238 s = +1;
239 }
240 divshape(idx) = s*dshape_cx(i)*shape_oy(j);
241 }
242 for (int j = 0; j <= pp1; j++)
243 for (int i = 0; i < pp1; i++)
244 {
245 int idx, s;
246 if ((idx = dof_map[o++]) < 0)
247 {
248 idx = -1 - idx, s = -1;
249 }
250 else
251 {
252 s = +1;
253 }
254 divshape(idx) = s*shape_ox(i)*dshape_cy(j);
255 }
256}
257
260 Vector &dofs) const
261{
262 MFEM_ASSERT(obasis1d.IsIntegratedType(), "Not integrated type");
264 Vector xk(vk, vc.GetVDim());
265
267 const int nqpt = ir.GetNPoints();
268
269 IntegrationPoint ip2d;
270
271 int o = 0;
272 for (int c = 0; c < 2; c++)
273 {
274 int im = (c == 0) ? order + 1 : order;
275 int jm = (c == 1) ? order + 1 : order;
276 for (int j = 0; j < jm; j++)
277 for (int i = 0; i < im; i++)
278 {
279 int idx = dof_map[o++];
280 if (idx < 0) { idx = -1 - idx; }
281 int ic = (c == 0) ? j : i;
282 const real_t h = cp[ic+1] - cp[ic];
283 real_t val = 0.0;
284 for (int k = 0; k < nqpt; k++)
285 {
286 const IntegrationPoint &ip1d = ir.IntPoint(k);
287 if (c == 0) { ip2d.Set2(cp[i], cp[j] + (h*ip1d.x)); }
288 else { ip2d.Set2(cp[i] + (h*ip1d.x), cp[j]); }
289 Trans.SetIntPoint(&ip2d);
290 vc.Eval(xk, Trans, ip2d);
291 // nk^t adj(J) xk
292 const real_t ipval = Trans.AdjugateJacobian().InnerProduct(vk,
293 nk + dof2nk[idx]*dim);
294 val += ip1d.weight*ipval;
295 }
296 dofs(idx) = val*h;
297 }
298 }
299}
300
302 Array<int> &face_map) const
303{
304 const int p = order;
305 const int pp1 = p + 1;
306 const int n_face_dofs = p;
307
308 std::vector<int> offsets;
309 std::vector<int> strides = {(face_id == 0 || face_id == 2) ? 1 : pp1};
310 switch (face_id)
311 {
312 case 0: offsets = {p*pp1}; break; // y = 0
313 case 1: offsets = {pp1 - 1}; break; // x = 1
314 case 2: offsets = {p*pp1 + p*(pp1 - 1)}; break; // y = 1
315 case 3: offsets = {0}; break; // x = 0
316 }
317
318 std::vector<int> n_dofs(dim - 1, p);
319 internal::FillFaceMap(n_face_dofs, offsets, strides, n_dofs, face_map);
320}
321
322
323const real_t RT_HexahedronElement::nk[18] =
324{ 0.,0.,-1., 0.,-1.,0., 1.,0.,0., 0.,1.,0., -1.,0.,0., 0.,0.,1. };
325
327 const int cb_type,
328 const int ob_type)
329 : VectorTensorFiniteElement(3, 3*(p + 1)*(p + 1)*(p + 2), p + 1, cb_type,
330 ob_type, H_DIV, DofMapType::L2_DOF_MAP),
331 dof2nk(dof),
332 cp(poly1d.ClosedPoints(p + 1, cb_type))
333{
334 if (obasis1d.IsIntegratedType()) { is_nodal = false; }
335
337
338 const real_t *op = poly1d.OpenPoints(p, ob_type);
339 const int dof3 = dof/3;
340
341#ifndef MFEM_THREAD_SAFE
342 shape_cx.SetSize(p + 2);
343 shape_ox.SetSize(p + 1);
344 shape_cy.SetSize(p + 2);
345 shape_oy.SetSize(p + 1);
346 shape_cz.SetSize(p + 2);
347 shape_oz.SetSize(p + 1);
348 dshape_cx.SetSize(p + 2);
349 dshape_cy.SetSize(p + 2);
350 dshape_cz.SetSize(p + 2);
351#endif
352
353 // faces
354 int o = 0;
355 for (int j = 0; j <= p; j++) // (3,2,1,0) -- bottom
356 for (int i = 0; i <= p; i++)
357 {
358 dof_map[2*dof3 + i + ((p - j) + 0*(p + 1))*(p + 1)] = o++;
359 }
360 for (int j = 0; j <= p; j++) // (0,1,5,4) -- front
361 for (int i = 0; i <= p; i++)
362 {
363 dof_map[1*dof3 + i + (0 + j*(p + 2))*(p + 1)] = o++;
364 }
365 for (int j = 0; j <= p; j++) // (1,2,6,5) -- right
366 for (int i = 0; i <= p; i++)
367 {
368 dof_map[0*dof3 + (p + 1) + (i + j*(p + 1))*(p + 2)] = o++;
369 }
370 for (int j = 0; j <= p; j++) // (2,3,7,6) -- back
371 for (int i = 0; i <= p; i++)
372 {
373 dof_map[1*dof3 + (p - i) + ((p + 1) + j*(p + 2))*(p + 1)] = o++;
374 }
375 for (int j = 0; j <= p; j++) // (3,0,4,7) -- left
376 for (int i = 0; i <= p; i++)
377 {
378 dof_map[0*dof3 + 0 + ((p - i) + j*(p + 1))*(p + 2)] = o++;
379 }
380 for (int j = 0; j <= p; j++) // (4,5,6,7) -- top
381 for (int i = 0; i <= p; i++)
382 {
383 dof_map[2*dof3 + i + (j + (p + 1)*(p + 1))*(p + 1)] = o++;
384 }
385
386 // interior
387 // x-components
388 for (int k = 0; k <= p; k++)
389 for (int j = 0; j <= p; j++)
390 for (int i = 1; i <= p; i++)
391 {
392 dof_map[0*dof3 + i + (j + k*(p + 1))*(p + 2)] = o++;
393 }
394 // y-components
395 for (int k = 0; k <= p; k++)
396 for (int j = 1; j <= p; j++)
397 for (int i = 0; i <= p; i++)
398 {
399 dof_map[1*dof3 + i + (j + k*(p + 2))*(p + 1)] = o++;
400 }
401 // z-components
402 for (int k = 1; k <= p; k++)
403 for (int j = 0; j <= p; j++)
404 for (int i = 0; i <= p; i++)
405 {
406 dof_map[2*dof3 + i + (j + k*(p + 1))*(p + 1)] = o++;
407 }
408
409 // dof orientations
410 // for odd p, do not change the orientations in the mid-planes
411 // {i = p/2 + 1}, {j = p/2 + 1}, {k = p/2 + 1} in the x, y, z-components
412 // respectively.
413 // x-components
414 for (int k = 0; k <= p; k++)
415 for (int j = 0; j <= p; j++)
416 for (int i = 0; i <= p/2; i++)
417 {
418 int idx = 0*dof3 + i + (j + k*(p + 1))*(p + 2);
419 dof_map[idx] = -1 - dof_map[idx];
420 }
421 // y-components
422 for (int k = 0; k <= p; k++)
423 for (int j = 0; j <= p/2; j++)
424 for (int i = 0; i <= p; i++)
425 {
426 int idx = 1*dof3 + i + (j + k*(p + 2))*(p + 1);
427 dof_map[idx] = -1 - dof_map[idx];
428 }
429 // z-components
430 for (int k = 0; k <= p/2; k++)
431 for (int j = 0; j <= p; j++)
432 for (int i = 0; i <= p; i++)
433 {
434 int idx = 2*dof3 + i + (j + k*(p + 1))*(p + 1);
435 dof_map[idx] = -1 - dof_map[idx];
436 }
437
438 o = 0;
439 // x-components
440 for (int k = 0; k <= p; k++)
441 for (int j = 0; j <= p; j++)
442 for (int i = 0; i <= p + 1; i++)
443 {
444 int idx;
445 if ((idx = dof_map[o++]) < 0)
446 {
447 idx = -1 - idx;
448 dof2nk[idx] = 4;
449 }
450 else
451 {
452 dof2nk[idx] = 2;
453 }
454 Nodes.IntPoint(idx).Set3(cp[i], op[j], op[k]);
455 }
456 // y-components
457 for (int k = 0; k <= p; k++)
458 for (int j = 0; j <= p + 1; j++)
459 for (int i = 0; i <= p; i++)
460 {
461 int idx;
462 if ((idx = dof_map[o++]) < 0)
463 {
464 idx = -1 - idx;
465 dof2nk[idx] = 1;
466 }
467 else
468 {
469 dof2nk[idx] = 3;
470 }
471 Nodes.IntPoint(idx).Set3(op[i], cp[j], op[k]);
472 }
473 // z-components
474 for (int k = 0; k <= p + 1; k++)
475 for (int j = 0; j <= p; j++)
476 for (int i = 0; i <= p; i++)
477 {
478 int idx;
479 if ((idx = dof_map[o++]) < 0)
480 {
481 idx = -1 - idx;
482 dof2nk[idx] = 0;
483 }
484 else
485 {
486 dof2nk[idx] = 5;
487 }
488 Nodes.IntPoint(idx).Set3(op[i], op[j], cp[k]);
489 }
490}
491
493 DenseMatrix &shape) const
494{
495 const int pp1 = order;
496
497#ifdef MFEM_THREAD_SAFE
498 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
499 Vector shape_cz(pp1 + 1), shape_oz(pp1);
500#endif
501
503 {
504#ifdef MFEM_THREAD_SAFE
505 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
506#endif
507 basis1d.Eval(ip.x, shape_cx, dshape_cx);
508 basis1d.Eval(ip.y, shape_cy, dshape_cy);
509 basis1d.Eval(ip.z, shape_cz, dshape_cz);
511 obasis1d.EvalIntegrated(dshape_cx, shape_ox);
512 obasis1d.EvalIntegrated(dshape_cy, shape_oy);
513 obasis1d.EvalIntegrated(dshape_cz, shape_oz);
514 }
515 else
516 {
517 basis1d.Eval(ip.x, shape_cx);
518 basis1d.Eval(ip.y, shape_cy);
519 basis1d.Eval(ip.z, shape_cz);
520 obasis1d.Eval(ip.x, shape_ox);
521 obasis1d.Eval(ip.y, shape_oy);
522 obasis1d.Eval(ip.z, shape_oz);
523 }
524
525 int o = 0;
526 // x-components
527 for (int k = 0; k < pp1; k++)
528 for (int j = 0; j < pp1; j++)
529 for (int i = 0; i <= pp1; i++)
530 {
531 int idx, s;
532 if ((idx = dof_map[o++]) < 0)
533 {
534 idx = -1 - idx, s = -1;
535 }
536 else
537 {
538 s = +1;
539 }
540 shape(idx,0) = s*shape_cx(i)*shape_oy(j)*shape_oz(k);
541 shape(idx,1) = 0.;
542 shape(idx,2) = 0.;
543 }
544 // y-components
545 for (int k = 0; k < pp1; k++)
546 for (int j = 0; j <= pp1; j++)
547 for (int i = 0; i < pp1; i++)
548 {
549 int idx, s;
550 if ((idx = dof_map[o++]) < 0)
551 {
552 idx = -1 - idx, s = -1;
553 }
554 else
555 {
556 s = +1;
557 }
558 shape(idx,0) = 0.;
559 shape(idx,1) = s*shape_ox(i)*shape_cy(j)*shape_oz(k);
560 shape(idx,2) = 0.;
561 }
562 // z-components
563 for (int k = 0; k <= pp1; k++)
564 for (int j = 0; j < pp1; j++)
565 for (int i = 0; i < pp1; i++)
566 {
567 int idx, s;
568 if ((idx = dof_map[o++]) < 0)
569 {
570 idx = -1 - idx, s = -1;
571 }
572 else
573 {
574 s = +1;
575 }
576 shape(idx,0) = 0.;
577 shape(idx,1) = 0.;
578 shape(idx,2) = s*shape_ox(i)*shape_oy(j)*shape_cz(k);
579 }
580}
581
583 Vector &divshape) const
584{
585 const int pp1 = order;
586
587#ifdef MFEM_THREAD_SAFE
588 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
589 Vector shape_cz(pp1 + 1), shape_oz(pp1);
590 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1), dshape_cz(pp1 + 1);
591#endif
592
593 basis1d.Eval(ip.x, shape_cx, dshape_cx);
594 basis1d.Eval(ip.y, shape_cy, dshape_cy);
595 basis1d.Eval(ip.z, shape_cz, dshape_cz);
597 {
599 obasis1d.EvalIntegrated(dshape_cx, shape_ox);
600 obasis1d.EvalIntegrated(dshape_cy, shape_oy);
601 obasis1d.EvalIntegrated(dshape_cz, shape_oz);
602 }
603 else
604 {
605 obasis1d.Eval(ip.x, shape_ox);
606 obasis1d.Eval(ip.y, shape_oy);
607 obasis1d.Eval(ip.z, shape_oz);
608 }
609
610 int o = 0;
611 // x-components
612 for (int k = 0; k < pp1; k++)
613 for (int j = 0; j < pp1; j++)
614 for (int i = 0; i <= pp1; i++)
615 {
616 int idx, s;
617 if ((idx = dof_map[o++]) < 0)
618 {
619 idx = -1 - idx, s = -1;
620 }
621 else
622 {
623 s = +1;
624 }
625 divshape(idx) = s*dshape_cx(i)*shape_oy(j)*shape_oz(k);
626 }
627 // y-components
628 for (int k = 0; k < pp1; k++)
629 for (int j = 0; j <= pp1; j++)
630 for (int i = 0; i < pp1; i++)
631 {
632 int idx, s;
633 if ((idx = dof_map[o++]) < 0)
634 {
635 idx = -1 - idx, s = -1;
636 }
637 else
638 {
639 s = +1;
640 }
641 divshape(idx) = s*shape_ox(i)*dshape_cy(j)*shape_oz(k);
642 }
643 // z-components
644 for (int k = 0; k <= pp1; k++)
645 for (int j = 0; j < pp1; j++)
646 for (int i = 0; i < pp1; i++)
647 {
648 int idx, s;
649 if ((idx = dof_map[o++]) < 0)
650 {
651 idx = -1 - idx, s = -1;
652 }
653 else
654 {
655 s = +1;
656 }
657 divshape(idx) = s*shape_ox(i)*shape_oy(j)*dshape_cz(k);
658 }
659}
660
663 Vector &dofs) const
664{
665 MFEM_ASSERT(obasis1d.IsIntegratedType(), "Not integrated type");
667 Vector xq(vq, vc.GetVDim());
668
670 const int nqpt = ir2d.GetNPoints();
671
672 IntegrationPoint ip3d;
673
674 int o = 0;
675 for (int c = 0; c < 3; c++)
676 {
677 int im = (c == 0) ? order + 1 : order;
678 int jm = (c == 1) ? order + 1 : order;
679 int km = (c == 2) ? order + 1 : order;
680 for (int k = 0; k < km; k++)
681 for (int j = 0; j < jm; j++)
682 for (int i = 0; i < im; i++)
683 {
684 int idx = dof_map[o++];
685 if (idx < 0) { idx = -1 - idx; }
686 int ic1, ic2;
687 if (c == 0) { ic1 = j; ic2 = k; }
688 else if (c == 1) { ic1 = i; ic2 = k; }
689 else { ic1 = i; ic2 = j; }
690 const real_t h1 = cp[ic1+1] - cp[ic1];
691 const real_t h2 = cp[ic2+1] - cp[ic2];
692 real_t val = 0.0;
693 for (int q = 0; q < nqpt; q++)
694 {
695 const IntegrationPoint &ip2d = ir2d.IntPoint(q);
696 if (c == 0) { ip3d.Set3(cp[i], cp[j] + h1*ip2d.x, cp[k] + h2*ip2d.y); }
697 else if (c == 1) { ip3d.Set3(cp[i] + h1*ip2d.x, cp[j], cp[k] + h2*ip2d.y); }
698 else { ip3d.Set3(cp[i] + h1*ip2d.x, cp[j] + h2*ip2d.y, cp[k]); }
699 Trans.SetIntPoint(&ip3d);
700 vc.Eval(xq, Trans, ip3d);
701 // nk^t adj(J) xq
702 const real_t ipval
703 = Trans.AdjugateJacobian().InnerProduct(vq, nk + dof2nk[idx]*dim);
704 val += ip2d.weight*ipval;
705 }
706 dofs(idx) = val*h1*h2;
707 }
708 }
709}
710
711void RT_HexahedronElement::GetFaceMap(const int face_id,
712 Array<int> &face_map) const
713{
714 const int p = order;
715 const int pp1 = p + 1;
716 int n_face_dofs = p*p;
717 std::vector<int> strides, offsets;
718 const int n_dof_per_dim = p*p*pp1;
719 const auto f = internal::GetFaceNormal3D(face_id);
720 const int face_normal = f.first, level = f.second;
721 if (face_normal == 0) // x-normal
722 {
723 offsets = {level ? pp1 - 1 : 0};
724 strides = {pp1, p*pp1};
725 }
726 else if (face_normal == 1) // y-normal
727 {
728 offsets = {n_dof_per_dim + (level ? p*(pp1 - 1) : 0)};
729 strides = {1, p*pp1};
730 }
731 else if (face_normal == 2) // z-normal
732 {
733 offsets = {2*n_dof_per_dim + (level ? p*p*(pp1 - 1) : 0)};
734 strides = {1, p};
735 }
736 std::vector<int> n_dofs = {p, p};
737 internal::FillFaceMap(n_face_dofs, offsets, strides, n_dofs, face_map);
738}
739
740
741const real_t RT_TriangleElement::nk[6] =
742{ 0., -1., 1., 1., -1., 0. };
743
744const real_t RT_TriangleElement::c = 1./3.;
745
747 : VectorFiniteElement(2, Geometry::TRIANGLE, (p + 1)*(p + 3), p + 1,
748 H_DIV, FunctionSpace::Pk),
749 dof2nk(dof)
750{
751 const real_t *iop = (p > 0) ? poly1d.OpenPoints(p - 1) : NULL;
752 const real_t *bop = poly1d.OpenPoints(p);
753
754#ifndef MFEM_THREAD_SAFE
755 shape_x.SetSize(p + 1);
756 shape_y.SetSize(p + 1);
757 shape_l.SetSize(p + 1);
758 dshape_x.SetSize(p + 1);
759 dshape_y.SetSize(p + 1);
760 dshape_l.SetSize(p + 1);
761 u.SetSize(dof, dim);
762 divu.SetSize(dof);
763#else
764 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
765#endif
766
767 // edges
768 int o = 0;
769 for (int i = 0; i <= p; i++) // (0,1)
770 {
771 Nodes.IntPoint(o).Set2(bop[i], 0.);
772 dof2nk[o++] = 0;
773 }
774 for (int i = 0; i <= p; i++) // (1,2)
775 {
776 Nodes.IntPoint(o).Set2(bop[p-i], bop[i]);
777 dof2nk[o++] = 1;
778 }
779 for (int i = 0; i <= p; i++) // (2,0)
780 {
781 Nodes.IntPoint(o).Set2(0., bop[p-i]);
782 dof2nk[o++] = 2;
783 }
784
785 // interior
786 for (int j = 0; j < p; j++)
787 for (int i = 0; i + j < p; i++)
788 {
789 real_t w = iop[i] + iop[j] + iop[p-1-i-j];
790 Nodes.IntPoint(o).Set2(iop[i]/w, iop[j]/w);
791 dof2nk[o++] = 0;
792 Nodes.IntPoint(o).Set2(iop[i]/w, iop[j]/w);
793 dof2nk[o++] = 2;
794 }
795
796 DenseMatrix T(dof);
797 for (int k = 0; k < dof; k++)
798 {
799 const IntegrationPoint &ip = Nodes.IntPoint(k);
800 poly1d.CalcBasis(p, ip.x, shape_x);
801 poly1d.CalcBasis(p, ip.y, shape_y);
802 poly1d.CalcBasis(p, 1. - ip.x - ip.y, shape_l);
803 const real_t *n_k = nk + 2*dof2nk[k];
804
805 o = 0;
806 for (int j = 0; j <= p; j++)
807 for (int i = 0; i + j <= p; i++)
808 {
809 real_t s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
810 T(o++, k) = s*n_k[0];
811 T(o++, k) = s*n_k[1];
812 }
813 for (int i = 0; i <= p; i++)
814 {
815 real_t s = shape_x(i)*shape_y(p-i);
816 T(o++, k) = s*((ip.x - c)*n_k[0] + (ip.y - c)*n_k[1]);
817 }
818 }
819
820 Ti.Factor(T);
821 // mfem::out << "RT_TriangleElement(" << p << ") : "; Ti.TestInversion();
822}
823
825 DenseMatrix &shape) const
826{
827 const int p = order - 1;
828
829#ifdef MFEM_THREAD_SAFE
830 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
831 DenseMatrix u(dof, dim);
832#endif
833
834 poly1d.CalcBasis(p, ip.x, shape_x);
835 poly1d.CalcBasis(p, ip.y, shape_y);
836 poly1d.CalcBasis(p, 1. - ip.x - ip.y, shape_l);
837
838 int o = 0;
839 for (int j = 0; j <= p; j++)
840 for (int i = 0; i + j <= p; i++)
841 {
842 real_t s = shape_x(i)*shape_y(j)*shape_l(p-i-j);
843 u(o,0) = s; u(o,1) = 0; o++;
844 u(o,0) = 0; u(o,1) = s; o++;
845 }
846 for (int i = 0; i <= p; i++)
847 {
848 real_t s = shape_x(i)*shape_y(p-i);
849 u(o,0) = (ip.x - c)*s;
850 u(o,1) = (ip.y - c)*s;
851 o++;
852 }
853
854 Ti.Mult(u, shape);
855}
856
858 Vector &divshape) const
859{
860 const int p = order - 1;
861
862#ifdef MFEM_THREAD_SAFE
863 Vector shape_x(p + 1), shape_y(p + 1), shape_l(p + 1);
864 Vector dshape_x(p + 1), dshape_y(p + 1), dshape_l(p + 1);
865 Vector divu(dof);
866#endif
867
868 poly1d.CalcBasis(p, ip.x, shape_x, dshape_x);
869 poly1d.CalcBasis(p, ip.y, shape_y, dshape_y);
870 poly1d.CalcBasis(p, 1. - ip.x - ip.y, shape_l, dshape_l);
871
872 int o = 0;
873 for (int j = 0; j <= p; j++)
874 for (int i = 0; i + j <= p; i++)
875 {
876 int k = p - i - j;
877 divu(o++) = (dshape_x(i)*shape_l(k) -
878 shape_x(i)*dshape_l(k))*shape_y(j);
879 divu(o++) = (dshape_y(j)*shape_l(k) -
880 shape_y(j)*dshape_l(k))*shape_x(i);
881 }
882 for (int i = 0; i <= p; i++)
883 {
884 int j = p - i;
885 divu(o++) = ((shape_x(i) + (ip.x - c)*dshape_x(i))*shape_y(j) +
886 (shape_y(j) + (ip.y - c)*dshape_y(j))*shape_x(i));
887 }
888
889 Ti.Mult(divu, divshape);
890}
891
892
893const real_t RT_TetrahedronElement::nk[12] =
894{ 1,1,1, -1,0,0, 0,-1,0, 0,0,-1 };
895// { .5,.5,.5, -.5,0,0, 0,-.5,0, 0,0,-.5}; // n_F |F|
896
897const real_t RT_TetrahedronElement::c = 1./4.;
898
900 : VectorFiniteElement(3, Geometry::TETRAHEDRON, (p + 1)*(p + 2)*(p + 4)/2,
901 p + 1, H_DIV, FunctionSpace::Pk),
902 dof2nk(dof)
903{
904 const real_t *iop = (p > 0) ? poly1d.OpenPoints(p - 1) : NULL;
905 const real_t *bop = poly1d.OpenPoints(p);
906
907#ifndef MFEM_THREAD_SAFE
908 shape_x.SetSize(p + 1);
909 shape_y.SetSize(p + 1);
910 shape_z.SetSize(p + 1);
911 shape_l.SetSize(p + 1);
912 dshape_x.SetSize(p + 1);
913 dshape_y.SetSize(p + 1);
914 dshape_z.SetSize(p + 1);
915 dshape_l.SetSize(p + 1);
916 u.SetSize(dof, dim);
917 divu.SetSize(dof);
918#else
919 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
920#endif
921
922 int o = 0;
923 // faces (see Mesh::GenerateFaces in mesh/mesh.cpp,
924 // the constructor of H1_TetrahedronElement)
925 for (int j = 0; j <= p; j++)
926 for (int i = 0; i + j <= p; i++) // (1,2,3)
927 {
928 real_t w = bop[i] + bop[j] + bop[p-i-j];
929 Nodes.IntPoint(o).Set3(bop[p-i-j]/w, bop[i]/w, bop[j]/w);
930 dof2nk[o++] = 0;
931 }
932 for (int j = 0; j <= p; j++)
933 for (int i = 0; i + j <= p; i++) // (0,3,2)
934 {
935 real_t w = bop[i] + bop[j] + bop[p-i-j];
936 Nodes.IntPoint(o).Set3(0., bop[j]/w, bop[i]/w);
937 dof2nk[o++] = 1;
938 }
939 for (int j = 0; j <= p; j++)
940 for (int i = 0; i + j <= p; i++) // (0,1,3)
941 {
942 real_t w = bop[i] + bop[j] + bop[p-i-j];
943 Nodes.IntPoint(o).Set3(bop[i]/w, 0., bop[j]/w);
944 dof2nk[o++] = 2;
945 }
946 for (int j = 0; j <= p; j++)
947 for (int i = 0; i + j <= p; i++) // (0,2,1)
948 {
949 real_t w = bop[i] + bop[j] + bop[p-i-j];
950 Nodes.IntPoint(o).Set3(bop[j]/w, bop[i]/w, 0.);
951 dof2nk[o++] = 3;
952 }
953
954 // interior
955 for (int k = 0; k < p; k++)
956 for (int j = 0; j + k < p; j++)
957 for (int i = 0; i + j + k < p; i++)
958 {
959 real_t w = iop[i] + iop[j] + iop[k] + iop[p-1-i-j-k];
960 Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
961 dof2nk[o++] = 1;
962 Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
963 dof2nk[o++] = 2;
964 Nodes.IntPoint(o).Set3(iop[i]/w, iop[j]/w, iop[k]/w);
965 dof2nk[o++] = 3;
966 }
967
968 DenseMatrix T(dof);
969 for (int m = 0; m < dof; m++)
970 {
971 const IntegrationPoint &ip = Nodes.IntPoint(m);
972 poly1d.CalcBasis(p, ip.x, shape_x);
973 poly1d.CalcBasis(p, ip.y, shape_y);
974 poly1d.CalcBasis(p, ip.z, shape_z);
975 poly1d.CalcBasis(p, 1. - ip.x - ip.y - ip.z, shape_l);
976 const real_t *nm = nk + 3*dof2nk[m];
977
978 o = 0;
979 for (int k = 0; k <= p; k++)
980 for (int j = 0; j + k <= p; j++)
981 for (int i = 0; i + j + k <= p; i++)
982 {
983 real_t s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
984 T(o++, m) = s * nm[0];
985 T(o++, m) = s * nm[1];
986 T(o++, m) = s * nm[2];
987 }
988 for (int j = 0; j <= p; j++)
989 for (int i = 0; i + j <= p; i++)
990 {
991 real_t s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
992 T(o++, m) = s*((ip.x - c)*nm[0] + (ip.y - c)*nm[1] +
993 (ip.z - c)*nm[2]);
994 }
995 }
996
997 Ti.Factor(T);
998 // mfem::out << "RT_TetrahedronElement(" << p << ") : "; Ti.TestInversion();
999}
1000
1002 DenseMatrix &shape) const
1003{
1004 const int p = order - 1;
1005
1006#ifdef MFEM_THREAD_SAFE
1007 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
1008 DenseMatrix u(dof, dim);
1009#endif
1010
1011 poly1d.CalcBasis(p, ip.x, shape_x);
1012 poly1d.CalcBasis(p, ip.y, shape_y);
1013 poly1d.CalcBasis(p, ip.z, shape_z);
1014 poly1d.CalcBasis(p, 1. - ip.x - ip.y - ip.z, shape_l);
1015
1016 int o = 0;
1017 for (int k = 0; k <= p; k++)
1018 for (int j = 0; j + k <= p; j++)
1019 for (int i = 0; i + j + k <= p; i++)
1020 {
1021 real_t s = shape_x(i)*shape_y(j)*shape_z(k)*shape_l(p-i-j-k);
1022 u(o,0) = s; u(o,1) = 0; u(o,2) = 0; o++;
1023 u(o,0) = 0; u(o,1) = s; u(o,2) = 0; o++;
1024 u(o,0) = 0; u(o,1) = 0; u(o,2) = s; o++;
1025 }
1026 for (int j = 0; j <= p; j++)
1027 for (int i = 0; i + j <= p; i++)
1028 {
1029 real_t s = shape_x(i)*shape_y(j)*shape_z(p-i-j);
1030 u(o,0) = (ip.x - c)*s; u(o,1) = (ip.y - c)*s; u(o,2) = (ip.z - c)*s;
1031 o++;
1032 }
1033
1034 Ti.Mult(u, shape);
1035}
1036
1038 Vector &divshape) const
1039{
1040 const int p = order - 1;
1041
1042#ifdef MFEM_THREAD_SAFE
1043 Vector shape_x(p + 1), shape_y(p + 1), shape_z(p + 1), shape_l(p + 1);
1044 Vector dshape_x(p + 1), dshape_y(p + 1), dshape_z(p + 1), dshape_l(p + 1);
1045 Vector divu(dof);
1046#endif
1047
1048 poly1d.CalcBasis(p, ip.x, shape_x, dshape_x);
1049 poly1d.CalcBasis(p, ip.y, shape_y, dshape_y);
1050 poly1d.CalcBasis(p, ip.z, shape_z, dshape_z);
1051 poly1d.CalcBasis(p, 1. - ip.x - ip.y - ip.z, shape_l, dshape_l);
1052
1053 int o = 0;
1054 for (int k = 0; k <= p; k++)
1055 for (int j = 0; j + k <= p; j++)
1056 for (int i = 0; i + j + k <= p; i++)
1057 {
1058 int l = p - i - j - k;
1059 divu(o++) = (dshape_x(i)*shape_l(l) -
1060 shape_x(i)*dshape_l(l))*shape_y(j)*shape_z(k);
1061 divu(o++) = (dshape_y(j)*shape_l(l) -
1062 shape_y(j)*dshape_l(l))*shape_x(i)*shape_z(k);
1063 divu(o++) = (dshape_z(k)*shape_l(l) -
1064 shape_z(k)*dshape_l(l))*shape_x(i)*shape_y(j);
1065 }
1066 for (int j = 0; j <= p; j++)
1067 for (int i = 0; i + j <= p; i++)
1068 {
1069 int k = p - i - j;
1070 divu(o++) =
1071 (shape_x(i) + (ip.x - c)*dshape_x(i))*shape_y(j)*shape_z(k) +
1072 (shape_y(j) + (ip.y - c)*dshape_y(j))*shape_x(i)*shape_z(k) +
1073 (shape_z(k) + (ip.z - c)*dshape_z(k))*shape_x(i)*shape_y(j);
1074 }
1075
1076 Ti.Mult(divu, divshape);
1077}
1078
1079const real_t RT_WedgeElement::nk[15] =
1080{ 0,0,-1, 0,0,1, 0,-1,0, 1,1,0, -1,0,0};
1081
1083 : VectorFiniteElement(3, Geometry::PRISM,
1084 (p + 2) * ((p + 1) * (p + 2)) / 2 +
1085 (p + 1) * (p + 1) * (p + 3), p + 1,
1086 H_DIV, FunctionSpace::Qk),
1087 dof2nk(dof),
1088 t_dof(dof),
1089 s_dof(dof),
1090 L2TriangleFE(p),
1091 RTTriangleFE(p),
1092 H1SegmentFE(p + 1),
1093 L2SegmentFE(p)
1094{
1095 MFEM_ASSERT(L2TriangleFE.GetDof() * H1SegmentFE.GetDof() +
1096 RTTriangleFE.GetDof() * L2SegmentFE.GetDof() == dof,
1097 "Mismatch in number of degrees of freedom "
1098 "when building RT_WedgeElement!");
1099
1100 const int pm1 = p - 1;
1101
1102#ifndef MFEM_THREAD_SAFE
1103 tl2_shape.SetSize(L2TriangleFE.GetDof());
1104 sh1_shape.SetSize(H1SegmentFE.GetDof());
1105 trt_shape.SetSize(RTTriangleFE.GetDof(), 2);
1106 sl2_shape.SetSize(L2SegmentFE.GetDof());
1107 sh1_dshape.SetSize(H1SegmentFE.GetDof(), 1);
1108 trt_dshape.SetSize(RTTriangleFE.GetDof());
1109#endif
1110
1111 const IntegrationRule &tl2_n = L2TriangleFE.GetNodes();
1112 const IntegrationRule &trt_n = RTTriangleFE.GetNodes();
1113 const IntegrationRule &sh1_n = H1SegmentFE.GetNodes();
1114 const IntegrationRule &sl2_n = L2SegmentFE.GetNodes();
1115
1116 // faces
1117 int o = 0;
1118 int l = 0;
1119 // (0,2,1) -- bottom
1120 for (int j = 0; j <= p; j++)
1121 for (int i = 0; i + j <= p; i++)
1122 {
1123 l = j + i * (2 * p + 3 - i) / 2;
1124 t_dof[o] = l; s_dof[o] = 0; dof2nk[o] = 0;
1125 const IntegrationPoint & t_ip = tl2_n.IntPoint(t_dof[o]);
1126 Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sh1_n.IntPoint(s_dof[o]).x);
1127 o++;
1128 }
1129 // (3,4,5) -- top
1130 l = 0;
1131 for (int j = 0; j <= p; j++)
1132 for (int i = 0; i + j <= p; i++)
1133 {
1134 t_dof[o] = l; s_dof[o] = 1; dof2nk[o] = 1; l++;
1135 const IntegrationPoint & t_ip = tl2_n.IntPoint(t_dof[o]);
1136 Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sh1_n.IntPoint(s_dof[o]).x);
1137 o++;
1138 }
1139 // (0, 1, 4, 3) -- xz plane
1140 for (int j = 0; j <= p; j++)
1141 for (int i = 0; i <= p; i++)
1142 {
1143 t_dof[o] = i; s_dof[o] = j; dof2nk[o] = 2;
1144 const IntegrationPoint & t_ip = trt_n.IntPoint(t_dof[o]);
1145 Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sl2_n.IntPoint(s_dof[o]).x);
1146 o++;
1147 }
1148 // (1, 2, 5, 4) -- (y-x)z plane
1149 for (int j = 0; j <= p; j++)
1150 for (int i = 0; i <= p; i++)
1151 {
1152 t_dof[o] = p + 1 + i; s_dof[o] = j; dof2nk[o] = 3;
1153 const IntegrationPoint & t_ip = trt_n.IntPoint(t_dof[o]);
1154 Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sl2_n.IntPoint(s_dof[o]).x);
1155 o++;
1156 }
1157 // (2, 0, 3, 5) -- yz plane
1158 for (int j = 0; j <= p; j++)
1159 for (int i = 0; i <= p; i++)
1160 {
1161 t_dof[o] = 2 * p + 2 + i; s_dof[o] = j; dof2nk[o] = 4;
1162 const IntegrationPoint & t_ip = trt_n.IntPoint(t_dof[o]);
1163 Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sl2_n.IntPoint(s_dof[o]).x);
1164 o++;
1165 }
1166
1167 // interior
1168 for (int k = 0; k < L2SegmentFE.GetDof(); k++)
1169 {
1170 l = 0;
1171 for (int j = 0; j <= pm1; j++)
1172 for (int i = 0; i + j <= pm1; i++)
1173 {
1174 t_dof[o] = 3 * (p + 1) + 2 * l; s_dof[o] = k;
1175 dof2nk[o] = 2;
1176 const IntegrationPoint & t_ip0 = trt_n.IntPoint(t_dof[o]);
1177 const IntegrationPoint & s_ip0 = sl2_n.IntPoint(s_dof[o]);
1178 Nodes.IntPoint(o).Set3(t_ip0.x, t_ip0.y, s_ip0.x);
1179 o++;
1180 t_dof[o] = 3 * (p + 1) + 2 * l + 1; s_dof[o] = k;
1181 dof2nk[o] = 4; l++;
1182 const IntegrationPoint & t_ip1 = trt_n.IntPoint(t_dof[o]);
1183 const IntegrationPoint & s_ip1 = sl2_n.IntPoint(s_dof[o]);
1184 Nodes.IntPoint(o).Set3(t_ip1.x, t_ip1.y, s_ip1.x);
1185 o++;
1186 }
1187 }
1188 for (int k = 2; k < H1SegmentFE.GetDof(); k++)
1189 {
1190 for (l = 0; l < L2TriangleFE.GetDof(); l++)
1191 {
1192 t_dof[o] = l; s_dof[o] = k; dof2nk[o] = 1;
1193 const IntegrationPoint & t_ip = tl2_n.IntPoint(t_dof[o]);
1194 Nodes.IntPoint(o).Set3(t_ip.x, t_ip.y, sh1_n.IntPoint(s_dof[o]).x);
1195 o++;
1196 }
1197 }
1198}
1199
1201 DenseMatrix &shape) const
1202{
1203#ifdef MFEM_THREAD_SAFE
1204 DenseMatrix trt_shape(RTTriangleFE.GetDof(), 2);
1205 Vector tl2_shape(L2TriangleFE.GetDof());
1206 Vector sh1_shape(H1SegmentFE.GetDof());
1207 Vector sl2_shape(L2SegmentFE.GetDof());
1208#endif
1209
1210 IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1211
1212 L2TriangleFE.CalcShape(ip, tl2_shape);
1213 RTTriangleFE.CalcVShape(ip, trt_shape);
1214 H1SegmentFE.CalcShape(ipz, sh1_shape);
1215 L2SegmentFE.CalcShape(ipz, sl2_shape);
1216
1217 for (int i=0; i<dof; i++)
1218 {
1219 if ( dof2nk[i] >= 2 )
1220 {
1221 shape(i, 0) = trt_shape(t_dof[i], 0) * sl2_shape[s_dof[i]];
1222 shape(i, 1) = trt_shape(t_dof[i], 1) * sl2_shape[s_dof[i]];
1223 shape(i, 2) = 0.0;
1224 }
1225 else
1226 {
1227 real_t s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1228 shape(i, 0) = 0.0;
1229 shape(i, 1) = 0.0;
1230 shape(i, 2) = s * tl2_shape[t_dof[i]] * sh1_shape(s_dof[i]);
1231 }
1232 }
1233}
1234
1236 Vector &divshape) const
1237{
1238#ifdef MFEM_THREAD_SAFE
1239 Vector trt_dshape(RTTriangleFE.GetDof());
1240 Vector tl2_shape(L2TriangleFE.GetDof());
1241 Vector sl2_shape(L2SegmentFE.GetDof());
1242 DenseMatrix sh1_dshape(H1SegmentFE.GetDof(), 1);
1243#endif
1244
1245 IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1246
1247 RTTriangleFE.CalcDivShape(ip, trt_dshape);
1248 L2TriangleFE.CalcShape(ip, tl2_shape);
1249
1250 L2SegmentFE.CalcShape(ipz, sl2_shape);
1251 H1SegmentFE.CalcDShape(ipz, sh1_dshape);
1252
1253 for (int i=0; i<dof; i++)
1254 {
1255 if ( dof2nk[i] >= 2 )
1256 {
1257 divshape(i) = trt_dshape(t_dof[i]) * sl2_shape(s_dof[i]);
1258 }
1259 else
1260 {
1261 real_t s = (dof2nk[i] == 0) ? -1.0 : 1.0;
1262 divshape(i) = s * tl2_shape(t_dof[i]) * sh1_dshape(s_dof[i], 0);
1263 }
1264 }
1265}
1266
1267const real_t RT_FuentesPyramidElement::nk[24] =
1268{
1269 0,0,-1, 0,-1,0, 1,0,1, 0,1,1, -1,0,0,
1270 M_SQRT2,0,M_SQRT1_2, 0,M_SQRT2,M_SQRT1_2, 0,0,1
1271};
1272
1274 : VectorFiniteElement(3, Geometry::PYRAMID, (p + 1)*(3*p*(p + 2) + 5),
1275 p + 1, H_DIV, FunctionSpace::Uk),
1276 dof2nk(dof)
1277{
1278 zmax = 0.0;
1279
1280 const real_t *iop = poly1d.OpenPoints(p);
1281 const real_t *icp = poly1d.ClosedPoints(p + 1);
1282 const real_t *bop = poly1d.OpenPoints(p);
1283
1284#ifndef MFEM_THREAD_SAFE
1285 tmp1_i.SetSize(p + 2);
1286 tmp1_ij.SetSize(p + 2, p + 2);
1287 tmp2_ij.SetSize(p + 2, dim);
1288 tmp3_ij.SetSize(p + 2, dim);
1289 tmp4_ij.SetSize(p + 1, p + 1);
1290 tmp1_ijk.SetSize(p + 1, p + 1, dim);
1291 tmp2_ijk.SetSize(p + 1, p + 1, dim);
1292 tmp3_ijk.SetSize(p + 1, p + 1, dim);
1293 tmp4_ijk.SetSize(p + 1, p + 2, dim);
1294 tmp5_ijk.SetSize(p + 1, p + 2, dim);
1295 tmp6_ijk.SetSize(p + 2, p + 2, dim);
1296 tmp7_ijk.SetSize(p + 2, p + 2, dim);
1297 u.SetSize(dof, dim);
1298 divu.SetSize(dof);
1299#else
1300 Vector tmp1_i(p + 2);
1301 DenseMatrix tmp1_ij(p + 2, p + 2);
1302 DenseMatrix tmp2_ij(p + 2, dim);
1303 DenseMatrix tmp3_ij(p + 2, dim);
1304 DenseTensor tmp1_ijk(p + 1, p + 1, dim);
1305 DenseTensor tmp2_ijk(p + 1, p + 1, dim);
1306 DenseTensor tmp3_ijk(p + 1, p + 1, dim);
1307 DenseTensor tmp4_ijk(p + 1, p + 2, dim);
1308 DenseTensor tmp5_ijk(p + 1, p + 2, dim);
1309 DenseTensor tmp6_ijk(p + 2, p + 2, dim);
1310 DenseTensor tmp7_ijk(p + 2, p + 2, dim);
1311 DenseMatrix u(dof, dim);
1312#endif
1313
1314 int o = 0;
1315
1316 // quadrilateral face
1317 for (int j = 0; j <= p; j++)
1318 for (int i = 0; i <= p; i++) // (3,2,1,0)
1319 {
1320 Nodes.IntPoint(o).Set3(bop[i], bop[p-j], 0.);
1321 dof2nk[o++] = 0;
1322 }
1323 // triangular faces
1324 for (int j = 0; j <= p; j++)
1325 for (int i = 0; i + j <= p; i++) // (0,1,4)
1326 {
1327 real_t w = bop[i] + bop[j] + bop[p-i-j];
1328 Nodes.IntPoint(o).Set3(bop[i]/w, 0., bop[j]/w);
1329 dof2nk[o++] = 1;
1330 }
1331 for (int j = 0; j <= p; j++)
1332 for (int i = 0; i + j <= p; i++) // (1,2,4)
1333 {
1334 real_t w = bop[i] + bop[j] + bop[p-i-j];
1335 Nodes.IntPoint(o).Set3(1.-bop[j]/w, bop[i]/w, bop[j]/w);
1336 dof2nk[o++] = 2;
1337 }
1338 for (int j = 0; j <= p; j++)
1339 for (int i = p - j; i >= 0; i--) // (2,3,4)
1340 {
1341 real_t w = bop[i] + bop[j] + bop[p-i-j];
1342 Nodes.IntPoint(o).Set3(bop[i]/w, 1.0-bop[j]/w, bop[j]/w);
1343 dof2nk[o++] = 3;
1344 }
1345 for (int j = 0; j <= p; j++)
1346 for (int i = p - j; i >= 0; i--) // (3,0,4)
1347 {
1348 real_t w = bop[i] + bop[j] + bop[p-i-j];
1349 Nodes.IntPoint(o).Set3(0., bop[i]/w, bop[j]/w);
1350 dof2nk[o++] = 4;
1351 }
1352
1353 // interior
1354 // x-components
1355 for (int k = 0; k <= p; k++)
1356 for (int j = 0; j <= p; j++)
1357 for (int i = 1; i <= p; i++)
1358 {
1359 real_t w = 1.0 - iop[k];
1360 Nodes.IntPoint(o).Set3(icp[i]*w, iop[j]*w, iop[k]);
1361 dof2nk[o++] = 5;
1362 }
1363 // y-components
1364 for (int k = 0; k <= p; k++)
1365 for (int j = 1; j <= p; j++)
1366 for (int i = 0; i <= p; i++)
1367 {
1368 real_t w = 1.0 - iop[k];
1369 Nodes.IntPoint(o).Set3(iop[i]*w, icp[j]*w, iop[k]);
1370 dof2nk[o++] = 6;
1371 }
1372 // z-components
1373 for (int k = 1; k <= p; k++)
1374 for (int j = 0; j <= p; j++)
1375 for (int i = 0; i <= p; i++)
1376 {
1377 real_t w = 1.0 - icp[k];
1378 Nodes.IntPoint(o).Set3(iop[i]*w, iop[j]*w, icp[k]);
1379 dof2nk[o++] = 7;
1380 }
1381
1382 DenseMatrix T(dof);
1383
1384 for (int m = 0; m < dof; m++)
1385 {
1386 const IntegrationPoint &ip = Nodes.IntPoint(m);
1387 const Vector nm({nk[3*dof2nk[m]], nk[3*dof2nk[m]+1], nk[3*dof2nk[m]+2]});
1388 calcBasis(order, ip, tmp1_i, tmp1_ij, tmp2_ij,
1389 tmp1_ijk, tmp2_ijk, tmp3_ijk, tmp4_ijk, tmp5_ijk, tmp6_ijk,
1390 tmp7_ijk,
1391 tmp3_ij, u);
1392 u.Mult(nm, T.GetColumn(m));
1393 }
1394
1395 Ti.Factor(T);
1396}
1397
1399 DenseMatrix &shape) const
1400{
1401#ifdef MFEM_THREAD_SAFE
1402 const int p = order - 1;
1403
1404 Vector tmp1_i(p + 2);
1405 DenseMatrix tmp1_ij(p + 2, p + 2);
1406 DenseMatrix tmp2_ij(p + 2, dim);
1407 DenseMatrix tmp3_ij(p + 2, dim);
1408 DenseTensor tmp1_ijk(p + 1, p + 1, dim);
1409 DenseTensor tmp2_ijk(p + 1, p + 1, dim);
1410 DenseTensor tmp3_ijk(p + 1, p + 1, dim);
1411 DenseTensor tmp4_ijk(p + 1, p + 2, dim);
1412 DenseTensor tmp5_ijk(p + 1, p + 2, dim);
1413 DenseTensor tmp6_ijk(p + 2, p + 2, dim);
1414 DenseTensor tmp7_ijk(p + 2, p + 2, dim);
1415 DenseMatrix u(dof, dim);
1416#endif
1417
1418 calcBasis(order, ip, tmp1_i, tmp1_ij, tmp2_ij,
1419 tmp1_ijk, tmp2_ijk, tmp3_ijk, tmp4_ijk, tmp5_ijk, tmp6_ijk,
1420 tmp7_ijk, tmp3_ij, u);
1421
1422 Ti.Mult(u, shape);
1423}
1424
1426 DenseMatrix &shape) const
1427{
1428#ifdef MFEM_THREAD_SAFE
1429 const int p = order - 1;
1430
1431 Vector tmp1_i(p + 2);
1432 DenseMatrix tmp1_ij(p + 2, p + 2);
1433 DenseMatrix tmp2_ij(p + 2, dim);
1434 DenseMatrix tmp3_ij(p + 2, dim);
1435 DenseTensor tmp1_ijk(p + 1, p + 1, dim);
1436 DenseTensor tmp2_ijk(p + 1, p + 1, dim);
1437 DenseTensor tmp3_ijk(p + 1, p + 1, dim);
1438 DenseTensor tmp4_ijk(p + 1, p + 2, dim);
1439 DenseTensor tmp5_ijk(p + 1, p + 2, dim);
1440 DenseTensor tmp6_ijk(p + 2, p + 2, dim);
1441 DenseTensor tmp7_ijk(p + 2, p + 2, dim);
1442#endif
1443
1444 calcBasis(order, ip, tmp1_i, tmp1_ij, tmp2_ij,
1445 tmp1_ijk, tmp2_ijk, tmp3_ijk, tmp4_ijk, tmp5_ijk, tmp6_ijk,
1446 tmp7_ijk, tmp3_ij, shape);
1447}
1448
1450 Vector &divshape) const
1451{
1452#ifdef MFEM_THREAD_SAFE
1453 const int p = order - 1;
1454
1455 Vector tmp1_i(p + 2);
1456 DenseMatrix tmp1_ij(p + 2, p + 2);
1457 DenseMatrix tmp2_ij(p + 2, dim);
1458 DenseMatrix tmp3_ij(p + 2, dim);
1459 DenseMatrix tmp4_ij(p + 1, p + 1);
1460 DenseTensor tmp1_ijk(p + 1, p + 1, dim);
1461 DenseTensor tmp2_ijk(p + 1, p + 1, dim);
1462 DenseTensor tmp3_ijk(p + 1, p + 1, dim);
1463 DenseTensor tmp4_ijk(p + 1, p + 2, dim);
1464 DenseTensor tmp5_ijk(p + 1, p + 2, dim);
1465 DenseTensor tmp6_ijk(p + 2, p + 2, dim);
1466 DenseTensor tmp7_ijk(p + 2, p + 2, dim);
1467 Vector divu(dof);
1468#endif
1469 divu = 0.0;
1470
1471 calcDivBasis(order, ip, tmp1_i, tmp1_ij, tmp2_ij,
1472 tmp1_ijk, tmp2_ijk, tmp3_ijk, tmp4_ij, tmp4_ijk, tmp5_ijk,
1473 tmp6_ijk, tmp7_ijk, tmp3_ij, divu);
1474
1475 Ti.Mult(divu, divshape);
1476}
1477
1479 Vector &dshape) const
1480{
1481#ifdef MFEM_THREAD_SAFE
1482 const int p = order - 1;
1483
1484 Vector tmp1_i(p + 2);
1485 DenseMatrix tmp1_ij(p + 2, p + 2);
1486 DenseMatrix tmp2_ij(p + 2, dim);
1487 DenseMatrix tmp3_ij(p + 2, dim);
1488 DenseMatrix tmp4_ij(p + 1, p + 1);
1489 DenseTensor tmp1_ijk(p + 1, p + 1, dim);
1490 DenseTensor tmp2_ijk(p + 1, p + 1, dim);
1491 DenseTensor tmp3_ijk(p + 1, p + 1, dim);
1492 DenseTensor tmp4_ijk(p + 1, p + 2, dim);
1493 DenseTensor tmp5_ijk(p + 1, p + 2, dim);
1494 DenseTensor tmp6_ijk(p + 2, p + 2, dim);
1495 DenseTensor tmp7_ijk(p + 2, p + 2, dim);
1496#endif
1497
1498 calcDivBasis(order, ip, tmp1_i, tmp1_ij, tmp2_ij,
1499 tmp1_ijk, tmp2_ijk, tmp3_ijk, tmp4_ij, tmp4_ijk, tmp5_ijk,
1500 tmp6_ijk, tmp7_ijk, tmp3_ij, dshape);
1501}
1502
1503void RT_FuentesPyramidElement::calcBasis(const int p,
1504 const IntegrationPoint &ip,
1505 Vector &phi_k,
1506 DenseMatrix &phi_ij,
1507 DenseMatrix &dphi_k,
1508 DenseTensor &VQ_ijk,
1509 DenseTensor &VT_ijk,
1510 DenseTensor &VTT_ijk,
1511 DenseTensor &E_ijk,
1512 DenseTensor &dE_ijk,
1513 DenseTensor &dphi_ijk,
1514 DenseTensor &VL_ijk,
1515 DenseMatrix &VR_ij,
1516 DenseMatrix &F) const
1517{
1518 real_t x = ip.x;
1519 real_t y = ip.y;
1520 real_t z = ip.z;
1521 Vector xy({x,y});
1522 real_t mu;
1523
1524 if (std::fabs(1.0 - z) < apex_tol)
1525 {
1526 z = 1.0 - apex_tol;
1527 y = 0.5 * (1.0 - z);
1528 x = 0.5 * (1.0 - z);
1529 xy(0) = x; xy(1) = y;
1530 }
1531 zmax = std::max(z, zmax);
1532
1533 F = 0.0;
1534
1535 int o = 0;
1536
1537 // Quadrilateral face
1538 if (z < 1.0)
1539 {
1540 V_Q(p, mu01(z, xy, 1), mu01_grad_mu01(z, xy, 1),
1541 mu01(z, xy, 2), mu01_grad_mu01(z, xy, 2),
1542 VQ_ijk);
1543
1544 const real_t muz3 = pow(mu0(z), 3);
1545
1546 for (int j=0; j<p; j++)
1547 for (int i=0; i<p; i++, o++)
1548 for (int k=0; k<3; k++)
1549 {
1550 F(o, k) = muz3 * VQ_ijk(i, j, k);
1551 }
1552 }
1553
1554 // Triangular faces
1555 if (z < 1.0)
1556 {
1557 Vector dmuz;
1558
1559 // (a,b) = (1,2), c = 0
1560 V_T(p, nu012(z, xy, 1), nu012_grad_nu012(z, xy, 1), VT_ijk);
1561 mu = mu0(z, xy, 2);
1562 dmuz.Destroy(); dmuz = grad_mu0(z, xy, 2);
1563 VT_T(p, nu012(z, xy, 1), nu01_grad_nu01(z, xy, 1),
1564 nu012_grad_nu012(z, xy, 1), mu, dmuz, VTT_ijk);
1565 for (int j=0; j<p; j++)
1566 for (int i=0; i+j<p; i++, o++)
1567 for (int k=0; k<3; k++)
1568 {
1569 F(o, k) = 0.5 * (mu * VT_ijk(i, j, k) + VTT_ijk(i, j, k));
1570 }
1571
1572 // (a,b) = (1,2), c = 1
1573 mu = mu1(z, xy, 2);
1574 dmuz.Destroy(); dmuz = grad_mu1(z, xy, 2);
1575 VT_T(p, nu012(z, xy, 1), nu01_grad_nu01(z, xy, 1),
1576 nu012_grad_nu012(z, xy, 1), mu, dmuz, VTT_ijk);
1577 for (int j=0; j<p; j++)
1578 for (int i=0; i+j<p; i++, o++)
1579 for (int k=0; k<3; k++)
1580 {
1581 F(o, k) = 0.5 * (mu * VT_ijk(i, j, k) + VTT_ijk(i, j, k));
1582 }
1583
1584 // (a,b) = (2,1), c = 0
1585 V_T(p, nu012(z, xy, 2), nu012_grad_nu012(z, xy, 2), VT_ijk);
1586 mu = mu0(z, xy, 1);
1587 dmuz.Destroy(); dmuz = grad_mu0(z, xy, 1);
1588 VT_T(p, nu012(z, xy, 2), nu01_grad_nu01(z, xy, 2),
1589 nu012_grad_nu012(z, xy, 2), mu, dmuz, VTT_ijk);
1590 for (int j=0; j<p; j++)
1591 for (int i=0; i+j<p; i++, o++)
1592 for (int k=0; k<3; k++)
1593 {
1594 F(o, k) = 0.5 * (mu * VT_ijk(i, j, k) + VTT_ijk(i, j, k));
1595 }
1596
1597 // (a,b) = (2,1), c = 1
1598 mu = mu1(z, xy, 1);
1599 dmuz.Destroy(); dmuz = grad_mu1(z, xy, 1);
1600 VT_T(p, nu012(z, xy, 2), nu01_grad_nu01(z, xy, 2),
1601 nu012_grad_nu012(z, xy, 2), mu, dmuz, VTT_ijk);
1602 for (int j=0; j<p; j++)
1603 for (int i=0; i+j<p; i++, o++)
1604 for (int k=0; k<3; k++)
1605 {
1606 F(o, k) = 0.5 * (mu * VT_ijk(i, j, k) + VTT_ijk(i, j, k));
1607 }
1608 }
1609
1610 // Interior
1611 // Family I
1612 if (z < 1.0 && p >= 2)
1613 {
1614 E_Q(p, mu01(z, xy, 1), grad_mu01(z, xy, 1),
1615 mu01(z, xy, 2), grad_mu01(z, xy, 2), E_ijk, dE_ijk);
1616 phi_E(p, mu01(z), grad_mu01(z), phi_k, dphi_k);
1617 const real_t muz = mu0(z);
1618 const Vector dmuz(grad_mu0(z));
1619
1620 Vector dmuphi(3), E(3), v(3);
1621
1622 for (int k=2; k<=p; k++)
1623 {
1624 dmuphi(0) = muz * dphi_k(k,0) + dmuz(0) * phi_k(k);
1625 dmuphi(1) = muz * dphi_k(k,1) + dmuz(1) * phi_k(k);
1626 dmuphi(2) = muz * dphi_k(k,2) + dmuz(2) * phi_k(k);
1627 for (int j=2; j<=p; j++)
1628 for (int i=0; i<p; i++, o++)
1629 {
1630 E(0) = E_ijk(i,j,0); E(1) = E_ijk(i,j,1); E(2) = E_ijk(i,j,2);
1631 dmuphi.cross3D(E, v);
1632 for (int l=0; l<3; l++)
1633 {
1634 F(o, l) = muz * phi_k(k) * dE_ijk(i,j,l) + v(l);
1635 }
1636 }
1637 }
1638 }
1639
1640 // Family II
1641 if (z < 1.0 && p >= 2)
1642 {
1643 E_Q(p, mu01(z, xy, 2), grad_mu01(z, xy, 2),
1644 mu01(z, xy, 1), grad_mu01(z, xy, 1), E_ijk, dE_ijk);
1645 // Re-using phi_E from Family I
1646 const real_t muz = mu0(z);
1647 const Vector dmuz(grad_mu0(z));
1648
1649 Vector dmuphi(3), E(3), v(3);
1650
1651 for (int k=2; k<=p; k++)
1652 {
1653 dmuphi(0) = muz * dphi_k(k,0) + dmuz(0) * phi_k(k);
1654 dmuphi(1) = muz * dphi_k(k,1) + dmuz(1) * phi_k(k);
1655 dmuphi(2) = muz * dphi_k(k,2) + dmuz(2) * phi_k(k);
1656 for (int j=2; j<=p; j++)
1657 for (int i=0; i<p; i++, o++)
1658 {
1659 E(0) = E_ijk(i,j,0); E(1) = E_ijk(i,j,1); E(2) = E_ijk(i,j,2);
1660 dmuphi.cross3D(E, v);
1661 for (int l=0; l<3; l++)
1662 {
1663 F(o, l) = muz * phi_k(k) * dE_ijk(i,j,l) + v(l);
1664 }
1665 }
1666 }
1667 }
1668 // Family III
1669 if (z < 1.0 && p >= 2)
1670 {
1671 phi_Q(p, mu01(z, xy, 2), grad_mu01(z, xy, 2),
1672 mu01(z, xy, 1), grad_mu01(z, xy, 1), phi_ij, dphi_ijk);
1673 const real_t muz = mu0(z);
1674 const Vector dmuz(grad_mu0(z));
1675
1676 for (int j=2; j<=p; j++)
1677 for (int i=2; i<=p; i++, o++)
1678 {
1679 const int n = std::max(i,j);
1680 const real_t nmu = n * pow(muz, n-1);
1681 F(o, 0) = nmu * (dphi_ijk(i,j,1) * dmuz(2) -
1682 dphi_ijk(i,j,2) * dmuz(1));
1683 F(o, 1) = nmu * (dphi_ijk(i,j,2) * dmuz(0) -
1684 dphi_ijk(i,j,0) * dmuz(2));
1685 F(o, 2) = nmu * (dphi_ijk(i,j,0) * dmuz(1) -
1686 dphi_ijk(i,j,1) * dmuz(0));
1687 }
1688 }
1689 // Family IV
1690 if (z < 1.0 && p >= 2)
1691 {
1692 // Re-using V_Q from Quadrilateral Face
1693 phi_E(p, mu01(z), phi_k);
1694
1695 const real_t muz2 = pow(mu0(z), 2);
1696
1697 for (int k=2; k<=p; k++)
1698 for (int j=0; j<p; j++)
1699 for (int i=0; i<p; i++, o++)
1700 for (int l=0; l<3; l++)
1701 {
1702 F(o, l) = muz2 * VQ_ijk(i, j, l) * phi_k(k);
1703 }
1704
1705 }
1706 // Family V
1707 if (z < 1.0 && p >= 2)
1708 {
1709 V_L(p, mu01(z, xy, 1), grad_mu01(z, xy, 1),
1710 mu01(z, xy, 2), grad_mu01(z, xy, 2), mu0(z), grad_mu0(z), VL_ijk);
1711
1712 const real_t muz = mu1(z);
1713
1714 for (int j=2; j<=p; j++)
1715 for (int i=2; i<=p; i++, o++)
1716 {
1717 const int n = std::max(i, j);
1718 const real_t muzi = pow(muz, n-1);
1719 for (int l=0; l<3; l++)
1720 {
1721 F(o, l) = muzi * VL_ijk(i, j, l);
1722 }
1723 }
1724 }
1725 // Family VI
1726 if (z < 1.0 && p >= 2)
1727 {
1728 V_R(p, mu01(z, xy, 1), grad_mu01(z, xy, 1),
1729 mu1(z, xy, 2), grad_mu1(z, xy, 2), mu0(z), grad_mu0(z), VR_ij);
1730
1731 const real_t muz = mu1(z);
1732
1733 for (int i=2; i<=p; i++, o++)
1734 {
1735 const real_t muzi = pow(muz, i-1);
1736 for (int l=0; l<3; l++)
1737 {
1738 F(o, l) = muzi * VR_ij(i, l);
1739 }
1740 }
1741 }
1742 // Family VII
1743 if (z < 1.0 && p >= 2)
1744 {
1745 V_R(p, mu01(z, xy, 2), grad_mu01(z, xy, 2),
1746 mu1(z, xy, 1), grad_mu1(z,xy,1), mu0(z), grad_mu0(z), VR_ij);
1747
1748 const real_t muz = mu1(z);
1749
1750 for (int i=2; i<=p; i++, o++)
1751 {
1752 const real_t muzi = pow(muz, i-1);
1753 for (int l=0; l<3; l++)
1754 {
1755 F(o, l) = muzi * VR_ij(i, l);
1756 }
1757 }
1758 }
1759}
1760
1761void RT_FuentesPyramidElement::calcDivBasis(const int p,
1762 const IntegrationPoint &ip,
1763 Vector &phi_k,
1764 DenseMatrix &phi_ij,
1765 DenseMatrix &dphi_k,
1766 DenseTensor &VQ_ijk,
1767 DenseTensor &VT_ijk,
1768 DenseTensor &VTT_ijk,
1769 DenseMatrix &dVTT_ij,
1770 DenseTensor &E_ijk,
1771 DenseTensor &dE_ijk,
1772 DenseTensor &dphi_ijk,
1773 DenseTensor &VL_ijk,
1774 DenseMatrix &VR_ij,
1775 Vector &dF) const
1776{
1777 real_t x = ip.x;
1778 real_t y = ip.y;
1779 real_t z = ip.z;
1780 Vector xy({x,y});
1781 real_t mu;
1782
1783 bool limz1 = false;
1784 if (std::fabs(1.0 - z) < apex_tol)
1785 {
1786 limz1 = true;
1787 z = 1.0 - apex_tol;
1788 y = 0.5 * (1.0 - z);
1789 x = 0.5 * (1.0 - z);
1790 xy(0) = x; xy(1) = y;
1791 }
1792 zmax = std::max(z, zmax);
1793
1794 dF = 0.0;
1795
1796 int o = 0;
1797
1798 // Quadrilateral face
1799 {
1800 V_Q(p, mu01(z, xy, 1), mu01_grad_mu01(z, xy, 1),
1801 mu01(z, xy, 2), mu01_grad_mu01(z, xy, 2),
1802 VQ_ijk);
1803
1804 const real_t muz2 = pow(mu0(z), 2);
1805 const Vector dmuz = grad_mu0(z);
1806
1807 const int o0 = o;
1808 for (int j=0; j<p; j++)
1809 for (int i=0; i<p; i++, o++)
1810 for (int k=0; k<3; k++)
1811 {
1812 dF(o) += 3.0 * muz2 * dmuz(k) * VQ_ijk(i, j, k);
1813 }
1814
1815 // Overwrite lowest order quadrilateral face DoF with known limiting
1816 // value
1817 if (limz1)
1818 {
1819 dF(o0) = -3.0;
1820 }
1821 }
1822
1823 // Triangular faces
1824 {
1825 Vector dmuz;
1826
1827 // (a,b) = (1,2), c = 0
1828 V_T(p, nu012(z, xy, 1), nu012_grad_nu012(z, xy, 1), VT_ijk);
1829 mu = mu0(z, xy, 2);
1830 dmuz.Destroy(); dmuz = grad_mu0(z, xy, 2);
1831 VT_T(p, nu012(z, xy, 1), nu01_grad_nu01(z, xy, 1),
1832 nu012_grad_nu012(z, xy, 1), grad_nu2(z, xy, 1), mu, dmuz,
1833 VTT_ijk, dVTT_ij);
1834 const int o1 = o;
1835 for (int j=0; j<p; j++)
1836 for (int i=0; i+j<p; i++, o++)
1837 {
1838 dF(o) = 0.5 * dVTT_ij(i, j);
1839 for (int k=0; k<3; k++)
1840 {
1841 dF(o) += 0.5 * dmuz(k) * VT_ijk(i, j, k);
1842 }
1843 }
1844
1845 // (a,b) = (1,2), c = 1
1846 mu = mu1(z, xy, 2);
1847 dmuz.Destroy(); dmuz = grad_mu1(z, xy, 2);
1848 VT_T(p, nu012(z, xy, 1), nu01_grad_nu01(z, xy, 1),
1849 nu012_grad_nu012(z, xy, 1), grad_nu2(z, xy, 1), mu, dmuz,
1850 VTT_ijk, dVTT_ij);
1851 const int o2 = o;
1852 for (int j=0; j<p; j++)
1853 for (int i=0; i+j<p; i++, o++)
1854 {
1855 dF(o) = 0.5 * dVTT_ij(i, j);
1856 for (int k=0; k<3; k++)
1857 {
1858 dF(o) += 0.5 * dmuz(k) * VT_ijk(i, j, k);
1859 }
1860 }
1861
1862 // (a,b) = (2,1), c = 0
1863 V_T(p, nu012(z, xy, 2), nu012_grad_nu012(z, xy, 2), VT_ijk);
1864 mu = mu0(z, xy, 1);
1865 dmuz.Destroy(); dmuz = grad_mu0(z, xy, 1);
1866 VT_T(p, nu012(z, xy, 2), nu01_grad_nu01(z, xy, 2),
1867 nu012_grad_nu012(z, xy, 2), grad_nu2(z, xy, 2), mu, dmuz,
1868 VTT_ijk, dVTT_ij);
1869 const int o3 = o;
1870 for (int j=0; j<p; j++)
1871 for (int i=0; i+j<p; i++, o++)
1872 {
1873 dF(o) = 0.5 * dVTT_ij(i, j);
1874 for (int k=0; k<3; k++)
1875 {
1876 dF(o) += 0.5 * dmuz(k) * VT_ijk(i, j, k);
1877 }
1878 }
1879
1880 // (a,b) = (2,1), c = 1
1881 mu = mu1(z, xy, 1);
1882 dmuz.Destroy(); dmuz = grad_mu1(z, xy, 1);
1883 VT_T(p, nu012(z, xy, 2), nu01_grad_nu01(z, xy, 2),
1884 nu012_grad_nu012(z, xy, 2), grad_nu2(z, xy, 2), mu, dmuz,
1885 VTT_ijk, dVTT_ij);
1886 const int o4 = o;
1887 for (int j=0; j<p; j++)
1888 for (int i=0; i+j<p; i++, o++)
1889 {
1890 dF(o) = 0.5 * dVTT_ij(i, j);
1891 for (int k=0; k<3; k++)
1892 {
1893 dF(o) += 0.5 * dmuz(k) * VT_ijk(i, j, k);
1894 }
1895 }
1896
1897 // Overwrite lowest order triangular face DoFs with known limiting values
1898 if (limz1)
1899 {
1900 dF(o1) = 1.5;
1901 dF(o2) = -1.5;
1902 dF(o3) = -1.5;
1903 dF(o4) = 1.5;
1904 }
1905 }
1906
1907 // Interior
1908 // Family I
1909 if (p >= 2)
1910 {
1911 // Divergence is zero so skip ahead
1912 o += (p-1) * (p-1) * p;
1913 }
1914
1915 // Family II
1916 if (p >= 2)
1917 {
1918 // Divergence is zero so skip ahead
1919 o += (p-1) * (p-1) * p;
1920 }
1921 // Family III
1922 if (p >= 2)
1923 {
1924 // Divergence is zero so skip ahead
1925 o += (p-1) * (p-1);
1926 }
1927 // Family IV
1928 if (p >= 2)
1929 {
1930 // Re-using V_Q from Quadrilateral Face
1931 phi_E(p, mu01(z), grad_mu01(z), phi_k, dphi_k);
1932
1933 const real_t muz2 = pow(mu0(z), 2);
1934 const Vector dmuz = grad_mu0(z);
1935
1936 for (int k=2; k<=p; k++)
1937 for (int j=0; j<p; j++)
1938 for (int i=0; i<p; i++, o++)
1939 for (int l=0; l<3; l++)
1940 {
1941 dF(o) += (muz2 * dphi_k(k, l) +
1942 2.0 * mu0(z) * phi_k(k) * dmuz(l)) * VQ_ijk(i, j, l);
1943 }
1944 }
1945 // Family V
1946 if (p >= 2)
1947 {
1948 V_L(p, mu01(z, xy, 1), grad_mu01(z, xy, 1),
1949 mu01(z, xy, 2), grad_mu01(z, xy, 2), mu0(z), grad_mu0(z), VL_ijk);
1950
1951 const real_t muz = mu1(z);
1952 const Vector dmuz = grad_mu1(z);
1953
1954 for (int j=2; j<=p; j++)
1955 for (int i=2; i<=p; i++, o++)
1956 {
1957 const int n = std::max(i, j);
1958 const real_t muzi = pow(muz, n-2);
1959 for (int l=0; l<3; l++)
1960 {
1961 dF(o) += (n-1) * muzi * dmuz(l) * VL_ijk(i, j, l);
1962 }
1963 }
1964 }
1965 // Family VI
1966 if (p >= 2)
1967 {
1968 V_R(p, mu01(z, xy, 1), grad_mu01(z, xy, 1),
1969 mu1(z, xy, 2), grad_mu1(z, xy, 2), mu0(z), grad_mu0(z), VR_ij);
1970
1971 const real_t muz = mu1(z);
1972 const Vector dmuz = grad_mu1(z);
1973
1974 for (int i=2; i<=p; i++, o++)
1975 {
1976 const real_t muzi = pow(muz, i-2);
1977 for (int l=0; l<3; l++)
1978 {
1979 dF(o) += (i-1) * muzi * dmuz(l) * VR_ij(i, l);
1980 }
1981 }
1982 }
1983 // Family VII
1984 if (p >= 2)
1985 {
1986 V_R(p, mu01(z, xy, 2), grad_mu01(z, xy, 2),
1987 mu1(z, xy, 1), grad_mu1(z,xy,1), mu0(z), grad_mu0(z), VR_ij);
1988
1989 const real_t muz = mu1(z);
1990 const Vector dmuz = grad_mu1(z);
1991
1992 for (int i=2; i<=p; i++, o++)
1993 {
1994 const real_t muzi = pow(muz, i-2);
1995 for (int l=0; l<3; l++)
1996 {
1997 dF(o) += (i-1) * muzi * dmuz(l) * VR_ij(i, l);
1998 }
1999 }
2000 }
2001}
2002
2003const real_t RT_R1D_SegmentElement::nk[9] = { 1.,0.,0., 0.,1.,0., 0.,0.,1. };
2004
2006 const int cb_type,
2007 const int ob_type)
2008 : VectorFiniteElement(1, Geometry::SEGMENT, 3 * p + 4, p + 1,
2009 H_DIV, FunctionSpace::Pk),
2010 dof2nk(dof),
2011 cbasis1d(poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
2012 obasis1d(poly1d.GetBasis(p, VerifyOpen(ob_type)))
2013{
2014 // Override default dimension for VectorFiniteElements
2015 vdim = 3;
2016
2017 const real_t *cp = poly1d.ClosedPoints(p + 1, cb_type);
2018 const real_t *op = poly1d.OpenPoints(p, ob_type);
2019
2020#ifndef MFEM_THREAD_SAFE
2021 shape_cx.SetSize(p + 2);
2022 shape_ox.SetSize(p + 1);
2023 dshape_cx.SetSize(p + 2);
2024#endif
2025
2026 dof_map.SetSize(dof);
2027
2028 int o = 0;
2029 // nodes
2030 // (0)
2031 Nodes.IntPoint(o).x = cp[0]; // x-directed
2032 dof_map[0] = o; dof2nk[o++] = 0;
2033
2034 // (1)
2035 Nodes.IntPoint(o).x = cp[p+1]; // x-directed
2036 dof_map[p+1] = o; dof2nk[o++] = 0;
2037
2038 // interior
2039 // x-components
2040 for (int i = 1; i <= p; i++)
2041 {
2042 Nodes.IntPoint(o).x = cp[i];
2043 dof_map[i] = o; dof2nk[o++] = 0;
2044 }
2045 // y-components
2046 for (int i = 0; i <= p; i++)
2047 {
2048 Nodes.IntPoint(o).x = op[i];
2049 dof_map[p+i+2] = o; dof2nk[o++] = 1;
2050 }
2051 // z-components
2052 for (int i = 0; i <= p; i++)
2053 {
2054 Nodes.IntPoint(o).x = op[i];
2055 dof_map[2*p+3+i] = o; dof2nk[o++] = 2;
2056 }
2057}
2058
2060 DenseMatrix &shape) const
2061{
2062 const int p = order;
2063
2064#ifdef MFEM_THREAD_SAFE
2065 Vector shape_cx(p + 1), shape_ox(p);
2066#endif
2067
2068 cbasis1d.Eval(ip.x, shape_cx);
2069 obasis1d.Eval(ip.x, shape_ox);
2070
2071 int o = 0;
2072 // x-components
2073 for (int i = 0; i <= p; i++)
2074 {
2075 int idx = dof_map[o++];
2076 shape(idx,0) = shape_cx(i);
2077 shape(idx,1) = 0.;
2078 shape(idx,2) = 0.;
2079 }
2080 // y-components
2081 for (int i = 0; i < p; i++)
2082 {
2083 int idx = dof_map[o++];
2084 shape(idx,0) = 0.;
2085 shape(idx,1) = shape_ox(i);
2086 shape(idx,2) = 0.;
2087 }
2088 // z-components
2089 for (int i = 0; i < p; i++)
2090 {
2091 int idx = dof_map[o++];
2092 shape(idx,0) = 0.;
2093 shape(idx,1) = 0.;
2094 shape(idx,2) = shape_ox(i);
2095 }
2096}
2097
2099 DenseMatrix &shape) const
2100{
2101 CalcVShape(Trans.GetIntPoint(), shape);
2102 const DenseMatrix & J = Trans.Jacobian();
2103 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
2104 "RT_R1D_SegmentElement cannot be embedded in "
2105 "2 or 3 dimensional spaces");
2106 for (int i=0; i<dof; i++)
2107 {
2108 shape(i, 0) *= J(0,0);
2109 }
2110 shape *= (1.0 / Trans.Weight());
2111}
2112
2114 Vector &divshape) const
2115{
2116 const int p = order;
2117
2118#ifdef MFEM_THREAD_SAFE
2119 Vector shape_cx(p + 1);
2120 Vector dshape_cx(p + 1);
2121#endif
2122
2123 cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
2124
2125 int o = 0;
2126 // x-components
2127 for (int i = 0; i <= p; i++)
2128 {
2129 int idx = dof_map[o++];
2130 divshape(idx) = dshape_cx(i);
2131 }
2132 // y-components
2133 for (int i = 0; i < p; i++)
2134 {
2135 int idx = dof_map[o++];
2136 divshape(idx) = 0.;
2137 }
2138 // z-components
2139 for (int i = 0; i < p; i++)
2140 {
2141 int idx = dof_map[o++];
2142 divshape(idx) = 0.;
2143 }
2144}
2145
2147 ElementTransformation &Trans,
2148 Vector &dofs) const
2149{
2150 real_t data[3];
2151 Vector vk1(data, 1);
2152 Vector vk3(data, 3);
2153
2154 real_t * nk_ptr = const_cast<real_t*>(nk);
2155
2156 for (int k = 0; k < dof; k++)
2157 {
2158 Trans.SetIntPoint(&Nodes.IntPoint(k));
2159
2160 vc.Eval(vk3, Trans, Nodes.IntPoint(k));
2161 // dof_k = nk^t adj(J) vk
2162 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
2163 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
2164
2165 dofs(k) = Trans.AdjugateJacobian().InnerProduct(vk1, n1) +
2166 Trans.Weight() * vk3(1) * n3(1) +
2167 Trans.Weight() * vk3(2) * n3(2);
2168 }
2169}
2170
2172 ElementTransformation &Trans,
2173 DenseMatrix &I) const
2174{
2175 if (fe.GetRangeType() == SCALAR)
2176 {
2178 Vector shape(fe.GetDof());
2179
2180 real_t * nk_ptr = const_cast<real_t*>(nk);
2181
2182 I.SetSize(dof, vdim*fe.GetDof());
2183 for (int k = 0; k < dof; k++)
2184 {
2185 const IntegrationPoint &ip = Nodes.IntPoint(k);
2186
2187 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
2188 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
2189
2190 fe.CalcShape(ip, shape);
2191 Trans.SetIntPoint(&ip);
2192 // Transform RT face normals from reference to physical space
2193 // vk = adj(J)^T nk
2194 Trans.AdjugateJacobian().MultTranspose(n1, vk);
2195 vk[1] = n3[1] * Trans.Weight();
2196 vk[2] = n3[2] * Trans.Weight();
2197 if (fe.GetMapType() == INTEGRAL)
2198 {
2199 real_t w = 1.0/Trans.Weight();
2200 for (int d = 0; d < 1; d++)
2201 {
2202 vk[d] *= w;
2203 }
2204 }
2205
2206 for (int j = 0; j < shape.Size(); j++)
2207 {
2208 real_t s = shape(j);
2209 if (fabs(s) < 1e-12)
2210 {
2211 s = 0.0;
2212 }
2213 // Project scalar basis function multiplied by each coordinate
2214 // direction onto the transformed face normals
2215 for (int d = 0; d < vdim; d++)
2216 {
2217 I(k,j+d*shape.Size()) = s*vk[d];
2218 }
2219 }
2220 }
2221 }
2222 else
2223 {
2226
2227 real_t * nk_ptr = const_cast<real_t*>(nk);
2228
2229 I.SetSize(dof, fe.GetDof());
2230 for (int k = 0; k < dof; k++)
2231 {
2232 const IntegrationPoint &ip = Nodes.IntPoint(k);
2233
2234 Vector n1(&nk_ptr[dof2nk[k] * 3], 1);
2235 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
2236
2237 Trans.SetIntPoint(&ip);
2238 // Transform RT face normals from reference to physical space
2239 // vk = adj(J)^T nk
2240 Trans.AdjugateJacobian().MultTranspose(n1, vk);
2241 // Compute fe basis functions in physical space
2242 fe.CalcVShape(Trans, vshape);
2243 // Project fe basis functions onto transformed face normals
2244 for (int j=0; j<vshape.Height(); j++)
2245 {
2246 I(k, j) = 0.0;
2247 I(k, j) += vshape(j, 0) * vk[0];
2248 if (vshape.Width() == 3)
2249 {
2250 I(k, j) += Trans.Weight() * vshape(j, 1) * n3(1);
2251 I(k, j) += Trans.Weight() * vshape(j, 2) * n3(2);
2252 }
2253 }
2254 }
2255 }
2256}
2257
2259 ElementTransformation &Trans,
2260 DenseMatrix &curl) const
2261{
2262 DenseMatrix curl_shape(fe.GetDof(), fe.GetRangeDim());
2263 Vector curl_k(fe.GetDof());
2264
2265 real_t * nk_ptr = const_cast<real_t*>(nk);
2266
2267 curl.SetSize(dof, fe.GetDof());
2268 for (int k = 0; k < dof; k++)
2269 {
2270 fe.CalcCurlShape(Nodes.IntPoint(k), curl_shape);
2271 curl_shape.Mult(nk_ptr + dof2nk[k] * 3, curl_k);
2272 for (int j = 0; j < curl_k.Size(); j++)
2273 {
2274 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
2275 }
2276 }
2277}
2278
2279const real_t RT_R2D_SegmentElement::nk[2] = { 0.,1.};
2280
2282 const int ob_type)
2283 : VectorFiniteElement(1, Geometry::SEGMENT, p + 1, p + 1,
2284 H_DIV, FunctionSpace::Pk),
2285 dof2nk(dof),
2286 obasis1d(poly1d.GetBasis(p, VerifyOpen(ob_type)))
2287{
2288 // Override default dimension for VectorFiniteElements
2289 vdim = 2;
2290
2291 const real_t *op = poly1d.OpenPoints(p, ob_type);
2292
2293#ifndef MFEM_THREAD_SAFE
2294 shape_ox.SetSize(p+1);
2295#endif
2296
2297 dof_map.SetSize(dof);
2298
2299 int o = 0;
2300 // interior
2301 // z-components
2302 for (int i = 0; i <= p; i++)
2303 {
2304 Nodes.IntPoint(o).x = op[i];
2305 dof_map[i] = o; dof2nk[o++] = 0;
2306 }
2307}
2308
2310 DenseMatrix &shape) const
2311{
2312 const int p = order;
2313
2314#ifdef MFEM_THREAD_SAFE
2315 Vector shape_ox(p);
2316#endif
2317
2318 obasis1d.Eval(ip.x, shape_ox);
2319
2320 int o = 0;
2321 // z-components
2322 for (int i = 0; i <= p; i++)
2323 {
2324 int idx = dof_map[o++];
2325 shape(idx,0) = shape_ox(i);
2326 shape(idx,1) = 0.;
2327 }
2328}
2329
2331 DenseMatrix &shape) const
2332{
2333 CalcVShape(Trans.GetIntPoint(), shape);
2334 const DenseMatrix & J = Trans.Jacobian();
2335 MFEM_ASSERT(J.Width() == 1 && J.Height() == 1,
2336 "RT_R2D_SegmentElement cannot be embedded in "
2337 "2 or 3 dimensional spaces");
2338 for (int i=0; i<dof; i++)
2339 {
2340 shape(i, 0) *= J(0,0);
2341 }
2342 shape *= (1.0 / Trans.Weight());
2343}
2344
2346 Vector &div_shape) const
2347{
2348 div_shape = 0.0;
2349}
2350
2351void RT_R2D_SegmentElement::LocalInterpolation(const VectorFiniteElement &cfe,
2352 ElementTransformation &Trans,
2353 DenseMatrix &I) const
2354{
2355 real_t vk[Geometry::MaxDim]; vk[1] = 0.0; vk[2] = 0.0;
2356 Vector xk(vk, dim);
2358 DenseMatrix vshape(cfe.GetDof(), vdim);
2359
2360 real_t * nk_ptr = const_cast<real_t*>(nk);
2361
2362 I.SetSize(dof, vshape.Height());
2363
2364 // assuming Trans is linear; this should be ok for all refinement types
2366 const DenseMatrix &adjJ = Trans.AdjugateJacobian();
2367 for (int k = 0; k < dof; k++)
2368 {
2369 Vector n2(&nk_ptr[dof2nk[k] * 2], 2);
2370
2371 Trans.Transform(Nodes.IntPoint(k), xk);
2372 ip.Set3(vk);
2373 cfe.CalcVShape(ip, vshape);
2374 // xk = |J| J^{-t} n_k
2375 adjJ.MultTranspose(n2, vk);
2376 // I_k = vshape_k.adj(J)^t.n_k, k=1,...,dof
2377 for (int j = 0; j < vshape.Height(); j++)
2378 {
2379 real_t Ikj = 0.;
2380 /*
2381 for (int i = 0; i < dim; i++)
2382 {
2383 Ikj += vshape(j, i) * vk[i];
2384 }
2385 */
2386 Ikj += Trans.Weight() * vshape(j, 1) * n2(1);
2387 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
2388 }
2389 }
2390}
2391
2393 const real_t *nk_fe)
2394 : VectorFiniteElement(2, G, Do, p + 1,
2395 H_DIV, FunctionSpace::Pk),
2396 nk(nk_fe),
2397 dof_map(dof),
2398 dof2nk(dof)
2399{
2400 // Override default dimension for VectorFiniteElements
2401 vdim = 3;
2402}
2403
2405 DenseMatrix &shape) const
2406{
2407 CalcVShape(Trans.GetIntPoint(), shape);
2408 const DenseMatrix & J = Trans.Jacobian();
2409 MFEM_ASSERT(J.Width() == 2 && J.Height() == 2,
2410 "RT_R2D_FiniteElement cannot be embedded in "
2411 "3 dimensional spaces");
2412 for (int i=0; i<dof; i++)
2413 {
2414 real_t sx = shape(i, 0);
2415 real_t sy = shape(i, 1);
2416 shape(i, 0) = sx * J(0, 0) + sy * J(0, 1);
2417 shape(i, 1) = sx * J(1, 0) + sy * J(1, 1);
2418 }
2419 shape *= (1.0 / Trans.Weight());
2420}
2421
2422void
2423RT_R2D_FiniteElement::LocalInterpolation(const VectorFiniteElement &cfe,
2424 ElementTransformation &Trans,
2425 DenseMatrix &I) const
2426{
2427 real_t vk[Geometry::MaxDim]; vk[2] = 0.0;
2428 Vector xk(vk, dim);
2430 DenseMatrix vshape(cfe.GetDof(), vdim);
2431
2432 real_t * nk_ptr = const_cast<real_t*>(nk);
2433
2434 I.SetSize(dof, vshape.Height());
2435
2436 // assuming Trans is linear; this should be ok for all refinement types
2438 const DenseMatrix &adjJ = Trans.AdjugateJacobian();
2439 for (int k = 0; k < dof; k++)
2440 {
2441 Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
2442 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
2443
2444 Trans.Transform(Nodes.IntPoint(k), xk);
2445 ip.Set3(vk);
2446 cfe.CalcVShape(ip, vshape);
2447 // xk = |J| J^{-t} n_k
2448 adjJ.MultTranspose(n2, vk);
2449 // I_k = vshape_k.adj(J)^t.n_k, k=1,...,dof
2450 for (int j = 0; j < vshape.Height(); j++)
2451 {
2452 real_t Ikj = 0.;
2453 for (int i = 0; i < dim; i++)
2454 {
2455 Ikj += vshape(j, i) * vk[i];
2456 }
2457 Ikj += Trans.Weight() * vshape(j, 2) * n3(2);
2458 I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
2459 }
2460 }
2461}
2462
2464 DenseMatrix &R) const
2465{
2466 real_t pt_data[Geometry::MaxDim];
2468 Vector pt(pt_data, dim);
2469
2470#ifdef MFEM_THREAD_SAFE
2472#endif
2473
2474 real_t * nk_ptr = const_cast<real_t*>(nk);
2475
2477 const DenseMatrix &J = Trans.Jacobian();
2478 const real_t weight = Trans.Weight();
2479 for (int j = 0; j < dof; j++)
2480 {
2481 Vector n2(&nk_ptr[dof2nk[j] * 3], 2);
2482 Vector n3(&nk_ptr[dof2nk[j] * 3], 3);
2483
2484 InvertLinearTrans(Trans, Nodes.IntPoint(j), pt);
2485 ip.Set(pt_data, dim);
2486 if (Geometries.CheckPoint(geom_type, ip)) // do we need an epsilon here?
2487 {
2488 CalcVShape(ip, vshape);
2489 J.MultTranspose(n2, pt_data);
2490 pt /= weight;
2491 for (int k = 0; k < dof; k++)
2492 {
2493 real_t R_jk = 0.0;
2494 for (int d = 0; d < dim; d++)
2495 {
2496 R_jk += vshape(k,d)*pt_data[d];
2497 }
2498 R_jk += vshape(k,2) * n3(2);
2499 R(j,k) = R_jk;
2500 }
2501 }
2502 else
2503 {
2504 // Set the whole row to avoid valgrind warnings in R.Threshold().
2505 R.SetRow(j, infinity());
2506 }
2507 }
2508 R.Threshold(1e-12);
2509}
2510
2512 ElementTransformation &Trans,
2513 Vector &dofs) const
2514{
2515 real_t data[3];
2516 Vector vk2(data, 2);
2517 Vector vk3(data, 3);
2518
2519 real_t * nk_ptr = const_cast<real_t*>(nk);
2520
2521 for (int k = 0; k < dof; k++)
2522 {
2523 Trans.SetIntPoint(&Nodes.IntPoint(k));
2524
2525 vc.Eval(vk3, Trans, Nodes.IntPoint(k));
2526 // dof_k = nk^t adj(J) vk
2527 Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
2528 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
2529
2530 dofs(k) = Trans.AdjugateJacobian().InnerProduct(vk2, n2) +
2531 Trans.Weight() * vk3(2) * n3(2);
2532 }
2533}
2534
2536 ElementTransformation &Trans,
2537 DenseMatrix &I) const
2538{
2539 if (fe.GetRangeType() == SCALAR)
2540 {
2542 Vector shape(fe.GetDof());
2543
2544 real_t * nk_ptr = const_cast<real_t*>(nk);
2545
2546 I.SetSize(dof, vdim*fe.GetDof());
2547 for (int k = 0; k < dof; k++)
2548 {
2549 const IntegrationPoint &ip = Nodes.IntPoint(k);
2550
2551 Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
2552 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
2553
2554 fe.CalcShape(ip, shape);
2555 Trans.SetIntPoint(&ip);
2556 // Transform RT face normals from reference to physical space
2557 // vk = adj(J)^T nk
2558 Trans.AdjugateJacobian().MultTranspose(n2, vk);
2559 vk[2] = n3[2] * Trans.Weight();
2560 if (fe.GetMapType() == INTEGRAL)
2561 {
2562 real_t w = 1.0/Trans.Weight();
2563 for (int d = 0; d < 2; d++)
2564 {
2565 vk[d] *= w;
2566 }
2567 }
2568
2569 for (int j = 0; j < shape.Size(); j++)
2570 {
2571 real_t s = shape(j);
2572 if (fabs(s) < 1e-12)
2573 {
2574 s = 0.0;
2575 }
2576 // Project scalar basis function multiplied by each coordinate
2577 // direction onto the transformed face normals
2578 for (int d = 0; d < vdim; d++)
2579 {
2580 I(k,j+d*shape.Size()) = s*vk[d];
2581 }
2582 }
2583 }
2584 }
2585 else
2586 {
2589
2590 real_t * nk_ptr = const_cast<real_t*>(nk);
2591
2592 I.SetSize(dof, fe.GetDof());
2593 for (int k = 0; k < dof; k++)
2594 {
2595 const IntegrationPoint &ip = Nodes.IntPoint(k);
2596
2597 Vector n2(&nk_ptr[dof2nk[k] * 3], 2);
2598 Vector n3(&nk_ptr[dof2nk[k] * 3], 3);
2599
2600 Trans.SetIntPoint(&ip);
2601 // Transform RT face normals from reference to physical space
2602 // vk = adj(J)^T nk
2603 Trans.AdjugateJacobian().MultTranspose(n2, vk);
2604 // Compute fe basis functions in physical space
2605 fe.CalcVShape(Trans, vshape);
2606 // Project fe basis functions onto transformed face normals
2607 for (int j=0; j<vshape.Height(); j++)
2608 {
2609 I(k, j) = 0.0;
2610 for (int i=0; i<2; i++)
2611 {
2612 I(k, j) += vshape(j, i) * vk[i];
2613 }
2614 if (vshape.Width() == 3)
2615 {
2616 I(k, j) += Trans.Weight() * vshape(j, 2) * n3(2);
2617 }
2618 }
2619 }
2620 }
2621}
2622
2624 ElementTransformation &Trans,
2625 DenseMatrix &curl) const
2626{
2627 DenseMatrix curl_shape(fe.GetDof(), fe.GetRangeDim());
2628 Vector curl_k(fe.GetDof());
2629
2630 real_t * nk_ptr = const_cast<real_t*>(nk);
2631
2632 curl.SetSize(dof, fe.GetDof());
2633 for (int k = 0; k < dof; k++)
2634 {
2635 fe.CalcCurlShape(Nodes.IntPoint(k), curl_shape);
2636 curl_shape.Mult(nk_ptr + dof2nk[k] * 3, curl_k);
2637 for (int j = 0; j < curl_k.Size(); j++)
2638 {
2639 curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
2640 }
2641 }
2642}
2643
2644const real_t RT_R2D_TriangleElement::nk_t[12] =
2645{ 0.,-1.,0., 1.,1.,0., -1.,0.,0., 0.,0.,1. };
2646
2648 : RT_R2D_FiniteElement(p, Geometry::TRIANGLE, ((p + 1)*(3 * p + 8))/2, nk_t),
2649 RT_FE(p),
2650 L2_FE(p)
2651{
2652 L2_FE.SetMapType(INTEGRAL);
2653
2654#ifndef MFEM_THREAD_SAFE
2655 rt_shape.SetSize(RT_FE.GetDof(), 2);
2656 l2_shape.SetSize(L2_FE.GetDof());
2657 rt_dshape.SetSize(RT_FE.GetDof());
2658#endif
2659
2660 int o = 0;
2661 int r = 0;
2662 int l = 0;
2663
2664 // Three edges
2665 for (int e=0; e<3; e++)
2666 {
2667 // Dofs in the plane
2668 for (int i=0; i<=p; i++)
2669 {
2670 dof_map[o] = r++; dof2nk[o++] = e;
2671 }
2672 }
2673
2674 // Interior dofs in the plane
2675 for (int j = 0; j < p; j++)
2676 for (int i = 0; i + j < p; i++)
2677 {
2678 dof_map[o] = r++; dof2nk[o++] = 0;
2679 dof_map[o] = r++; dof2nk[o++] = 2;
2680 }
2681
2682 // Interior z-directed dofs
2683 for (int j = 0; j <= p; j++)
2684 for (int i = 0; i + j <= p; i++)
2685 {
2686 dof_map[o] = -1 - l++; dof2nk[o++] = 3;
2687 }
2688
2689 MFEM_VERIFY(r == RT_FE.GetDof(),
2690 "RT_R2D_Triangle incorrect number of RT dofs.");
2691 MFEM_VERIFY(l == L2_FE.GetDof(),
2692 "RT_R2D_Triangle incorrect number of L2 dofs.");
2693 MFEM_VERIFY(o == GetDof(),
2694 "RT_R2D_Triangle incorrect number of dofs.");
2695
2696 const IntegrationRule & rt_Nodes = RT_FE.GetNodes();
2697 const IntegrationRule & l2_Nodes = L2_FE.GetNodes();
2698
2699 for (int i=0; i<dof; i++)
2700 {
2701 int idx = dof_map[i];
2702 if (idx >= 0)
2703 {
2704 const IntegrationPoint & ip = rt_Nodes.IntPoint(idx);
2705 Nodes.IntPoint(i).Set2(ip.x, ip.y);
2706 }
2707 else
2708 {
2709 const IntegrationPoint & ip = l2_Nodes.IntPoint(-idx-1);
2710 Nodes.IntPoint(i).Set2(ip.x, ip.y);
2711 }
2712 }
2713}
2714
2716 DenseMatrix &shape) const
2717{
2718#ifdef MFEM_THREAD_SAFE
2719 DenseMatrix rt_shape(RT_FE.GetDof(), 2);
2720 Vector l2_shape(L2_FE.GetDof());
2721#endif
2722
2723 RT_FE.CalcVShape(ip, rt_shape);
2724 L2_FE.CalcShape(ip, l2_shape);
2725
2726 for (int i=0; i<dof; i++)
2727 {
2728 int idx = dof_map[i];
2729 if (idx >= 0)
2730 {
2731 shape(i, 0) = rt_shape(idx, 0);
2732 shape(i, 1) = rt_shape(idx, 1);
2733 shape(i, 2) = 0.0;
2734 }
2735 else
2736 {
2737 shape(i, 0) = 0.0;
2738 shape(i, 1) = 0.0;
2739 shape(i, 2) = l2_shape(-idx-1);
2740 }
2741 }
2742}
2743
2745 Vector &div_shape) const
2746{
2747#ifdef MFEM_THREAD_SAFE
2748 Vector rt_dshape(RT_FE.GetDof());
2749#endif
2750
2751 RT_FE.CalcDivShape(ip, rt_dshape);
2752
2753 for (int i=0; i<dof; i++)
2754 {
2755 int idx = dof_map[i];
2756 if (idx >= 0)
2757 {
2758 div_shape(i) = rt_dshape(idx);
2759 }
2760 else
2761 {
2762 div_shape(i) = 0.0;
2763 }
2764 }
2765}
2766
2767const real_t RT_R2D_QuadrilateralElement::nk_q[15] =
2768{ 0., -1., 0., 1., 0., 0., 0., 1., 0., -1., 0., 0., 0., 0., 1. };
2769
2771 const int cb_type,
2772 const int ob_type)
2773 : RT_R2D_FiniteElement(p, Geometry::SQUARE, (3*p + 5)*(p + 1), nk_q),
2774 cbasis1d(poly1d.GetBasis(p + 1, VerifyClosed(cb_type))),
2775 obasis1d(poly1d.GetBasis(p, VerifyOpen(ob_type)))
2776{
2777 const real_t *cp = poly1d.ClosedPoints(p + 1, cb_type);
2778 const real_t *op = poly1d.OpenPoints(p, ob_type);
2779 const int dofx = (p + 1)*(p + 2);
2780 const int dofy = (p + 1)*(p + 2);
2781 const int dofxy = dofx + dofy;
2782
2783#ifndef MFEM_THREAD_SAFE
2784 shape_cx.SetSize(p + 2);
2785 shape_ox.SetSize(p + 1);
2786 shape_cy.SetSize(p + 2);
2787 shape_oy.SetSize(p + 1);
2788 dshape_cx.SetSize(p + 2);
2789 dshape_cy.SetSize(p + 2);
2790#endif
2791
2792 // edges
2793 int o = 0;
2794 for (int i = 0; i <= p; i++) // (0,1)
2795 {
2796 dof_map[dofx + i + 0*(p + 1)] = o++;
2797 }
2798 for (int i = 0; i <= p; i++) // (1,2)
2799 {
2800 dof_map[(p + 1) + i*(p + 2)] = o++;
2801 }
2802 for (int i = 0; i <= p; i++) // (2,3)
2803 {
2804 dof_map[dofx + (p - i) + (p + 1)*(p + 1)] = o++;
2805 }
2806 for (int i = 0; i <= p; i++) // (3,0)
2807 {
2808 dof_map[0 + (p - i)*(p + 2)] = o++;
2809 }
2810
2811 // interior
2812 for (int j = 0; j <= p; j++) // x-components
2813 for (int i = 1; i <= p; i++)
2814 {
2815 dof_map[i + j*(p + 2)] = o++;
2816 }
2817 for (int j = 1; j <= p; j++) // y-components
2818 for (int i = 0; i <= p; i++)
2819 {
2820 dof_map[dofx + i + j*(p + 1)] = o++;
2821 }
2822 for (int j = 0; j <= p; j++) // z-components
2823 for (int i = 0; i <= p; i++)
2824 {
2825 dof_map[dofxy + i + j*(p + 1)] = o++;
2826 }
2827
2828 // dof orientations
2829 // x-components
2830 for (int j = 0; j <= p; j++)
2831 for (int i = 0; i <= p/2; i++)
2832 {
2833 int idx = i + j*(p + 2);
2834 dof_map[idx] = -1 - dof_map[idx];
2835 }
2836 if (p%2 == 1)
2837 for (int j = p/2 + 1; j <= p; j++)
2838 {
2839 int idx = (p/2 + 1) + j*(p + 2);
2840 dof_map[idx] = -1 - dof_map[idx];
2841 }
2842 // y-components
2843 for (int j = 0; j <= p/2; j++)
2844 for (int i = 0; i <= p; i++)
2845 {
2846 int idx = dofx + i + j*(p + 1);
2847 dof_map[idx] = -1 - dof_map[idx];
2848 }
2849 if (p%2 == 1)
2850 for (int i = 0; i <= p/2; i++)
2851 {
2852 int idx = dofx + i + (p/2 + 1)*(p + 1);
2853 dof_map[idx] = -1 - dof_map[idx];
2854 }
2855
2856 o = 0;
2857 for (int j = 0; j <= p; j++)
2858 for (int i = 0; i <= p + 1; i++)
2859 {
2860 int idx;
2861 if ((idx = dof_map[o++]) < 0)
2862 {
2863 idx = -1 - idx;
2864 dof2nk[idx] = 3;
2865 }
2866 else
2867 {
2868 dof2nk[idx] = 1;
2869 }
2870 Nodes.IntPoint(idx).Set2(cp[i], op[j]);
2871 }
2872 for (int j = 0; j <= p + 1; j++)
2873 for (int i = 0; i <= p; i++)
2874 {
2875 int idx;
2876 if ((idx = dof_map[o++]) < 0)
2877 {
2878 idx = -1 - idx;
2879 dof2nk[idx] = 0;
2880 }
2881 else
2882 {
2883 dof2nk[idx] = 2;
2884 }
2885 Nodes.IntPoint(idx).Set2(op[i], cp[j]);
2886 }
2887 for (int j = 0; j <= p; j++)
2888 for (int i = 0; i <= p; i++)
2889 {
2890 int idx = dof_map[o++];
2891 dof2nk[idx] = 4;
2892 Nodes.IntPoint(idx).Set2(op[i], op[j]);
2893 }
2894}
2895
2897 DenseMatrix &shape) const
2898{
2899 const int pp1 = order;
2900
2901#ifdef MFEM_THREAD_SAFE
2902 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2903#endif
2904
2905 cbasis1d.Eval(ip.x, shape_cx);
2906 obasis1d.Eval(ip.x, shape_ox);
2907 cbasis1d.Eval(ip.y, shape_cy);
2908 obasis1d.Eval(ip.y, shape_oy);
2909
2910 int o = 0;
2911 for (int j = 0; j < pp1; j++)
2912 for (int i = 0; i <= pp1; i++)
2913 {
2914 int idx, s;
2915 if ((idx = dof_map[o++]) < 0)
2916 {
2917 idx = -1 - idx, s = -1;
2918 }
2919 else
2920 {
2921 s = +1;
2922 }
2923 shape(idx,0) = s*shape_cx(i)*shape_oy(j);
2924 shape(idx,1) = 0.;
2925 shape(idx,2) = 0.;
2926 }
2927 for (int j = 0; j <= pp1; j++)
2928 for (int i = 0; i < pp1; i++)
2929 {
2930 int idx, s;
2931 if ((idx = dof_map[o++]) < 0)
2932 {
2933 idx = -1 - idx, s = -1;
2934 }
2935 else
2936 {
2937 s = +1;
2938 }
2939 shape(idx,0) = 0.;
2940 shape(idx,1) = s*shape_ox(i)*shape_cy(j);
2941 shape(idx,2) = 0.;
2942 }
2943 for (int j = 0; j < pp1; j++)
2944 for (int i = 0; i < pp1; i++)
2945 {
2946 int idx = dof_map[o++];
2947 shape(idx,0) = 0.;
2948 shape(idx,1) = 0.;
2949 shape(idx,2) = shape_ox(i)*shape_oy(j);
2950 }
2951}
2952
2954 Vector &divshape) const
2955{
2956 const int pp1 = order;
2957
2958#ifdef MFEM_THREAD_SAFE
2959 Vector shape_cx(pp1 + 1), shape_ox(pp1), shape_cy(pp1 + 1), shape_oy(pp1);
2960 Vector dshape_cx(pp1 + 1), dshape_cy(pp1 + 1);
2961#endif
2962
2963 cbasis1d.Eval(ip.x, shape_cx, dshape_cx);
2964 obasis1d.Eval(ip.x, shape_ox);
2965 cbasis1d.Eval(ip.y, shape_cy, dshape_cy);
2966 obasis1d.Eval(ip.y, shape_oy);
2967
2968 int o = 0;
2969 for (int j = 0; j < pp1; j++)
2970 for (int i = 0; i <= pp1; i++)
2971 {
2972 int idx, s;
2973 if ((idx = dof_map[o++]) < 0)
2974 {
2975 idx = -1 - idx, s = -1;
2976 }
2977 else
2978 {
2979 s = +1;
2980 }
2981 divshape(idx) = s*dshape_cx(i)*shape_oy(j);
2982 }
2983 for (int j = 0; j <= pp1; j++)
2984 for (int i = 0; i < pp1; i++)
2985 {
2986 int idx, s;
2987 if ((idx = dof_map[o++]) < 0)
2988 {
2989 idx = -1 - idx, s = -1;
2990 }
2991 else
2992 {
2993 s = +1;
2994 }
2995 divshape(idx) = s*shape_ox(i)*dshape_cy(j);
2996 }
2997 for (int j = 0; j < pp1; j++)
2998 for (int i = 0; i < pp1; i++)
2999 {
3000 int idx = dof_map[o++];
3001 divshape(idx) = 0.;
3002 }
3003}
3004
3005}
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
void Factor()
Factor the current DenseMatrix, *a.
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication with the inverse of dense matrix.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication.
Definition densemat.cpp:125
void MultTranspose(const real_t *x, real_t *y) const
Multiply a vector with the transpose matrix.
Definition densemat.cpp:172
void Threshold(real_t eps)
Replace small entries, abs(a_ij) <= eps, with zero.
void SetRow(int r, const real_t *row)
real_t InnerProduct(const real_t *x, const real_t *y) const
Compute y^t A x.
Definition densemat.cpp:301
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:108
void GetColumn(int c, Vector &col) const
Rank 3 tensor (array of matrices)
void SetSize(int i, int j, int k, MemoryType mt_=MemoryType::PRESERVE)
const DenseMatrix & AdjugateJacobian()
Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoin...
Definition eltrans.hpp:148
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition eltrans.hpp:111
real_t Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint....
Definition eltrans.hpp:144
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint,...
Definition eltrans.hpp:132
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition eltrans.hpp:106
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Abstract class for all finite elements.
Definition fe_base.hpp:244
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_base.cpp:40
int dof
Number of degrees of freedom.
Definition fe_base.hpp:253
int GetRangeDim() const
Returns the vector dimension for vector-valued finite elements, which is also the dimension of the in...
Definition fe_base.hpp:325
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in reference space at the given p...
Definition fe_base.cpp:65
IntegrationRule Nodes
Definition fe_base.hpp:256
int vdim
Vector dimension of vector-valued basis functions.
Definition fe_base.hpp:247
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition fe_base.hpp:360
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition fe_base.hpp:351
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:400
Geometry::Type geom_type
Geometry::Type of the reference element.
Definition fe_base.hpp:249
DenseMatrix vshape
Definition fe_base.hpp:258
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:334
int order
Order/degree of the shape functions.
Definition fe_base.hpp:254
int dim
Dimension of reference space.
Definition fe_base.hpp:246
static DenseMatrix grad_mu01(real_t z)
void phi_Q(int p, Vector s, Vector t, DenseMatrix &u) const
void V_R(int p, Vector s, const DenseMatrix &grad_s, real_t mu, Vector grad_mu, real_t t, Vector grad_t, DenseMatrix &u) const
static constexpr real_t apex_tol
static void phi_E(int p, real_t s0, real_t s1, real_t *u)
static Vector nu012(real_t z, Vector xy, unsigned int ab)
static Vector nu01_grad_nu01(real_t z, Vector xy, unsigned int ab)
void VT_T(int p, Vector s, Vector sds, Vector sdsxds, real_t mu, Vector grad_mu, DenseTensor &u) const
static Vector mu01(real_t z)
void V_L(int p, Vector sx, const DenseMatrix &grad_sx, Vector sy, const DenseMatrix &grad_sy, real_t t, Vector grad_t, DenseTensor &u) const
static Vector grad_nu2(real_t z, const Vector xy, unsigned int ab)
static Vector grad_mu1(real_t z)
void V_Q(int p, Vector s, Vector ds, Vector t, Vector dt, DenseTensor &u) const
static real_t mu0(real_t z)
void V_T(int p, Vector s, Vector sdsxds, DenseTensor &u) const
static Vector grad_mu0(real_t z)
void E_Q(int p, Vector s, Vector ds, Vector t, DenseTensor &u) const
static Vector mu01_grad_mu01(real_t z, Vector xy, unsigned int ab)
static real_t mu1(real_t z)
static Vector nu012_grad_nu012(real_t z, Vector xy, unsigned int ab)
Describes the function space on each element.
Definition fe_base.hpp:226
static const int MaxDim
Definition geom.hpp:47
const IntegrationPoint & GetCenter(int GeomType) const
Return the center of the given Geometry::Type, GeomType.
Definition geom.hpp:75
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
Definition geom.cpp:435
void CalcShape(const IntegrationPoint &ip, Vector &shape) const override
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_h1.cpp:40
void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const override
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Definition fe_h1.cpp:59
Class for integration point with weight.
Definition intrules.hpp:35
void Set(const real_t *p, const int dim)
Definition intrules.hpp:46
void Set2(const real_t x1, const real_t x2)
Definition intrules.hpp:89
void Set3(const real_t x1, const real_t x2, const real_t x3)
Definition intrules.hpp:79
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void CalcShape(const IntegrationPoint &ip, Vector &shape) const override
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_l2.cpp:39
void CalcShape(const IntegrationPoint &ip, Vector &shape) const override
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
Definition fe_l2.cpp:615
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
void Eval(const real_t x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
Definition fe_base.cpp:1777
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
Definition fe_base.hpp:1047
void ScaleIntegrated(bool scale_integrated_)
Set whether the "integrated" basis should be scaled by the subcell sizes. Has no effect for non-integ...
Definition fe_base.cpp:2034
void EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the "integrated" basis type using pre-computed closed basis derivatives.
Definition fe_base.cpp:2009
const real_t * ClosedPoints(const int p, const int btype=BasisType::GaussLobatto, bool on_device=false)
Get coordinates of a closed (GaussLegendre) set of points if degree p.
Definition fe_base.hpp:1121
const real_t * OpenPoints(const int p, const int btype=BasisType::GaussLegendre, bool on_device=false)
Get coordinates of an open (GaussLegendre) set of points if degree p.
Definition fe_base.hpp:1113
static void CalcBasis(const int p, const real_t x, real_t *u)
Evaluate the values of a hierarchical 1D basis at point x hierarchical = k-th basis function is degre...
Definition fe_base.hpp:1140
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:1398
RT_FuentesPyramidElement(const int p)
Definition fe_rt.cpp:1273
void CalcRawDivShape(const IntegrationPoint &ip, Vector &dshape) const
Definition fe_rt.cpp:1478
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:1449
void CalcRawVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition fe_rt.cpp:1425
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:492
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition fe_rt.cpp:661
void GetFaceMap(const int face_id, Array< int > &face_map) const override
Return the mapping from lexicographically ordered face DOFs to lexicographically ordered element DOFs...
Definition fe_rt.cpp:711
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:582
RT_HexahedronElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_HexahedronElement of order p and closed and open BasisType cb_type and ob_type.
Definition fe_rt.cpp:326
void ProjectIntegrated(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition fe_rt.cpp:258
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:203
void GetFaceMap(const int face_id, Array< int > &face_map) const override
Return the mapping from lexicographic face DOFs to lexicographic element DOFs for the given local fac...
Definition fe_rt.cpp:301
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:142
RT_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_type.
Definition fe_rt.cpp:26
void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const override
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
Definition fe_rt.cpp:2146
RT_R1D_SegmentElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_R1D_SegmentElement of order p and closed and open BasisType cb_type and ob_type.
Definition fe_rt.cpp:2005
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:2059
void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const override
Compute the discrete curl matrix from the given FiniteElement onto 'this' FiniteElement....
Definition fe_rt.cpp:2258
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:2113
void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const override
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs.
Definition fe_rt.cpp:2463
void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const override
Given a vector coefficient and a transformation, compute its projection (approximation) in the local ...
Definition fe_rt.cpp:2511
void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const override
Compute the discrete curl matrix from the given FiniteElement onto 'this' FiniteElement....
Definition fe_rt.cpp:2623
void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in physical space at the point ...
Definition fe_rt.cpp:2404
RT_R2D_FiniteElement(int p, Geometry::Type G, int Do, const real_t *nk_fe)
Definition fe_rt.cpp:2392
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:2896
RT_R2D_QuadrilateralElement(const int p, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Construct the RT_QuadrilateralElement of order p and closed and open BasisType cb_type and ob_type.
Definition fe_rt.cpp:2770
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:2953
void CalcDivShape(const IntegrationPoint &ip, Vector &div_shape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:2345
RT_R2D_SegmentElement(const int p, const int ob_type=BasisType::GaussLegendre)
Construct the RT_R2D_SegmentElement of order p and open BasisType ob_type.
Definition fe_rt.cpp:2281
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:2309
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:2715
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:2744
RT_R2D_TriangleElement(const int p)
Construct the RT_R2D_TriangleElement of order p.
Definition fe_rt.cpp:2647
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:1037
RT_TetrahedronElement(const int p)
Construct the RT_TetrahedronElement of order p.
Definition fe_rt.cpp:899
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:1001
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:824
RT_TriangleElement(const int p)
Construct the RT_TriangleElement of order p.
Definition fe_rt.cpp:746
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:857
void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const override
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition fe_rt.cpp:1235
RT_WedgeElement(const int p)
Definition fe_rt.cpp:1082
void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const override
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition fe_rt.cpp:1200
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition fe_base.hpp:687
Poly_1D::Basis & basis1d
Definition fe_base.hpp:1251
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
Intermediate class for finite elements whose basis functions return vector values.
Definition fe_base.hpp:813
Vector data type.
Definition vector.hpp:82
int Size() const
Returns the size of the vector.
Definition vector.hpp:226
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
real_t mu
Definition ex25.cpp:140
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition vector.hpp:47
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
Geometry Geometries
Definition fe.cpp:49
float real_t
Definition config.hpp:43
Poly_1D poly1d
Definition fe.cpp:28
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
Definition fe_base.cpp:748
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492
real_t p(const Vector &x, real_t t)