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