MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
fe_pos.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// H1 Finite Element classes utilizing the Bernstein basis
13
14#include "fe_pos.hpp"
15#include "face_map_utils.hpp"
16#include "../bilininteg.hpp"
17#include "../lininteg.hpp"
18#include "../coefficient.hpp"
19
20namespace mfem
21{
22
23using namespace std;
24
26 Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
27{
28 for (int i = 0; i < dof; i++)
29 {
30 const IntegrationPoint &ip = Nodes.IntPoint(i);
31 Trans.SetIntPoint(&ip);
32 dofs(i) = coeff.Eval(Trans, ip);
33 }
34}
35
37 VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
38{
39 MFEM_ASSERT(dofs.Size() == vc.GetVDim()*dof, "");
40 Vector x(vc.GetVDim());
41
42 for (int i = 0; i < dof; i++)
43 {
44 const IntegrationPoint &ip = Nodes.IntPoint(i);
45 Trans.SetIntPoint(&ip);
46 vc.Eval (x, Trans, ip);
47 for (int j = 0; j < x.Size(); j++)
48 {
49 dofs(dof*j+i) = x(j);
50 }
51 }
52}
53
55 const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
56{
57 const NodalFiniteElement *nfe =
58 dynamic_cast<const NodalFiniteElement *>(&fe);
59
60 if (nfe && dof == nfe->GetDof())
61 {
62 nfe->Project(*this, Trans, I);
63 I.Invert();
64 }
65 else
66 {
67 // local L2 projection
68 DenseMatrix pos_mass, mixed_mass;
69 MassIntegrator mass_integ;
70
71 mass_integ.AssembleElementMatrix(*this, Trans, pos_mass);
72 mass_integ.AssembleElementMatrix2(fe, *this, Trans, mixed_mass);
73
74 DenseMatrixInverse pos_mass_inv(pos_mass);
75 I.SetSize(dof, fe.GetDof());
76 pos_mass_inv.Mult(mixed_mass, I);
77 }
78}
79
80
82 const int dims, const int p, const DofMapType dmtype)
83 : PositiveFiniteElement(dims, GetTensorProductGeometry(dims),
84 Pow(p + 1, dims), p,
85 dims > 1 ? FunctionSpace::Qk : FunctionSpace::Pk),
86 TensorBasisElement(dims, p, BasisType::Positive, dmtype) { }
87
89 Array<int> &face_map) const
90{
91 internal::GetTensorFaceMap(dim, order, face_id, face_map);
92}
93
94
96 : PositiveFiniteElement(2, Geometry::SQUARE, 9, 2, FunctionSpace::Qk)
97{
98 Nodes.IntPoint(0).x = 0.0;
99 Nodes.IntPoint(0).y = 0.0;
100 Nodes.IntPoint(1).x = 1.0;
101 Nodes.IntPoint(1).y = 0.0;
102 Nodes.IntPoint(2).x = 1.0;
103 Nodes.IntPoint(2).y = 1.0;
104 Nodes.IntPoint(3).x = 0.0;
105 Nodes.IntPoint(3).y = 1.0;
106 Nodes.IntPoint(4).x = 0.5;
107 Nodes.IntPoint(4).y = 0.0;
108 Nodes.IntPoint(5).x = 1.0;
109 Nodes.IntPoint(5).y = 0.5;
110 Nodes.IntPoint(6).x = 0.5;
111 Nodes.IntPoint(6).y = 1.0;
112 Nodes.IntPoint(7).x = 0.0;
113 Nodes.IntPoint(7).y = 0.5;
114 Nodes.IntPoint(8).x = 0.5;
115 Nodes.IntPoint(8).y = 0.5;
116}
117
119 Vector &shape) const
120{
121 real_t x = ip.x, y = ip.y;
122 real_t l1x, l2x, l3x, l1y, l2y, l3y;
123
124 l1x = (1. - x) * (1. - x);
125 l2x = 2. * x * (1. - x);
126 l3x = x * x;
127 l1y = (1. - y) * (1. - y);
128 l2y = 2. * y * (1. - y);
129 l3y = y * y;
130
131 shape(0) = l1x * l1y;
132 shape(4) = l2x * l1y;
133 shape(1) = l3x * l1y;
134 shape(7) = l1x * l2y;
135 shape(8) = l2x * l2y;
136 shape(5) = l3x * l2y;
137 shape(3) = l1x * l3y;
138 shape(6) = l2x * l3y;
139 shape(2) = l3x * l3y;
140}
141
143 DenseMatrix &dshape) const
144{
145 real_t x = ip.x, y = ip.y;
146 real_t l1x, l2x, l3x, l1y, l2y, l3y;
147 real_t d1x, d2x, d3x, d1y, d2y, d3y;
148
149 l1x = (1. - x) * (1. - x);
150 l2x = 2. * x * (1. - x);
151 l3x = x * x;
152 l1y = (1. - y) * (1. - y);
153 l2y = 2. * y * (1. - y);
154 l3y = y * y;
155
156 d1x = 2. * x - 2.;
157 d2x = 2. - 4. * x;
158 d3x = 2. * x;
159 d1y = 2. * y - 2.;
160 d2y = 2. - 4. * y;
161 d3y = 2. * y;
162
163 dshape(0,0) = d1x * l1y;
164 dshape(0,1) = l1x * d1y;
165
166 dshape(4,0) = d2x * l1y;
167 dshape(4,1) = l2x * d1y;
168
169 dshape(1,0) = d3x * l1y;
170 dshape(1,1) = l3x * d1y;
171
172 dshape(7,0) = d1x * l2y;
173 dshape(7,1) = l1x * d2y;
174
175 dshape(8,0) = d2x * l2y;
176 dshape(8,1) = l2x * d2y;
177
178 dshape(5,0) = d3x * l2y;
179 dshape(5,1) = l3x * d2y;
180
181 dshape(3,0) = d1x * l3y;
182 dshape(3,1) = l1x * d3y;
183
184 dshape(6,0) = d2x * l3y;
185 dshape(6,1) = l2x * d3y;
186
187 dshape(2,0) = d3x * l3y;
188 dshape(2,1) = l3x * d3y;
189}
190
192 ElementTransformation &Trans, DenseMatrix &I) const
193{
194 real_t s[9];
195 IntegrationPoint tr_ip;
196 Vector xx(&tr_ip.x, 2), shape(s, 9);
197
198 for (int i = 0; i < 9; i++)
199 {
200 Trans.Transform(Nodes.IntPoint(i), xx);
201 CalcShape(tr_ip, shape);
202 for (int j = 0; j < 9; j++)
203 if (fabs(I(i,j) = s[j]) < 1.0e-12)
204 {
205 I(i,j) = 0.0;
206 }
207 }
208 for (int i = 0; i < 9; i++)
209 {
210 real_t *d = &I(0,i);
211 d[4] = 2. * d[4] - 0.5 * (d[0] + d[1]);
212 d[5] = 2. * d[5] - 0.5 * (d[1] + d[2]);
213 d[6] = 2. * d[6] - 0.5 * (d[2] + d[3]);
214 d[7] = 2. * d[7] - 0.5 * (d[3] + d[0]);
215 d[8] = 4. * d[8] - 0.5 * (d[4] + d[5] + d[6] + d[7]) -
216 0.25 * (d[0] + d[1] + d[2] + d[3]);
217 }
218}
219
221 Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
222{
223 real_t *d = dofs.GetData();
224
225 for (int i = 0; i < 9; i++)
226 {
227 const IntegrationPoint &ip = Nodes.IntPoint(i);
228 Trans.SetIntPoint(&ip);
229 d[i] = coeff.Eval(Trans, ip);
230 }
231 d[4] = 2. * d[4] - 0.5 * (d[0] + d[1]);
232 d[5] = 2. * d[5] - 0.5 * (d[1] + d[2]);
233 d[6] = 2. * d[6] - 0.5 * (d[2] + d[3]);
234 d[7] = 2. * d[7] - 0.5 * (d[3] + d[0]);
235 d[8] = 4. * d[8] - 0.5 * (d[4] + d[5] + d[6] + d[7]) -
236 0.25 * (d[0] + d[1] + d[2] + d[3]);
237}
238
241 Vector &dofs) const
242{
243 real_t v[3];
244 Vector x (v, vc.GetVDim());
245
246 for (int i = 0; i < 9; i++)
247 {
248 const IntegrationPoint &ip = Nodes.IntPoint(i);
249 Trans.SetIntPoint(&ip);
250 vc.Eval (x, Trans, ip);
251 for (int j = 0; j < x.Size(); j++)
252 {
253 dofs(9*j+i) = v[j];
254 }
255 }
256 for (int j = 0; j < x.Size(); j++)
257 {
258 real_t *d = &dofs(9*j);
259
260 d[4] = 2. * d[4] - 0.5 * (d[0] + d[1]);
261 d[5] = 2. * d[5] - 0.5 * (d[1] + d[2]);
262 d[6] = 2. * d[6] - 0.5 * (d[2] + d[3]);
263 d[7] = 2. * d[7] - 0.5 * (d[3] + d[0]);
264 d[8] = 4. * d[8] - 0.5 * (d[4] + d[5] + d[6] + d[7]) -
265 0.25 * (d[0] + d[1] + d[2] + d[3]);
266 }
267}
268
269
271 : PositiveFiniteElement(1, Geometry::SEGMENT, 3, 2)
272{
273 Nodes.IntPoint(0).x = 0.0;
274 Nodes.IntPoint(1).x = 1.0;
275 Nodes.IntPoint(2).x = 0.5;
276}
277
279 Vector &shape) const
280{
281 const real_t x = ip.x, x1 = 1. - x;
282
283 shape(0) = x1 * x1;
284 shape(1) = x * x;
285 shape(2) = 2. * x * x1;
286}
287
289 DenseMatrix &dshape) const
290{
291 const real_t x = ip.x;
292
293 dshape(0,0) = 2. * x - 2.;
294 dshape(1,0) = 2. * x;
295 dshape(2,0) = 2. - 4. * x;
296}
297
298
300 : PositiveTensorFiniteElement(1, p, H1_DOF_MAP)
301{
302#ifndef MFEM_THREAD_SAFE
303 // thread private versions; see class header.
304 shape_x.SetSize(p+1);
305 dshape_x.SetSize(p+1);
306#endif
307
308 // Endpoints need to be first in the list, so reorder them.
309 Nodes.IntPoint(0).x = 0.0;
310 Nodes.IntPoint(1).x = 1.0;
311 for (int i = 1; i < p; i++)
312 {
313 Nodes.IntPoint(i+1).x = real_t(i)/p;
314 }
315}
316
318 Vector &shape) const
319{
320 const int p = order;
321
322#ifdef MFEM_THREAD_SAFE
323 Vector shape_x(p+1);
324#endif
325
326 Poly_1D::CalcBernstein(p, ip.x, shape_x.GetData() );
327
328 // Endpoints need to be first in the list, so reorder them.
329 shape(0) = shape_x(0);
330 shape(1) = shape_x(p);
331 for (int i = 1; i < p; i++)
332 {
333 shape(i+1) = shape_x(i);
334 }
335}
336
338 DenseMatrix &dshape) const
339{
340 const int p = order;
341
342#ifdef MFEM_THREAD_SAFE
343 Vector shape_x(p+1), dshape_x(p+1);
344#endif
345
346 Poly_1D::CalcBernstein(p, ip.x, shape_x.GetData(), dshape_x.GetData() );
347
348 // Endpoints need to be first in the list, so reorder them.
349 dshape(0,0) = dshape_x(0);
350 dshape(1,0) = dshape_x(p);
351 for (int i = 1; i < p; i++)
352 {
353 dshape(i+1,0) = dshape_x(i);
354 }
355}
356
357void H1Pos_SegmentElement::ProjectDelta(int vertex, Vector &dofs) const
358{
359 dofs = 0.0;
360 dofs[vertex] = 1.0;
361}
362
363
365 : PositiveTensorFiniteElement(2, p, H1_DOF_MAP)
366{
367#ifndef MFEM_THREAD_SAFE
368 const int p1 = p + 1;
369
370 shape_x.SetSize(p1);
371 shape_y.SetSize(p1);
372 dshape_x.SetSize(p1);
373 dshape_y.SetSize(p1);
374#endif
375
376 int o = 0;
377 for (int j = 0; j <= p; j++)
378 for (int i = 0; i <= p; i++)
379 {
380 Nodes.IntPoint(dof_map[o++]).Set2(real_t(i)/p, real_t(j)/p);
381 }
382}
383
385 Vector &shape) const
386{
387 const int p = order;
388
389#ifdef MFEM_THREAD_SAFE
390 Vector shape_x(p+1), shape_y(p+1);
391#endif
392
393 Poly_1D::CalcBernstein(p, ip.x, shape_x);
394 Poly_1D::CalcBernstein(p, ip.y, shape_y);
395
396 // Reorder so that vertices are at the beginning of the list
397 for (int o = 0, j = 0; j <= p; j++)
398 for (int i = 0; i <= p; i++)
399 {
400 shape(dof_map[o++]) = shape_x(i)*shape_y(j);
401 }
402}
403
405 DenseMatrix &dshape) const
406{
407 const int p = order;
408
409#ifdef MFEM_THREAD_SAFE
410 Vector shape_x(p+1), shape_y(p+1), dshape_x(p+1), dshape_y(p+1);
411#endif
412
413 Poly_1D::CalcBernstein(p, ip.x, shape_x, dshape_x);
414 Poly_1D::CalcBernstein(p, ip.y, shape_y, dshape_y);
415
416 // Reorder so that vertices are at the beginning of the list
417 for (int o = 0, j = 0; j <= p; j++)
418 for (int i = 0; i <= p; i++)
419 {
420 dshape(dof_map[o],0) = dshape_x(i)* shape_y(j);
421 dshape(dof_map[o],1) = shape_x(i)*dshape_y(j); o++;
422 }
423}
424
426{
427 dofs = 0.0;
428 dofs[vertex] = 1.0;
429}
430
431
433 : PositiveTensorFiniteElement(3, p, H1_DOF_MAP)
434{
435#ifndef MFEM_THREAD_SAFE
436 const int p1 = p + 1;
437
438 shape_x.SetSize(p1);
439 shape_y.SetSize(p1);
440 shape_z.SetSize(p1);
441 dshape_x.SetSize(p1);
442 dshape_y.SetSize(p1);
443 dshape_z.SetSize(p1);
444#endif
445
446 int o = 0;
447 for (int k = 0; k <= p; k++)
448 for (int j = 0; j <= p; j++)
449 for (int i = 0; i <= p; i++)
450 Nodes.IntPoint(dof_map[o++]).Set3(real_t(i)/p, real_t(j)/p,
451 real_t(k)/p);
452}
453
455 Vector &shape) const
456{
457 const int p = order;
458
459#ifdef MFEM_THREAD_SAFE
460 Vector shape_x(p+1), shape_y(p+1), shape_z(p+1);
461#endif
462
463 Poly_1D::CalcBernstein(p, ip.x, shape_x.GetData() );
464 Poly_1D::CalcBernstein(p, ip.y, shape_y.GetData() );
465 Poly_1D::CalcBernstein(p, ip.z, shape_z.GetData() );
466
467 for (int o = 0, k = 0; k <= p; k++)
468 for (int j = 0; j <= p; j++)
469 for (int i = 0; i <= p; i++)
470 {
471 shape(dof_map[o++]) = shape_x(i)*shape_y(j)*shape_z(k);
472 }
473}
474
476 DenseMatrix &dshape) const
477{
478 const int p = order;
479
480#ifdef MFEM_THREAD_SAFE
481 Vector shape_x(p+1), shape_y(p+1), shape_z(p+1);
482 Vector dshape_x(p+1), dshape_y(p+1), dshape_z(p+1);
483#endif
484
485 Poly_1D::CalcBernstein(p, ip.x, shape_x.GetData(), dshape_x.GetData() );
486 Poly_1D::CalcBernstein(p, ip.y, shape_y.GetData(), dshape_y.GetData() );
487 Poly_1D::CalcBernstein(p, ip.z, shape_z.GetData(), dshape_z.GetData() );
488
489 for (int o = 0, k = 0; k <= p; k++)
490 for (int j = 0; j <= p; j++)
491 for (int i = 0; i <= p; i++)
492 {
493 dshape(dof_map[o],0) = dshape_x(i)* shape_y(j)* shape_z(k);
494 dshape(dof_map[o],1) = shape_x(i)*dshape_y(j)* shape_z(k);
495 dshape(dof_map[o],2) = shape_x(i)* shape_y(j)*dshape_z(k); o++;
496 }
497}
498
500{
501 dofs = 0.0;
502 dofs[vertex] = 1.0;
503}
504
505
507 : PositiveFiniteElement(2, Geometry::TRIANGLE, ((p + 1)*(p + 2))/2, p,
508 FunctionSpace::Pk)
509{
510#ifndef MFEM_THREAD_SAFE
512 dshape_1d.SetSize(p + 1);
514#endif
516
517 struct Index
518 {
519 int p2p3;
520 Index(int p) { p2p3 = 2*p + 3; }
521 int operator()(int i, int j) { return ((p2p3-j)*j)/2+i; }
522 };
523 Index idx(p);
524
525 // vertices
526 dof_map[idx(0,0)] = 0;
527 Nodes.IntPoint(0).Set2(0., 0.);
528 dof_map[idx(p,0)] = 1;
529 Nodes.IntPoint(1).Set2(1., 0.);
530 dof_map[idx(0,p)] = 2;
531 Nodes.IntPoint(2).Set2(0., 1.);
532
533 // edges
534 int o = 3;
535 for (int i = 1; i < p; i++)
536 {
537 dof_map[idx(i,0)] = o;
538 Nodes.IntPoint(o++).Set2(real_t(i)/p, 0.);
539 }
540 for (int i = 1; i < p; i++)
541 {
542 dof_map[idx(p-i,i)] = o;
543 Nodes.IntPoint(o++).Set2(real_t(p-i)/p, real_t(i)/p);
544 }
545 for (int i = 1; i < p; i++)
546 {
547 dof_map[idx(0,p-i)] = o;
548 Nodes.IntPoint(o++).Set2(0., real_t(p-i)/p);
549 }
550
551 // interior
552 for (int j = 1; j < p; j++)
553 for (int i = 1; i + j < p; i++)
554 {
555 dof_map[idx(i,j)] = o;
556 Nodes.IntPoint(o++).Set2(real_t(i)/p, real_t(j)/p);
557 }
558}
559
560// static method
562 const int p, const real_t l1, const real_t l2, real_t *shape)
563{
564 const real_t l3 = 1. - l1 - l2;
565
566 // The (i,j) basis function is given by: T(i,j,p-i-j) l1^i l2^j l3^{p-i-j},
567 // where T(i,j,k) = (i+j+k)! / (i! j! k!)
568 // Another expression is given by the terms of the expansion:
569 // (l1 + l2 + l3)^p =
570 // \sum_{j=0}^p \binom{p}{j} l2^j
571 // \sum_{i=0}^{p-j} \binom{p-j}{i} l1^i l3^{p-j-i}
572 const int *bp = Poly_1D::Binom(p);
573 real_t z = 1.;
574 for (int o = 0, j = 0; j <= p; j++)
575 {
576 Poly_1D::CalcBinomTerms(p - j, l1, l3, &shape[o]);
577 real_t s = bp[j]*z;
578 for (int i = 0; i <= p - j; i++)
579 {
580 shape[o++] *= s;
581 }
582 z *= l2;
583 }
584}
585
586// static method
588 const int p, const real_t l1, const real_t l2,
589 real_t *dshape_1d, real_t *dshape)
590{
591 const int dof = ((p + 1)*(p + 2))/2;
592 const real_t l3 = 1. - l1 - l2;
593
594 const int *bp = Poly_1D::Binom(p);
595 real_t z = 1.;
596 for (int o = 0, j = 0; j <= p; j++)
597 {
599 real_t s = bp[j]*z;
600 for (int i = 0; i <= p - j; i++)
601 {
602 dshape[o++] = s*dshape_1d[i];
603 }
604 z *= l2;
605 }
606 z = 1.;
607 for (int i = 0; i <= p; i++)
608 {
610 real_t s = bp[i]*z;
611 for (int o = i, j = 0; j <= p - i; j++)
612 {
613 dshape[dof + o] = s*dshape_1d[j];
614 o += p + 1 - j;
615 }
616 z *= l1;
617 }
618}
619
621 Vector &shape) const
622{
623#ifdef MFEM_THREAD_SAFE
625#endif
626 CalcShape(order, ip.x, ip.y, m_shape.GetData());
627 for (int i = 0; i < dof; i++)
628 {
629 shape(dof_map[i]) = m_shape(i);
630 }
631}
632
634 DenseMatrix &dshape) const
635{
636#ifdef MFEM_THREAD_SAFE
639#endif
641 for (int d = 0; d < 2; d++)
642 {
643 for (int i = 0; i < dof; i++)
644 {
645 dshape(dof_map[i],d) = m_dshape(i,d);
646 }
647 }
648}
649
650
652 : PositiveFiniteElement(3, Geometry::TETRAHEDRON,
653 ((p + 1)*(p + 2)*(p + 3))/6, p, FunctionSpace::Pk)
654{
655#ifndef MFEM_THREAD_SAFE
657 dshape_1d.SetSize(p + 1);
659#endif
661
662 struct Index
663 {
664 int p, dof;
665 int tri(int k) { return (k*(k + 1))/2; }
666 int tet(int k) { return (k*(k + 1)*(k + 2))/6; }
667 Index(int p_) { p = p_; dof = tet(p + 1); }
668 int operator()(int i, int j, int k)
669 { return dof - tet(p - k) - tri(p + 1 - k - j) + i; }
670 };
671 Index idx(p);
672
673 // vertices
674 dof_map[idx(0,0,0)] = 0;
675 Nodes.IntPoint(0).Set3(0., 0., 0.);
676 dof_map[idx(p,0,0)] = 1;
677 Nodes.IntPoint(1).Set3(1., 0., 0.);
678 dof_map[idx(0,p,0)] = 2;
679 Nodes.IntPoint(2).Set3(0., 1., 0.);
680 dof_map[idx(0,0,p)] = 3;
681 Nodes.IntPoint(3).Set3(0., 0., 1.);
682
683 // edges (see Tetrahedron::edges in mesh/tetrahedron.cpp)
684 int o = 4;
685 for (int i = 1; i < p; i++) // (0,1)
686 {
687 dof_map[idx(i,0,0)] = o;
688 Nodes.IntPoint(o++).Set3(real_t(i)/p, 0., 0.);
689 }
690 for (int i = 1; i < p; i++) // (0,2)
691 {
692 dof_map[idx(0,i,0)] = o;
693 Nodes.IntPoint(o++).Set3(0., real_t(i)/p, 0.);
694 }
695 for (int i = 1; i < p; i++) // (0,3)
696 {
697 dof_map[idx(0,0,i)] = o;
698 Nodes.IntPoint(o++).Set3(0., 0., real_t(i)/p);
699 }
700 for (int i = 1; i < p; i++) // (1,2)
701 {
702 dof_map[idx(p-i,i,0)] = o;
703 Nodes.IntPoint(o++).Set3(real_t(p-i)/p, real_t(i)/p, 0.);
704 }
705 for (int i = 1; i < p; i++) // (1,3)
706 {
707 dof_map[idx(p-i,0,i)] = o;
708 Nodes.IntPoint(o++).Set3(real_t(p-i)/p, 0., real_t(i)/p);
709 }
710 for (int i = 1; i < p; i++) // (2,3)
711 {
712 dof_map[idx(0,p-i,i)] = o;
713 Nodes.IntPoint(o++).Set3(0., real_t(p-i)/p, real_t(i)/p);
714 }
715
716 // faces (see Mesh::GenerateFaces in mesh/mesh.cpp)
717 for (int j = 1; j < p; j++)
718 for (int i = 1; i + j < p; i++) // (1,2,3)
719 {
720 dof_map[idx(p-i-j,i,j)] = o;
721 Nodes.IntPoint(o++).Set3(real_t(p-i-j)/p, real_t(i)/p, real_t(j)/p);
722 }
723 for (int j = 1; j < p; j++)
724 for (int i = 1; i + j < p; i++) // (0,3,2)
725 {
726 dof_map[idx(0,j,i)] = o;
727 Nodes.IntPoint(o++).Set3(0., real_t(j)/p, real_t(i)/p);
728 }
729 for (int j = 1; j < p; j++)
730 for (int i = 1; i + j < p; i++) // (0,1,3)
731 {
732 dof_map[idx(i,0,j)] = o;
733 Nodes.IntPoint(o++).Set3(real_t(i)/p, 0., real_t(j)/p);
734 }
735 for (int j = 1; j < p; j++)
736 for (int i = 1; i + j < p; i++) // (0,2,1)
737 {
738 dof_map[idx(j,i,0)] = o;
739 Nodes.IntPoint(o++).Set3(real_t(j)/p, real_t(i)/p, 0.);
740 }
741
742 // interior
743 for (int k = 1; k < p; k++)
744 for (int j = 1; j + k < p; j++)
745 for (int i = 1; i + j + k < p; i++)
746 {
747 dof_map[idx(i,j,k)] = o;
748 Nodes.IntPoint(o++).Set3(real_t(i)/p, real_t(j)/p, real_t(k)/p);
749 }
750}
751
752// static method
754 const int p, const real_t l1, const real_t l2, const real_t l3,
755 real_t *shape)
756{
757 const real_t l4 = 1. - l1 - l2 - l3;
758
759 // The basis functions are the terms in the expansion:
760 // (l1 + l2 + l3 + l4)^p =
761 // \sum_{k=0}^p \binom{p}{k} l3^k
762 // \sum_{j=0}^{p-k} \binom{p-k}{j} l2^j
763 // \sum_{i=0}^{p-k-j} \binom{p-k-j}{i} l1^i l4^{p-k-j-i}
764 const int *bp = Poly_1D::Binom(p);
765 real_t l3k = 1.;
766 for (int o = 0, k = 0; k <= p; k++)
767 {
768 const int *bpk = Poly_1D::Binom(p - k);
769 const real_t ek = bp[k]*l3k;
770 real_t l2j = 1.;
771 for (int j = 0; j <= p - k; j++)
772 {
773 Poly_1D::CalcBinomTerms(p - k - j, l1, l4, &shape[o]);
774 real_t ekj = ek*bpk[j]*l2j;
775 for (int i = 0; i <= p - k - j; i++)
776 {
777 shape[o++] *= ekj;
778 }
779 l2j *= l2;
780 }
781 l3k *= l3;
782 }
783}
784
785// static method
787 const int p, const real_t l1, const real_t l2, const real_t l3,
788 real_t *dshape_1d, real_t *dshape)
789{
790 const int dof = ((p + 1)*(p + 2)*(p + 3))/6;
791 const real_t l4 = 1. - l1 - l2 - l3;
792
793 // For the x derivatives, differentiate the terms of the expression:
794 // \sum_{k=0}^p \binom{p}{k} l3^k
795 // \sum_{j=0}^{p-k} \binom{p-k}{j} l2^j
796 // \sum_{i=0}^{p-k-j} \binom{p-k-j}{i} l1^i l4^{p-k-j-i}
797 const int *bp = Poly_1D::Binom(p);
798 real_t l3k = 1.;
799 for (int o = 0, k = 0; k <= p; k++)
800 {
801 const int *bpk = Poly_1D::Binom(p - k);
802 const real_t ek = bp[k]*l3k;
803 real_t l2j = 1.;
804 for (int j = 0; j <= p - k; j++)
805 {
806 Poly_1D::CalcDBinomTerms(p - k - j, l1, l4, dshape_1d);
807 real_t ekj = ek*bpk[j]*l2j;
808 for (int i = 0; i <= p - k - j; i++)
809 {
810 dshape[o++] = dshape_1d[i]*ekj;
811 }
812 l2j *= l2;
813 }
814 l3k *= l3;
815 }
816 // For the y derivatives, differentiate the terms of the expression:
817 // \sum_{k=0}^p \binom{p}{k} l3^k
818 // \sum_{i=0}^{p-k} \binom{p-k}{i} l1^i
819 // \sum_{j=0}^{p-k-i} \binom{p-k-i}{j} l2^j l4^{p-k-j-i}
820 l3k = 1.;
821 for (int ok = 0, k = 0; k <= p; k++)
822 {
823 const int *bpk = Poly_1D::Binom(p - k);
824 const real_t ek = bp[k]*l3k;
825 real_t l1i = 1.;
826 for (int i = 0; i <= p - k; i++)
827 {
828 Poly_1D::CalcDBinomTerms(p - k - i, l2, l4, dshape_1d);
829 real_t eki = ek*bpk[i]*l1i;
830 int o = ok + i;
831 for (int j = 0; j <= p - k - i; j++)
832 {
833 dshape[dof + o] = dshape_1d[j]*eki;
834 o += p - k - j + 1;
835 }
836 l1i *= l1;
837 }
838 l3k *= l3;
839 ok += ((p - k + 2)*(p - k + 1))/2;
840 }
841 // For the z derivatives, differentiate the terms of the expression:
842 // \sum_{j=0}^p \binom{p}{j} l2^j
843 // \sum_{i=0}^{p-j} \binom{p-j}{i} l1^i
844 // \sum_{k=0}^{p-j-i} \binom{p-j-i}{k} l3^k l4^{p-k-j-i}
845 real_t l2j = 1.;
846 for (int j = 0; j <= p; j++)
847 {
848 const int *bpj = Poly_1D::Binom(p - j);
849 const real_t ej = bp[j]*l2j;
850 real_t l1i = 1.;
851 for (int i = 0; i <= p - j; i++)
852 {
853 Poly_1D::CalcDBinomTerms(p - j - i, l3, l4, dshape_1d);
854 real_t eji = ej*bpj[i]*l1i;
855 int m = ((p + 2)*(p + 1))/2;
856 int n = ((p - j + 2)*(p - j + 1))/2;
857 for (int o = i, k = 0; k <= p - j - i; k++)
858 {
859 // m = ((p - k + 2)*(p - k + 1))/2;
860 // n = ((p - k - j + 2)*(p - k - j + 1))/2;
861 o += m;
862 dshape[2*dof + o - n] = dshape_1d[k]*eji;
863 m -= p - k + 1;
864 n -= p - k - j + 1;
865 }
866 l1i *= l1;
867 }
868 l2j *= l2;
869 }
870}
871
873 Vector &shape) const
874{
875#ifdef MFEM_THREAD_SAFE
877#endif
878 CalcShape(order, ip.x, ip.y, ip.z, m_shape.GetData());
879 for (int i = 0; i < dof; i++)
880 {
881 shape(dof_map[i]) = m_shape(i);
882 }
883}
884
886 DenseMatrix &dshape) const
887{
888#ifdef MFEM_THREAD_SAFE
891#endif
892 CalcDShape(order, ip.x, ip.y, ip.z, dshape_1d.GetData(), m_dshape.Data());
893 for (int d = 0; d < 3; d++)
894 {
895 for (int i = 0; i < dof; i++)
896 {
897 dshape(dof_map[i],d) = m_dshape(i,d);
898 }
899 }
900}
901
902
905 ((p + 1)*(p + 1)*(p + 2))/2, p, FunctionSpace::Qk),
906 TriangleFE(p),
907 SegmentFE(p)
908{
909#ifndef MFEM_THREAD_SAFE
914#endif
915
918
919 // Nodal DoFs
920 t_dof[0] = 0; s_dof[0] = 0;
921 t_dof[1] = 1; s_dof[1] = 0;
922 t_dof[2] = 2; s_dof[2] = 0;
923 t_dof[3] = 0; s_dof[3] = 1;
924 t_dof[4] = 1; s_dof[4] = 1;
925 t_dof[5] = 2; s_dof[5] = 1;
926
927 // Edge DoFs
928 int ne = p-1;
929 for (int i=1; i<p; i++)
930 {
931 t_dof[5 + 0 * ne + i] = 2 + 0 * ne + i; s_dof[5 + 0 * ne + i] = 0;
932 t_dof[5 + 1 * ne + i] = 2 + 1 * ne + i; s_dof[5 + 1 * ne + i] = 0;
933 t_dof[5 + 2 * ne + i] = 2 + 2 * ne + i; s_dof[5 + 2 * ne + i] = 0;
934 t_dof[5 + 3 * ne + i] = 2 + 0 * ne + i; s_dof[5 + 3 * ne + i] = 1;
935 t_dof[5 + 4 * ne + i] = 2 + 1 * ne + i; s_dof[5 + 4 * ne + i] = 1;
936 t_dof[5 + 5 * ne + i] = 2 + 2 * ne + i; s_dof[5 + 5 * ne + i] = 1;
937 t_dof[5 + 6 * ne + i] = 0; s_dof[5 + 6 * ne + i] = i + 1;
938 t_dof[5 + 7 * ne + i] = 1; s_dof[5 + 7 * ne + i] = i + 1;
939 t_dof[5 + 8 * ne + i] = 2; s_dof[5 + 8 * ne + i] = i + 1;
940 }
941
942 // Triangular Face DoFs
943 int k=0;
944 int nt = (p-1)*(p-2)/2;
945 for (int j=1; j<p; j++)
946 {
947 for (int i=1; i<j; i++)
948 {
949 t_dof[6 + 9 * ne + k] = 3 * p + k; s_dof[6 + 9 * ne + k] = 0;
950 t_dof[6 + 9 * ne + nt + k] = 3 * p + k; s_dof[6 + 9 * ne + nt + k] = 1;
951 k++;
952 }
953 }
954
955 // Quadrilateral Face DoFs
956 k=0;
957 int nq = (p-1)*(p-1);
958 for (int j=1; j<p; j++)
959 {
960 for (int i=1; i<p; i++)
961 {
962 t_dof[6 + 9 * ne + 2 * nt + 0 * nq + k] = 2 + 0 * ne + i;
963 t_dof[6 + 9 * ne + 2 * nt + 1 * nq + k] = 2 + 1 * ne + i;
964 t_dof[6 + 9 * ne + 2 * nt + 2 * nq + k] = 2 + 2 * ne + i;
965
966 s_dof[6 + 9 * ne + 2 * nt + 0 * nq + k] = 1 + j;
967 s_dof[6 + 9 * ne + 2 * nt + 1 * nq + k] = 1 + j;
968 s_dof[6 + 9 * ne + 2 * nt + 2 * nq + k] = 1 + j;
969
970 k++;
971 }
972 }
973
974 // Interior DoFs
975 int m=0;
976 for (k=1; k<p; k++)
977 {
978 int l=0;
979 for (int j=1; j<p; j++)
980 {
981 for (int i=1; i<j; i++)
982 {
983 t_dof[6 + 9 * ne + 2 * nt + 3 * nq + m] = 3 * p + l;
984 s_dof[6 + 9 * ne + 2 * nt + 3 * nq + m] = 1 + k;
985 l++; m++;
986 }
987 }
988 }
989
990 // Define Nodes
991 const IntegrationRule & t_Nodes = TriangleFE.GetNodes();
992 const IntegrationRule & s_Nodes = SegmentFE.GetNodes();
993 for (int i=0; i<dof; i++)
994 {
995 Nodes.IntPoint(i).x = t_Nodes.IntPoint(t_dof[i]).x;
996 Nodes.IntPoint(i).y = t_Nodes.IntPoint(t_dof[i]).y;
997 Nodes.IntPoint(i).z = s_Nodes.IntPoint(s_dof[i]).x;
998 }
999}
1000
1002 Vector &shape) const
1003{
1004#ifdef MFEM_THREAD_SAFE
1007#endif
1008
1009 IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1010
1013
1014 for (int i=0; i<dof; i++)
1015 {
1016 shape[i] = t_shape[t_dof[i]] * s_shape[s_dof[i]];
1017 }
1018}
1019
1021 DenseMatrix &dshape) const
1022{
1023#ifdef MFEM_THREAD_SAFE
1028#endif
1029
1030 IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1031
1036
1037 for (int i=0; i<dof; i++)
1038 {
1039 dshape(i, 0) = t_dshape(t_dof[i],0) * s_shape[s_dof[i]];
1040 dshape(i, 1) = t_dshape(t_dof[i],1) * s_shape[s_dof[i]];
1041 dshape(i, 2) = t_shape[t_dof[i]] * s_dshape(s_dof[i],0);
1042 }
1043}
1044
1046 : PositiveTensorFiniteElement(1, p, L2_DOF_MAP)
1047{
1048#ifndef MFEM_THREAD_SAFE
1049 shape_x.SetSize(p + 1);
1050 dshape_x.SetDataAndSize(NULL, p + 1);
1051#endif
1052
1053 if (p == 0)
1054 {
1055 Nodes.IntPoint(0).x = 0.5;
1056 }
1057 else
1058 {
1059 for (int i = 0; i <= p; i++)
1060 {
1061 Nodes.IntPoint(i).x = real_t(i)/p;
1062 }
1063 }
1064}
1065
1067 Vector &shape) const
1068{
1069 Poly_1D::CalcBernstein(order, ip.x, shape);
1070}
1071
1073 DenseMatrix &dshape) const
1074{
1075#ifdef MFEM_THREAD_SAFE
1076 Vector shape_x(dof), dshape_x(dshape.Data(), dof);
1077#else
1078 dshape_x.SetData(dshape.Data());
1079#endif
1080 Poly_1D::CalcBernstein(order, ip.x, shape_x, dshape_x);
1081}
1082
1083void L2Pos_SegmentElement::ProjectDelta(int vertex, Vector &dofs) const
1084{
1085 dofs = 0.0;
1086 dofs[vertex*order] = 1.0;
1087}
1088
1089
1091 : PositiveTensorFiniteElement(2, p, L2_DOF_MAP)
1092{
1093#ifndef MFEM_THREAD_SAFE
1094 shape_x.SetSize(p + 1);
1095 shape_y.SetSize(p + 1);
1096 dshape_x.SetSize(p + 1);
1097 dshape_y.SetSize(p + 1);
1098#endif
1099
1100 if (p == 0)
1101 {
1102 Nodes.IntPoint(0).Set2(0.5, 0.5);
1103 }
1104 else
1105 {
1106 for (int o = 0, j = 0; j <= p; j++)
1107 for (int i = 0; i <= p; i++)
1108 {
1109 Nodes.IntPoint(o++).Set2(real_t(i)/p, real_t(j)/p);
1110 }
1111 }
1112}
1113
1115 Vector &shape) const
1116{
1117 const int p = order;
1118
1119#ifdef MFEM_THREAD_SAFE
1120 Vector shape_x(p+1), shape_y(p+1);
1121#endif
1122
1123 Poly_1D::CalcBernstein(p, ip.x, shape_x);
1124 Poly_1D::CalcBernstein(p, ip.y, shape_y);
1125
1126 for (int o = 0, j = 0; j <= p; j++)
1127 for (int i = 0; i <= p; i++)
1128 {
1129 shape(o++) = shape_x(i)*shape_y(j);
1130 }
1131}
1132
1134 DenseMatrix &dshape) const
1135{
1136 const int p = order;
1137
1138#ifdef MFEM_THREAD_SAFE
1139 Vector shape_x(p+1), shape_y(p+1), dshape_x(p+1), dshape_y(p+1);
1140#endif
1141
1142 Poly_1D::CalcBernstein(p, ip.x, shape_x, dshape_x);
1143 Poly_1D::CalcBernstein(p, ip.y, shape_y, dshape_y);
1144
1145 for (int o = 0, j = 0; j <= p; j++)
1146 for (int i = 0; i <= p; i++)
1147 {
1148 dshape(o,0) = dshape_x(i)* shape_y(j);
1149 dshape(o,1) = shape_x(i)*dshape_y(j); o++;
1150 }
1151}
1152
1154{
1155 const int p = order;
1156
1157 dofs = 0.0;
1158 switch (vertex)
1159 {
1160 case 0: dofs[0] = 1.0; break;
1161 case 1: dofs[p] = 1.0; break;
1162 case 2: dofs[p*(p + 2)] = 1.0; break;
1163 case 3: dofs[p*(p + 1)] = 1.0; break;
1164 }
1165}
1166
1167
1169 : PositiveTensorFiniteElement(3, p, L2_DOF_MAP)
1170{
1171#ifndef MFEM_THREAD_SAFE
1172 shape_x.SetSize(p + 1);
1173 shape_y.SetSize(p + 1);
1174 shape_z.SetSize(p + 1);
1175 dshape_x.SetSize(p + 1);
1176 dshape_y.SetSize(p + 1);
1177 dshape_z.SetSize(p + 1);
1178#endif
1179
1180 if (p == 0)
1181 {
1182 Nodes.IntPoint(0).Set3(0.5, 0.5, 0.5);
1183 }
1184 else
1185 {
1186 for (int o = 0, k = 0; k <= p; k++)
1187 for (int j = 0; j <= p; j++)
1188 for (int i = 0; i <= p; i++)
1189 {
1190 Nodes.IntPoint(o++).Set3(real_t(i)/p, real_t(j)/p, real_t(k)/p);
1191 }
1192 }
1193}
1194
1196 Vector &shape) const
1197{
1198 const int p = order;
1199
1200#ifdef MFEM_THREAD_SAFE
1201 Vector shape_x(p+1), shape_y(p+1), shape_z(p+1);
1202#endif
1203
1204 Poly_1D::CalcBernstein(p, ip.x, shape_x);
1205 Poly_1D::CalcBernstein(p, ip.y, shape_y);
1206 Poly_1D::CalcBernstein(p, ip.z, shape_z);
1207
1208 for (int o = 0, k = 0; k <= p; k++)
1209 for (int j = 0; j <= p; j++)
1210 for (int i = 0; i <= p; i++)
1211 {
1212 shape(o++) = shape_x(i)*shape_y(j)*shape_z(k);
1213 }
1214}
1215
1217 DenseMatrix &dshape) const
1218{
1219 const int p = order;
1220
1221#ifdef MFEM_THREAD_SAFE
1222 Vector shape_x(p+1), shape_y(p+1), shape_z(p+1);
1223 Vector dshape_x(p+1), dshape_y(p+1), dshape_z(p+1);
1224#endif
1225
1226 Poly_1D::CalcBernstein(p, ip.x, shape_x, dshape_x);
1227 Poly_1D::CalcBernstein(p, ip.y, shape_y, dshape_y);
1228 Poly_1D::CalcBernstein(p, ip.z, shape_z, dshape_z);
1229
1230 for (int o = 0, k = 0; k <= p; k++)
1231 for (int j = 0; j <= p; j++)
1232 for (int i = 0; i <= p; i++)
1233 {
1234 dshape(o,0) = dshape_x(i)* shape_y(j)* shape_z(k);
1235 dshape(o,1) = shape_x(i)*dshape_y(j)* shape_z(k);
1236 dshape(o,2) = shape_x(i)* shape_y(j)*dshape_z(k); o++;
1237 }
1238}
1239
1241{
1242 const int p = order;
1243
1244 dofs = 0.0;
1245 switch (vertex)
1246 {
1247 case 0: dofs[0] = 1.0; break;
1248 case 1: dofs[p] = 1.0; break;
1249 case 2: dofs[p*(p + 2)] = 1.0; break;
1250 case 3: dofs[p*(p + 1)] = 1.0; break;
1251 case 4: dofs[p*(p + 1)*(p + 1)] = 1.0; break;
1252 case 5: dofs[p + p*(p + 1)*(p + 1)] = 1.0; break;
1253 case 6: dofs[dof - 1] = 1.0; break;
1254 case 7: dofs[dof - p - 1] = 1.0; break;
1255 }
1256}
1257
1258
1260 : PositiveFiniteElement(2, Geometry::TRIANGLE, ((p + 1)*(p + 2))/2, p,
1261 FunctionSpace::Pk)
1262{
1263#ifndef MFEM_THREAD_SAFE
1264 dshape_1d.SetSize(p + 1);
1265#endif
1266
1267 if (p == 0)
1268 {
1269 Nodes.IntPoint(0).Set2(1./3, 1./3);
1270 }
1271 else
1272 {
1273 for (int o = 0, j = 0; j <= p; j++)
1274 for (int i = 0; i + j <= p; i++)
1275 {
1276 Nodes.IntPoint(o++).Set2(real_t(i)/p, real_t(j)/p);
1277 }
1278 }
1279}
1280
1282 Vector &shape) const
1283{
1285}
1286
1288 DenseMatrix &dshape) const
1289{
1290#ifdef MFEM_THREAD_SAFE
1291 Vector dshape_1d(order + 1);
1292#endif
1293
1294 H1Pos_TriangleElement::CalcDShape(order, ip.x, ip.y, dshape_1d.GetData(),
1295 dshape.Data());
1296}
1297
1298void L2Pos_TriangleElement::ProjectDelta(int vertex, Vector &dofs) const
1299{
1300 dofs = 0.0;
1301 switch (vertex)
1302 {
1303 case 0: dofs[0] = 1.0; break;
1304 case 1: dofs[order] = 1.0; break;
1305 case 2: dofs[dof-1] = 1.0; break;
1306 }
1307}
1308
1309
1311 : PositiveFiniteElement(3, Geometry::TETRAHEDRON,
1312 ((p + 1)*(p + 2)*(p + 3))/6, p, FunctionSpace::Pk)
1313{
1314#ifndef MFEM_THREAD_SAFE
1315 dshape_1d.SetSize(p + 1);
1316#endif
1317
1318 if (p == 0)
1319 {
1320 Nodes.IntPoint(0).Set3(0.25, 0.25, 0.25);
1321 }
1322 else
1323 {
1324 for (int o = 0, k = 0; k <= p; k++)
1325 for (int j = 0; j + k <= p; j++)
1326 for (int i = 0; i + j + k <= p; i++)
1327 {
1328 Nodes.IntPoint(o++).Set3(real_t(i)/p, real_t(j)/p, real_t(k)/p);
1329 }
1330 }
1331}
1332
1334 Vector &shape) const
1335{
1337 shape.GetData());
1338}
1339
1341 DenseMatrix &dshape) const
1342{
1343#ifdef MFEM_THREAD_SAFE
1344 Vector dshape_1d(order + 1);
1345#endif
1346
1348 dshape_1d.GetData(), dshape.Data());
1349}
1350
1352{
1353 dofs = 0.0;
1354 switch (vertex)
1355 {
1356 case 0: dofs[0] = 1.0; break;
1357 case 1: dofs[order] = 1.0; break;
1358 case 2: dofs[(order*(order+3))/2] = 1.0; break;
1359 case 3: dofs[dof-1] = 1.0; break;
1360 }
1361}
1362
1363
1365 : PositiveFiniteElement(3, Geometry::PRISM,
1366 ((p + 1)*(p + 1)*(p + 2))/2, p, FunctionSpace::Qk),
1367 TriangleFE(p),
1368 SegmentFE(p)
1369{
1370#ifndef MFEM_THREAD_SAFE
1375#endif
1376
1377 t_dof.SetSize(dof);
1378 s_dof.SetSize(dof);
1379
1380 // Interior DoFs
1381 int m=0;
1382 for (int k=0; k<=p; k++)
1383 {
1384 int l=0;
1385 for (int j=0; j<=p; j++)
1386 {
1387 for (int i=0; i<=j; i++)
1388 {
1389 t_dof[m] = l;
1390 s_dof[m] = k;
1391 l++; m++;
1392 }
1393 }
1394 }
1395
1396 // Define Nodes
1397 const IntegrationRule & t_Nodes = TriangleFE.GetNodes();
1398 const IntegrationRule & s_Nodes = SegmentFE.GetNodes();
1399 for (int i=0; i<dof; i++)
1400 {
1401 Nodes.IntPoint(i).x = t_Nodes.IntPoint(t_dof[i]).x;
1402 Nodes.IntPoint(i).y = t_Nodes.IntPoint(t_dof[i]).y;
1403 Nodes.IntPoint(i).z = s_Nodes.IntPoint(s_dof[i]).x;
1404 }
1405}
1406
1408 Vector &shape) const
1409{
1410#ifdef MFEM_THREAD_SAFE
1413#endif
1414
1415 IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1416
1419
1420 for (int i=0; i<dof; i++)
1421 {
1422 shape[i] = t_shape[t_dof[i]] * s_shape[s_dof[i]];
1423 }
1424}
1425
1427 DenseMatrix &dshape) const
1428{
1429#ifdef MFEM_THREAD_SAFE
1434#endif
1435
1436 IntegrationPoint ipz; ipz.x = ip.z; ipz.y = 0.0; ipz.z = 0.0;
1437
1442
1443 for (int i=0; i<dof; i++)
1444 {
1445 dshape(i, 0) = t_dshape(t_dof[i],0) * s_shape[s_dof[i]];
1446 dshape(i, 1) = t_dshape(t_dof[i],1) * s_shape[s_dof[i]];
1447 dshape(i, 2) = t_shape[t_dof[i]] * s_dshape(s_dof[i],0);
1448 }
1449}
1450
1451}
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
Possible basis types. Note that not all elements can use all BasisType(s).
Definition fe_base.hpp:26
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
Definition fe_pos.cpp:220
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_pos.cpp:118
BiQuadPos2DFiniteElement()
Construct the BiQuadPos2DFiniteElement.
Definition fe_pos.cpp:95
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_pos.cpp:142
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Return the local interpolation matrix I (Dof x Dof) where the fine element is the image of the base g...
Definition fe_pos.cpp:191
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
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
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
void Invert()
Replaces the current matrix with its inverse.
Definition densemat.cpp:726
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
int dof
Number of degrees of freedom.
Definition fe_base.hpp:248
IntegrationRule Nodes
Definition fe_base.hpp:251
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:395
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
H1Pos_HexahedronElement(const int p)
Construct the H1Pos_HexahedronElement of order p.
Definition fe_pos.cpp:432
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_pos.cpp:475
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_pos.cpp:454
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:499
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_pos.cpp:404
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_pos.cpp:384
H1Pos_QuadrilateralElement(const int p)
Construct the H1Pos_QuadrilateralElement of order p.
Definition fe_pos.cpp:364
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:425
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:357
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_pos.cpp:317
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_pos.cpp:337
H1Pos_SegmentElement(const int p)
Construct the H1Pos_SegmentElement of order p.
Definition fe_pos.cpp:299
H1Pos_TetrahedronElement(const int p)
Construct the H1Pos_TetrahedronElement of order p.
Definition fe_pos.cpp:651
static void CalcShape(const int p, const real_t x, const real_t y, const real_t z, real_t *shape)
Definition fe_pos.cpp:753
static void CalcDShape(const int p, const real_t x, const real_t y, const real_t z, real_t *dshape_1d, real_t *dshape)
Definition fe_pos.cpp:786
H1Pos_TriangleElement(const int p)
Construct the H1Pos_TriangleElement of order p.
Definition fe_pos.cpp:506
static void CalcDShape(const int p, const real_t x, const real_t y, real_t *dshape_1d, real_t *dshape)
Definition fe_pos.cpp:587
static void CalcShape(const int p, const real_t x, const real_t y, real_t *shape)
Definition fe_pos.cpp:561
H1Pos_TriangleElement TriangleFE
Definition fe_pos.hpp:246
H1Pos_SegmentElement SegmentFE
Definition fe_pos.hpp:247
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_pos.cpp:1020
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_pos.cpp:1001
H1Pos_WedgeElement(const int p)
Construct the H1Pos_WedgeElement of order p.
Definition fe_pos.cpp:903
Class for integration point with weight.
Definition intrules.hpp:35
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
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
L2Pos_HexahedronElement(const int p)
Construct the L2Pos_HexahedronElement of order p.
Definition fe_pos.cpp:1168
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_pos.cpp:1195
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_pos.cpp:1216
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:1240
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_pos.cpp:1114
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:1153
L2Pos_QuadrilateralElement(const int p)
Construct the L2Pos_QuadrilateralElement of order p.
Definition fe_pos.cpp:1090
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_pos.cpp:1133
L2Pos_SegmentElement(const int p)
Construct the L2Pos_SegmentElement of order p.
Definition fe_pos.cpp:1045
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_pos.cpp:1066
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_pos.cpp:1072
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:1083
L2Pos_TetrahedronElement(const int p)
Construct the L2Pos_TetrahedronElement of order p.
Definition fe_pos.cpp:1310
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_pos.cpp:1333
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_pos.cpp:1340
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:1351
L2Pos_TriangleElement(const int p)
Construct the L2Pos_TriangleElement of order p.
Definition fe_pos.cpp:1259
virtual void ProjectDelta(int vertex, Vector &dofs) const
Project a delta function centered on the given vertex in the local finite dimensional space represent...
Definition fe_pos.cpp:1298
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_pos.cpp:1287
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_pos.cpp:1281
L2Pos_TriangleElement TriangleFE
Definition fe_pos.hpp:360
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_pos.cpp:1426
L2Pos_WedgeElement(const int p)
Construct the L2Pos_WedgeElement of order p.
Definition fe_pos.cpp:1364
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_pos.cpp:1407
L2Pos_SegmentElement SegmentFE
Definition fe_pos.hpp:361
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Class for standard nodal finite elements.
Definition fe_base.hpp:715
void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const override
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
Definition fe_base.cpp:796
static void CalcDBinomTerms(const int p, const real_t x, const real_t y, real_t *d)
Compute the derivatives (w.r.t. x) of the terms in the expansion of the binomial (x + y)^p assuming t...
Definition fe_base.cpp:2141
static const int * Binom(const int p)
Get a pointer to an array containing the binomial coefficients "pchoose k" for k=0,...
Definition fe_base.cpp:2028
static void CalcBernstein(const int p, const real_t x, real_t *u)
Compute the values of the Bernstein basis functions of order p at coordinate x and store the results ...
Definition fe_base.hpp:1164
static void CalcBinomTerms(const int p, const real_t x, const real_t y, real_t *u)
Compute the p terms in the expansion of the binomial (x + y)^p and store them in the already allocate...
Definition fe_base.cpp:2077
Class for finite elements utilizing the always positive Bernstein basis.
Definition fe_pos.hpp:23
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Given a coefficient and a transformation, compute its projection (approximation) in the local finite ...
Definition fe_pos.cpp:25
PositiveTensorFiniteElement(const int dims, const int p, const DofMapType dmtype)
Definition fe_pos.cpp:81
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_pos.cpp:88
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_pos.cpp:278
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_pos.cpp:288
QuadPos1DFiniteElement()
Construct the QuadPos1DFiniteElement.
Definition fe_pos.cpp:270
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 ...
Vector data type.
Definition vector.hpp:80
void SetDataAndSize(real_t *d, int s)
Set the Vector data and size.
Definition vector.hpp:175
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 * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
void SetData(real_t *d)
Definition vector.hpp:168
Linear1DFiniteElement SegmentFE
Definition segment.cpp:52
float real_t
Definition config.hpp:43
MFEM_EXPORT Linear2DFiniteElement TriangleFE
Definition fe.cpp:32
real_t p(const Vector &x, real_t t)
RefCoord s[3]