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