MFEM  v4.6.0
Finite element discretization library
fe_base.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 // Finite Element Base classes
13 
14 #include "fe_base.hpp"
15 #include "face_map_utils.hpp"
16 #include "../coefficient.hpp"
17 
18 namespace mfem
19 {
20 
21 using namespace std;
22 
24  int Do, int O, int F)
25  : Nodes(Do)
26 {
27  dim = D ; geom_type = G ; dof = Do ; order = O ; func_space = F;
28  vdim = 0 ; cdim = 0;
30  map_type = VALUE;
31  deriv_type = NONE;
34  for (int i = 0; i < Geometry::MaxDim; i++) { orders[i] = -1; }
35 #ifndef MFEM_THREAD_SAFE
37 #endif
38 }
39 
41  const IntegrationPoint &ip, DenseMatrix &shape) const
42 {
43  MFEM_ABORT("method is not implemented for this class");
44 }
45 
47  ElementTransformation &Trans, DenseMatrix &shape) const
48 {
49  MFEM_ABORT("method is not implemented for this class");
50 }
51 
53  const IntegrationPoint &ip, Vector &divshape) const
54 {
55  MFEM_ABORT("method is not implemented for this class");
56 }
57 
59  ElementTransformation &Trans, Vector &div_shape) const
60 {
61  CalcDivShape(Trans.GetIntPoint(), div_shape);
62  div_shape *= (1.0 / Trans.Weight());
63 }
64 
66  DenseMatrix &curl_shape) const
67 {
68  MFEM_ABORT("method is not implemented for this class");
69 }
70 
72  DenseMatrix &curl_shape) const
73 {
74  switch (dim)
75  {
76  case 3:
77  {
78 #ifdef MFEM_THREAD_SAFE
80 #endif
82  MultABt(vshape, Trans.Jacobian(), curl_shape);
83  curl_shape *= (1.0 / Trans.Weight());
84  break;
85  }
86  case 2:
87  // This is valid for both 2x2 and 3x2 Jacobians
88  CalcCurlShape(Trans.GetIntPoint(), curl_shape);
89  curl_shape *= (1.0 / Trans.Weight());
90  break;
91  default:
92  MFEM_ABORT("Invalid dimension, Dim = " << dim);
93  }
94 }
95 
96 void FiniteElement::GetFaceDofs(int face, int **dofs, int *ndofs) const
97 {
98  MFEM_ABORT("method is not overloaded");
99 }
100 
102  DenseMatrix &h) const
103 {
104  MFEM_ABORT("method is not overloaded");
105 }
106 
108  DenseMatrix &I) const
109 {
110  MFEM_ABORT("method is not overloaded");
111 }
112 
114  DenseMatrix &) const
115 {
116  MFEM_ABORT("method is not overloaded");
117 }
118 
120  ElementTransformation &Trans,
121  DenseMatrix &I) const
122 {
123  MFEM_ABORT("method is not overloaded");
124 }
125 
127  Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
128 {
129  MFEM_ABORT("method is not overloaded");
130 }
131 
133  VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
134 {
135  MFEM_ABORT("method is not overloaded");
136 }
137 
139  Vector &dofs) const
140 {
141  mfem_error("FiniteElement::ProjectFromNodes() (vector) is not overloaded!");
142 }
143 
145  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
146 {
147  MFEM_ABORT("method is not overloaded");
148 }
149 
150 void FiniteElement::ProjectDelta(int vertex, Vector &dofs) const
151 {
152  MFEM_ABORT("method is not implemented for this element");
153 }
154 
156  const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
157 {
158  MFEM_ABORT("method is not implemented for this element");
159 }
160 
162  const FiniteElement &fe, ElementTransformation &Trans,
163  DenseMatrix &grad) const
164 {
165  MFEM_ABORT("method is not implemented for this element");
166 }
167 
169  const FiniteElement &fe, ElementTransformation &Trans,
170  DenseMatrix &curl) const
171 {
172  MFEM_ABORT("method is not implemented for this element");
173 }
174 
176  const FiniteElement &fe, ElementTransformation &Trans,
177  DenseMatrix &div) const
178 {
179  MFEM_ABORT("method is not implemented for this element");
180 }
181 
183  Vector &shape) const
184 {
185  CalcShape(Trans.GetIntPoint(), shape);
186  if (map_type == INTEGRAL)
187  {
188  shape /= Trans.Weight();
189  }
190 }
191 
193  DenseMatrix &dshape) const
194 {
195  MFEM_ASSERT(map_type == VALUE, "");
196 #ifdef MFEM_THREAD_SAFE
198 #endif
199  CalcDShape(Trans.GetIntPoint(), vshape);
200  Mult(vshape, Trans.InverseJacobian(), dshape);
201 }
202 
204  Vector &Laplacian) const
205 {
206  MFEM_ASSERT(map_type == VALUE, "");
207 
208  // Simpler routine if mapping is affine
209  if (Trans.Hessian().FNorm2() < 1e-20)
210  {
211  CalcPhysLinLaplacian(Trans, Laplacian);
212  return;
213  }
214 
215  // Compute full Hessian first if non-affine
216  int size = (dim*(dim+1))/2;
217  DenseMatrix hess(dof, size);
218  CalcPhysHessian(Trans,hess);
219 
220  if (dim == 3)
221  {
222  for (int nd = 0; nd < dof; nd++)
223  {
224  Laplacian[nd] = hess(nd,0) + hess(nd,4) + hess(nd,5);
225  }
226  }
227  else if (dim == 2)
228  {
229  for (int nd = 0; nd < dof; nd++)
230  {
231  Laplacian[nd] = hess(nd,0) + hess(nd,2);
232  }
233  }
234  else
235  {
236  for (int nd = 0; nd < dof; nd++)
237  {
238  Laplacian[nd] = hess(nd,0);
239  }
240  }
241 }
242 
243 // Assume a linear mapping
245  Vector &Laplacian) const
246 {
247  MFEM_ASSERT(map_type == VALUE, "");
248  int size = (dim*(dim+1))/2;
249  DenseMatrix hess(dof, size);
250  DenseMatrix Gij(dim,dim);
251  Vector scale(size);
252 
253  CalcHessian(Trans.GetIntPoint(), hess);
254  MultAAt(Trans.InverseJacobian(), Gij);
255 
256  if (dim == 3)
257  {
258  scale[0] = Gij(0,0);
259  scale[1] = 2*Gij(0,1);
260  scale[2] = 2*Gij(0,2);
261 
262  scale[3] = 2*Gij(1,2);
263  scale[4] = Gij(2,2);
264 
265  scale[5] = Gij(1,1);
266  }
267  else if (dim == 2)
268  {
269  scale[0] = Gij(0,0);
270  scale[1] = 2*Gij(0,1);
271  scale[2] = Gij(1,1);
272  }
273  else
274  {
275  scale[0] = Gij(0,0);
276  }
277 
278  for (int nd = 0; nd < dof; nd++)
279  {
280  Laplacian[nd] = 0.0;
281  for (int ii = 0; ii < size; ii++)
282  {
283  Laplacian[nd] += hess(nd,ii)*scale[ii];
284  }
285  }
286 }
287 
289  DenseMatrix& Hessian) const
290 {
291  MFEM_ASSERT(map_type == VALUE, "");
292 
293  // Roll 2-Tensors in vectors and 4-Tensor in Matrix, exploiting symmetry
294  Array<int> map(dim*dim);
295  if (dim == 3)
296  {
297  map[0] = 0;
298  map[1] = 1;
299  map[2] = 2;
300 
301  map[3] = 1;
302  map[4] = 5;
303  map[5] = 3;
304 
305  map[6] = 2;
306  map[7] = 3;
307  map[8] = 4;
308  }
309  else if (dim == 2)
310  {
311  map[0] = 0;
312  map[1] = 1;
313 
314  map[2] = 1;
315  map[3] = 2;
316  }
317  else
318  {
319  map[0] = 0;
320  }
321 
322  // Hessian in ref coords
323  int size = (dim*(dim+1))/2;
324  DenseMatrix hess(dof, size);
325  CalcHessian(Trans.GetIntPoint(), hess);
326 
327  // Gradient in physical coords
328  if (Trans.Hessian().FNorm2() > 1e-10)
329  {
330  DenseMatrix grad(dof, dim);
331  CalcPhysDShape(Trans, grad);
332  DenseMatrix gmap(dof, size);
333  Mult(grad,Trans.Hessian(),gmap);
334  hess -= gmap;
335  }
336 
337  // LHM
338  DenseMatrix lhm(size,size);
339  DenseMatrix invJ = Trans.Jacobian();
340  lhm = 0.0;
341  for (int i = 0; i < dim; i++)
342  {
343  for (int j = 0; j < dim; j++)
344  {
345  for (int k = 0; k < dim; k++)
346  {
347  for (int l = 0; l < dim; l++)
348  {
349  lhm(map[i*dim+j],map[k*dim+l]) += invJ(i,k)*invJ(j,l);
350  }
351  }
352  }
353  }
354  // Correct multiplicity
355  Vector mult(size);
356  mult = 0.0;
357  for (int i = 0; i < dim*dim; i++) { mult[map[i]]++; }
358  lhm.InvRightScaling(mult);
359 
360  // Hessian in physical coords
361  lhm.Invert();
362  Mult( hess, lhm, Hessian);
363 }
364 
366  DofToQuad::Mode mode) const
367 {
368  MFEM_VERIFY(mode == DofToQuad::FULL, "invalid mode requested");
369 
370  for (int i = 0; i < dof2quad_array.Size(); i++)
371  {
372  const DofToQuad &d2q = *dof2quad_array[i];
373  if (d2q.IntRule == &ir && d2q.mode == mode) { return d2q; }
374  }
375 
376 #ifdef MFEM_THREAD_SAFE
378 #endif
379 
380  DofToQuad *d2q = new DofToQuad;
381  const int nqpt = ir.GetNPoints();
382  d2q->FE = this;
383  d2q->IntRule = &ir;
384  d2q->mode = mode;
385  d2q->ndof = dof;
386  d2q->nqpt = nqpt;
387  if (range_type == SCALAR)
388  {
389  d2q->B.SetSize(nqpt*dof);
390  d2q->Bt.SetSize(dof*nqpt);
391 
392  Vector shape;
393  vshape.GetColumnReference(0, shape);
394  for (int i = 0; i < nqpt; i++)
395  {
396  const IntegrationPoint &ip = ir.IntPoint(i);
397  CalcShape(ip, shape);
398  for (int j = 0; j < dof; j++)
399  {
400  d2q->B[i+nqpt*j] = d2q->Bt[j+dof*i] = shape(j);
401  }
402  }
403  }
404  else if (range_type == VECTOR)
405  {
406  d2q->B.SetSize(nqpt*dim*dof);
407  d2q->Bt.SetSize(dof*nqpt*dim);
408 
409  for (int i = 0; i < nqpt; i++)
410  {
411  const IntegrationPoint &ip = ir.IntPoint(i);
412  CalcVShape(ip, vshape);
413  for (int d = 0; d < dim; d++)
414  {
415  for (int j = 0; j < dof; j++)
416  {
417  d2q->B[i+nqpt*(d+dim*j)] = d2q->Bt[j+dof*(i+nqpt*d)] = vshape(j, d);
418  }
419  }
420  }
421  }
422  else
423  {
424  // Skip B and Bt for unknown range type
425  }
426  switch (deriv_type)
427  {
428  case GRAD:
429  {
430  d2q->G.SetSize(nqpt*dim*dof);
431  d2q->Gt.SetSize(dof*nqpt*dim);
432 
433  for (int i = 0; i < nqpt; i++)
434  {
435  const IntegrationPoint &ip = ir.IntPoint(i);
436  CalcDShape(ip, vshape);
437  for (int d = 0; d < dim; d++)
438  {
439  for (int j = 0; j < dof; j++)
440  {
441  d2q->G[i+nqpt*(d+dim*j)] = d2q->Gt[j+dof*(i+nqpt*d)] = vshape(j, d);
442  }
443  }
444  }
445  break;
446  }
447  case DIV:
448  {
449  d2q->G.SetSize(nqpt*dof);
450  d2q->Gt.SetSize(dof*nqpt);
451 
452  Vector divshape;
453  vshape.GetColumnReference(0, divshape);
454  for (int i = 0; i < nqpt; i++)
455  {
456  const IntegrationPoint &ip = ir.IntPoint(i);
457  CalcDivShape(ip, divshape);
458  for (int j = 0; j < dof; j++)
459  {
460  d2q->G[i+nqpt*j] = d2q->Gt[j+dof*i] = divshape(j);
461  }
462  }
463  break;
464  }
465  case CURL:
466  {
467  d2q->G.SetSize(nqpt*cdim*dof);
468  d2q->Gt.SetSize(dof*nqpt*cdim);
469 
470  DenseMatrix curlshape(vshape.GetData(), dof, cdim); // cdim <= dim
471  for (int i = 0; i < nqpt; i++)
472  {
473  const IntegrationPoint &ip = ir.IntPoint(i);
474  CalcCurlShape(ip, curlshape);
475  for (int d = 0; d < cdim; d++)
476  {
477  for (int j = 0; j < dof; j++)
478  {
479  d2q->G[i+nqpt*(d+cdim*j)] = d2q->Gt[j+dof*(i+nqpt*d)] = curlshape(j, d);
480  }
481  }
482  }
483  break;
484  }
485  case NONE:
486  default:
487  // Skip G and Gt for unknown derivative type
488  break;
489  }
490  dof2quad_array.Append(d2q);
491  return *d2q;
492 }
493 
494 void FiniteElement::GetFaceMap(const int face_id,
495  Array<int> &face_map) const
496 {
497  MFEM_ABORT("method is not implemented for this element");
498 }
499 
501 {
502  for (int i = 0; i < dof2quad_array.Size(); i++)
503  {
504  delete dof2quad_array[i];
505  }
506 }
507 
508 
511  const ScalarFiniteElement &fine_fe) const
512 {
513  double v[Geometry::MaxDim];
514  Vector vv(v, dim);
515  IntegrationPoint f_ip;
516 
517 #ifdef MFEM_THREAD_SAFE
518  Vector shape(dof);
519 #else
520  Vector shape;
521  vshape.GetColumnReference(0, shape);
522 #endif
523 
524  MFEM_ASSERT(map_type == fine_fe.GetMapType(), "");
525 
526  I.SetSize(fine_fe.dof, dof);
527  for (int i = 0; i < fine_fe.dof; i++)
528  {
529  Trans.Transform(fine_fe.Nodes.IntPoint(i), vv);
530  f_ip.Set(v, dim);
531  CalcShape(f_ip, shape);
532  for (int j = 0; j < dof; j++)
533  {
534  if (fabs(I(i,j) = shape(j)) < 1.0e-12)
535  {
536  I(i,j) = 0.0;
537  }
538  }
539  }
540  if (map_type == INTEGRAL)
541  {
542  // assuming Trans is linear; this should be ok for all refinement types
544  I *= Trans.Weight();
545  }
546 }
547 
550  const ScalarFiniteElement &fine_fe) const
551 {
552  // General "interpolation", defined by L2 projection
553 
554  double v[Geometry::MaxDim];
555  Vector vv(v, dim);
556  IntegrationPoint f_ip;
557 
558  const int fs = fine_fe.GetDof(), cs = this->GetDof();
559  I.SetSize(fs, cs);
560  Vector fine_shape(fs), coarse_shape(cs);
561  DenseMatrix fine_mass(fs), fine_coarse_mass(fs, cs); // initialized with 0
562  const int ir_order =
563  std::max(GetOrder(), fine_fe.GetOrder()) + fine_fe.GetOrder();
564  const IntegrationRule &ir = IntRules.Get(fine_fe.GetGeomType(), ir_order);
565 
566  for (int i = 0; i < ir.GetNPoints(); i++)
567  {
568  const IntegrationPoint &ip = ir.IntPoint(i);
569  fine_fe.CalcShape(ip, fine_shape);
570  Trans.Transform(ip, vv);
571  f_ip.Set(v, dim);
572  this->CalcShape(f_ip, coarse_shape);
573 
574  AddMult_a_VVt(ip.weight, fine_shape, fine_mass);
575  AddMult_a_VWt(ip.weight, fine_shape, coarse_shape, fine_coarse_mass);
576  }
577 
578  DenseMatrixInverse fine_mass_inv(fine_mass);
579  fine_mass_inv.Mult(fine_coarse_mass, I);
580 
581  if (map_type == INTEGRAL)
582  {
583  // assuming Trans is linear; this should be ok for all refinement types
585  I *= Trans.Weight();
586  }
587 }
588 
591  const ScalarFiniteElement &coarse_fe) const
592 {
593  // General "restriction", defined by L2 projection
594  double v[Geometry::MaxDim];
595  Vector vv(v, dim);
596 
597  const int cs = coarse_fe.GetDof(), fs = this->GetDof();
598  R.SetSize(cs, fs);
599  Vector fine_shape(fs), coarse_shape(cs);
600  DenseMatrix coarse_mass(cs), coarse_fine_mass(cs, fs); // initialized with 0
601  const int ir_order = GetOrder() + coarse_fe.GetOrder();
602  const IntegrationRule &ir = IntRules.Get(coarse_fe.GetGeomType(), ir_order);
603 
604  // integrate coarse_mass in the coarse space
605  for (int i = 0; i < ir.GetNPoints(); i++)
606  {
607  const IntegrationPoint &c_ip = ir.IntPoint(i);
608  coarse_fe.CalcShape(c_ip, coarse_shape);
609  AddMult_a_VVt(c_ip.weight, coarse_shape, coarse_mass);
610  }
611 
612  // integrate coarse_fine_mass in the fine space
614  for (int i = 0; i < ir.GetNPoints(); i++)
615  {
616  const IntegrationPoint &f_ip = ir.IntPoint(i);
617  this->CalcShape(f_ip, fine_shape);
618  Trans.Transform(f_ip, vv);
619 
620  IntegrationPoint c_ip;
621  c_ip.Set(v, dim);
622  coarse_fe.CalcShape(c_ip, coarse_shape);
623  AddMult_a_VWt(f_ip.weight*Trans.Weight(), coarse_shape, fine_shape,
624  coarse_fine_mass);
625  }
626 
627  DenseMatrixInverse coarse_mass_inv(coarse_mass);
628  coarse_mass_inv.Mult(coarse_fine_mass, R);
629 
630  if (map_type == INTEGRAL)
631  {
632  // assuming Trans is linear; this should be ok for all refinement types
634  R *= 1.0 / Trans.Weight();
635  }
636 }
637 
639  const FiniteElement &fe, ElementTransformation &Trans,
640  DenseMatrix &curl) const
641 {
642  DenseMatrix curl_shape(fe.GetDof(), 1);
643 
644  curl.SetSize(dof, fe.GetDof());
645  for (int i = 0; i < dof; i++)
646  {
647  fe.CalcCurlShape(Nodes.IntPoint(i), curl_shape);
648 
649  double w = 1.0;
651  {
652  Trans.SetIntPoint(&Nodes.IntPoint(i));
653  w /= Trans.Weight();
654  }
655  for (int j = 0; j < fe.GetDof(); j++)
656  {
657  curl(i,j) = w * curl_shape(j,0);
658  }
659  }
660 }
661 
663  const IntegrationPoint &pt, Vector &x)
664 {
665  // invert a linear transform with one Newton step
666  IntegrationPoint p0;
667  p0.Set3(0, 0, 0);
668  trans.Transform(p0, x);
669 
670  double store[3];
671  Vector v(store, x.Size());
672  pt.Get(store, x.Size());
673  v -= x;
674 
675  trans.InverseJacobian().Mult(v, x);
676 }
677 
679  DenseMatrix &R) const
680 {
681  IntegrationPoint ipt;
682  Vector pt(&ipt.x, dim);
683 
684 #ifdef MFEM_THREAD_SAFE
685  Vector shape(dof);
686 #else
687  Vector shape;
688  vshape.GetColumnReference(0, shape);
689 #endif
690 
691  Trans.SetIntPoint(&Nodes[0]);
692 
693  for (int j = 0; j < dof; j++)
694  {
695  InvertLinearTrans(Trans, Nodes[j], pt);
696  if (Geometries.CheckPoint(geom_type, ipt)) // do we need an epsilon here?
697  {
698  CalcShape(ipt, shape);
699  R.SetRow(j, shape);
700  }
701  else
702  {
703  // Set the whole row to avoid valgrind warnings in R.Threshold().
704  R.SetRow(j, infinity());
705  }
706  }
707  R.Threshold(1e-12);
708 }
709 
711  Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
712 {
713  for (int i = 0; i < dof; i++)
714  {
715  const IntegrationPoint &ip = Nodes.IntPoint(i);
716  // some coefficients expect that Trans.IntPoint is the same
717  // as the second argument of Eval
718  Trans.SetIntPoint(&ip);
719  dofs(i) = coeff.Eval(Trans, ip);
720  if (map_type == INTEGRAL)
721  {
722  dofs(i) *= Trans.Weight();
723  }
724  }
725 }
726 
728  VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
729 {
730  MFEM_ASSERT(dofs.Size() == vc.GetVDim()*dof, "");
731  Vector x(vc.GetVDim());
732 
733  for (int i = 0; i < dof; i++)
734  {
735  const IntegrationPoint &ip = Nodes.IntPoint(i);
736  Trans.SetIntPoint(&ip);
737  vc.Eval (x, Trans, ip);
738  if (map_type == INTEGRAL)
739  {
740  x *= Trans.Weight();
741  }
742  for (int j = 0; j < x.Size(); j++)
743  {
744  dofs(dof*j+i) = x(j);
745  }
746  }
747 }
748 
750  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
751 {
752  // (mc.height x mc.width) @ DOFs -> (dof x mc.width x mc.height) in dofs
753  MFEM_ASSERT(dofs.Size() == mc.GetHeight()*mc.GetWidth()*dof, "");
754  DenseMatrix MQ(mc.GetHeight(), mc.GetWidth());
755 
756  for (int k = 0; k < dof; k++)
757  {
758  T.SetIntPoint(&Nodes.IntPoint(k));
759  mc.Eval(MQ, T, Nodes.IntPoint(k));
760  if (map_type == INTEGRAL) { MQ *= T.Weight(); }
761  for (int r = 0; r < MQ.Height(); r++)
762  {
763  for (int d = 0; d < MQ.Width(); d++)
764  {
765  dofs(k+dof*(d+MQ.Width()*r)) = MQ(r,d);
766  }
767  }
768  }
769 }
770 
772  const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
773 {
774  if (fe.GetRangeType() == SCALAR)
775  {
776  Vector shape(fe.GetDof());
777 
778  I.SetSize(dof, fe.GetDof());
779  if (map_type == fe.GetMapType())
780  {
781  for (int k = 0; k < dof; k++)
782  {
783  fe.CalcShape(Nodes.IntPoint(k), shape);
784  for (int j = 0; j < shape.Size(); j++)
785  {
786  I(k,j) = (fabs(shape(j)) < 1e-12) ? 0.0 : shape(j);
787  }
788  }
789  }
790  else
791  {
792  for (int k = 0; k < dof; k++)
793  {
794  Trans.SetIntPoint(&Nodes.IntPoint(k));
795  fe.CalcPhysShape(Trans, shape);
796  if (map_type == INTEGRAL)
797  {
798  shape *= Trans.Weight();
799  }
800  for (int j = 0; j < shape.Size(); j++)
801  {
802  I(k,j) = (fabs(shape(j)) < 1e-12) ? 0.0 : shape(j);
803  }
804  }
805  }
806  }
807  else
808  {
809  DenseMatrix vshape(fe.GetDof(), std::max(Trans.GetSpaceDim(),
810  fe.GetVDim()));
811 
812  I.SetSize(vshape.Width()*dof, fe.GetDof());
813  for (int k = 0; k < dof; k++)
814  {
815  Trans.SetIntPoint(&Nodes.IntPoint(k));
816  fe.CalcVShape(Trans, vshape);
817  if (map_type == INTEGRAL)
818  {
819  vshape *= Trans.Weight();
820  }
821  for (int j = 0; j < vshape.Height(); j++)
822  for (int d = 0; d < vshape.Width(); d++)
823  {
824  I(k+d*dof,j) = vshape(j,d);
825  }
826  }
827  }
828 }
829 
831  const FiniteElement &fe, ElementTransformation &Trans,
832  DenseMatrix &grad) const
833 {
834  MFEM_ASSERT(fe.GetMapType() == VALUE, "");
835  MFEM_ASSERT(Trans.GetSpaceDim() == dim, "")
836 
837  DenseMatrix dshape(fe.GetDof(), dim), grad_k(fe.GetDof(), dim), Jinv(dim);
838 
839  grad.SetSize(dim*dof, fe.GetDof());
840  for (int k = 0; k < dof; k++)
841  {
842  const IntegrationPoint &ip = Nodes.IntPoint(k);
843  fe.CalcDShape(ip, dshape);
844  Trans.SetIntPoint(&ip);
845  CalcInverse(Trans.Jacobian(), Jinv);
846  Mult(dshape, Jinv, grad_k);
847  if (map_type == INTEGRAL)
848  {
849  grad_k *= Trans.Weight();
850  }
851  for (int j = 0; j < grad_k.Height(); j++)
852  for (int d = 0; d < dim; d++)
853  {
854  grad(k+d*dof,j) = grad_k(j,d);
855  }
856  }
857 }
858 
860  const FiniteElement &fe, ElementTransformation &Trans,
861  DenseMatrix &div) const
862 {
863  double detJ;
864  Vector div_shape(fe.GetDof());
865 
866  div.SetSize(dof, fe.GetDof());
867  for (int k = 0; k < dof; k++)
868  {
869  const IntegrationPoint &ip = Nodes.IntPoint(k);
870  fe.CalcDivShape(ip, div_shape);
871  if (map_type == VALUE)
872  {
873  Trans.SetIntPoint(&ip);
874  detJ = Trans.Weight();
875  for (int j = 0; j < div_shape.Size(); j++)
876  {
877  div(k,j) = (fabs(div_shape(j)) < 1e-12) ? 0.0 : div_shape(j)/detJ;
878  }
879  }
880  else
881  {
882  for (int j = 0; j < div_shape.Size(); j++)
883  {
884  div(k,j) = (fabs(div_shape(j)) < 1e-12) ? 0.0 : div_shape(j);
885  }
886  }
887  }
888 }
889 
890 
892  int Do, int O, int M, int F)
893  : FiniteElement(D, G, Do, O, F)
894 {
895  range_type = VECTOR;
896  map_type = M;
897  SetDerivMembers();
898  is_nodal = true;
899  vdim = dim;
900  if (map_type == H_CURL)
901  {
902  cdim = (dim == 3) ? 3 : 1;
903  }
904 }
905 
906 void VectorFiniteElement::CalcShape(
907  const IntegrationPoint &ip, Vector &shape ) const
908 {
909  mfem_error("Error: Cannot use scalar CalcShape(...) function with\n"
910  " VectorFiniteElements!");
911 }
912 
913 void VectorFiniteElement::CalcDShape(
914  const IntegrationPoint &ip, DenseMatrix &dshape ) const
915 {
916  mfem_error("Error: Cannot use scalar CalcDShape(...) function with\n"
917  " VectorFiniteElements!");
918 }
919 
921 {
922  switch (map_type)
923  {
924  case H_DIV:
925  deriv_type = DIV;
928  break;
929  case H_CURL:
930  switch (dim)
931  {
932  case 3: // curl: 3D H_CURL -> 3D H_DIV
933  deriv_type = CURL;
936  break;
937  case 2:
938  // curl: 2D H_CURL -> INTEGRAL
939  deriv_type = CURL;
942  break;
943  case 1:
944  deriv_type = NONE;
947  break;
948  default:
949  MFEM_ABORT("Invalid dimension, Dim = " << dim);
950  }
951  break;
952  default:
953  MFEM_ABORT("Invalid MapType = " << map_type);
954  }
955 }
956 
958  ElementTransformation &Trans, DenseMatrix &shape) const
959 {
960  MFEM_ASSERT(map_type == H_DIV, "");
961 #ifdef MFEM_THREAD_SAFE
963 #endif
964  CalcVShape(Trans.GetIntPoint(), vshape);
965  MultABt(vshape, Trans.Jacobian(), shape);
966  shape *= (1.0 / Trans.Weight());
967 }
968 
970  ElementTransformation &Trans, DenseMatrix &shape) const
971 {
972  MFEM_ASSERT(map_type == H_CURL, "");
973 #ifdef MFEM_THREAD_SAFE
975 #endif
976  CalcVShape(Trans.GetIntPoint(), vshape);
977  Mult(vshape, Trans.InverseJacobian(), shape);
978 }
979 
981  const double *nk, const Array<int> &d2n,
982  VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
983 {
984  double vk[Geometry::MaxDim];
985  const int sdim = Trans.GetSpaceDim();
986  MFEM_ASSERT(vc.GetVDim() == sdim, "");
987  Vector xk(vk, sdim);
988  const bool square_J = (dim == sdim);
989 
990  for (int k = 0; k < dof; k++)
991  {
992  Trans.SetIntPoint(&Nodes.IntPoint(k));
993  vc.Eval(xk, Trans, Nodes.IntPoint(k));
994  // dof_k = nk^t adj(J) xk
995  dofs(k) = Trans.AdjugateJacobian().InnerProduct(vk, nk + d2n[k]*dim);
996  if (!square_J) { dofs(k) /= Trans.Weight(); }
997  }
998 }
999 
1001  const double *nk, const Array<int> &d2n,
1002  Vector &vc, ElementTransformation &Trans, Vector &dofs) const
1003 {
1004  const int sdim = Trans.GetSpaceDim();
1005  const bool square_J = (dim == sdim);
1006 
1007  for (int k = 0; k < dof; k++)
1008  {
1009  Trans.SetIntPoint(&Nodes.IntPoint(k));
1010  // dof_k = nk^t adj(J) xk
1011  dofs(k) = Trans.AdjugateJacobian().InnerProduct(
1012  &vc[k*sdim], nk + d2n[k]*dim);
1013  if (!square_J) { dofs(k) /= Trans.Weight(); }
1014  }
1015 }
1016 
1018  const double *nk, const Array<int> &d2n,
1019  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
1020 {
1021  // project the rows of the matrix coefficient in an RT space
1022 
1023  const int sdim = T.GetSpaceDim();
1024  MFEM_ASSERT(mc.GetWidth() == sdim, "");
1025  const bool square_J = (dim == sdim);
1026  DenseMatrix MQ(mc.GetHeight(), mc.GetWidth());
1027  Vector nk_phys(sdim), dofs_k(MQ.Height());
1028  MFEM_ASSERT(dofs.Size() == dof*MQ.Height(), "");
1029 
1030  for (int k = 0; k < dof; k++)
1031  {
1032  T.SetIntPoint(&Nodes.IntPoint(k));
1033  mc.Eval(MQ, T, Nodes.IntPoint(k));
1034  // nk_phys = adj(J)^t nk
1035  T.AdjugateJacobian().MultTranspose(nk + d2n[k]*dim, nk_phys);
1036  if (!square_J) { nk_phys /= T.Weight(); }
1037  MQ.Mult(nk_phys, dofs_k);
1038  for (int r = 0; r < MQ.Height(); r++)
1039  {
1040  dofs(k+dof*r) = dofs_k(r);
1041  }
1042  }
1043 }
1044 
1046  const double *nk, const Array<int> &d2n, const FiniteElement &fe,
1047  ElementTransformation &Trans, DenseMatrix &I) const
1048 {
1049  if (fe.GetRangeType() == SCALAR)
1050  {
1051  double vk[Geometry::MaxDim];
1052  Vector shape(fe.GetDof());
1053  int sdim = Trans.GetSpaceDim();
1054 
1055  I.SetSize(dof, sdim*fe.GetDof());
1056  for (int k = 0; k < dof; k++)
1057  {
1058  const IntegrationPoint &ip = Nodes.IntPoint(k);
1059 
1060  fe.CalcShape(ip, shape);
1061  Trans.SetIntPoint(&ip);
1062  // Transform RT face normals from reference to physical space
1063  // vk = adj(J)^T nk
1064  Trans.AdjugateJacobian().MultTranspose(nk + d2n[k]*dim, vk);
1065  if (fe.GetMapType() == INTEGRAL)
1066  {
1067  double w = 1.0/Trans.Weight();
1068  for (int d = 0; d < dim; d++)
1069  {
1070  vk[d] *= w;
1071  }
1072  }
1073 
1074  for (int j = 0; j < shape.Size(); j++)
1075  {
1076  double s = shape(j);
1077  if (fabs(s) < 1e-12)
1078  {
1079  s = 0.0;
1080  }
1081  // Project scalar basis function multiplied by each coordinate
1082  // direction onto the transformed face normals
1083  for (int d = 0; d < sdim; d++)
1084  {
1085  I(k,j+d*shape.Size()) = s*vk[d];
1086  }
1087  }
1088  }
1089  }
1090  else
1091  {
1092  int sdim = Trans.GetSpaceDim();
1093  double vk[Geometry::MaxDim];
1094  DenseMatrix vshape(fe.GetDof(), sdim);
1095  Vector vshapenk(fe.GetDof());
1096  const bool square_J = (dim == sdim);
1097 
1098  I.SetSize(dof, fe.GetDof());
1099  for (int k = 0; k < dof; k++)
1100  {
1101  const IntegrationPoint &ip = Nodes.IntPoint(k);
1102 
1103  Trans.SetIntPoint(&ip);
1104  // Transform RT face normals from reference to physical space
1105  // vk = adj(J)^T nk
1106  Trans.AdjugateJacobian().MultTranspose(nk + d2n[k]*dim, vk);
1107  // Compute fe basis functions in physical space
1108  fe.CalcVShape(Trans, vshape);
1109  // Project fe basis functions onto transformed face normals
1110  vshape.Mult(vk, vshapenk);
1111  if (!square_J) { vshapenk /= Trans.Weight(); }
1112  for (int j=0; j<vshapenk.Size(); j++)
1113  {
1114  I(k,j) = vshapenk(j);
1115  }
1116  }
1117  }
1118 }
1119 
1121  const double *nk, const Array<int> &d2n, const FiniteElement &fe,
1122  ElementTransformation &Trans, DenseMatrix &grad) const
1123 {
1124  if (dim != 2)
1125  {
1126  mfem_error("VectorFiniteElement::ProjectGrad_RT works only in 2D!");
1127  }
1128 
1129  DenseMatrix dshape(fe.GetDof(), fe.GetDim());
1130  Vector grad_k(fe.GetDof());
1131  double tk[2];
1132 
1133  grad.SetSize(dof, fe.GetDof());
1134  for (int k = 0; k < dof; k++)
1135  {
1136  fe.CalcDShape(Nodes.IntPoint(k), dshape);
1137  tk[0] = nk[d2n[k]*dim+1];
1138  tk[1] = -nk[d2n[k]*dim];
1139  dshape.Mult(tk, grad_k);
1140  for (int j = 0; j < grad_k.Size(); j++)
1141  {
1142  grad(k,j) = (fabs(grad_k(j)) < 1e-12) ? 0.0 : grad_k(j);
1143  }
1144  }
1145 }
1146 
1148  const double *tk, const Array<int> &d2t, const FiniteElement &fe,
1149  ElementTransformation &Trans, DenseMatrix &curl) const
1150 {
1151 #ifdef MFEM_THREAD_SAFE
1154  DenseMatrix JtJ(dim, dim);
1155 #else
1156  curlshape.SetSize(fe.GetDof(), dim);
1157  curlshape_J.SetSize(fe.GetDof(), dim);
1158  JtJ.SetSize(dim, dim);
1159 #endif
1160 
1161  Vector curl_k(fe.GetDof());
1162 
1163  curl.SetSize(dof, fe.GetDof());
1164  for (int k = 0; k < dof; k++)
1165  {
1166  const IntegrationPoint &ip = Nodes.IntPoint(k);
1167 
1168  // calculate J^t * J / |J|
1169  Trans.SetIntPoint(&ip);
1170  MultAtB(Trans.Jacobian(), Trans.Jacobian(), JtJ);
1171  JtJ *= 1.0 / Trans.Weight();
1172 
1173  // transform curl of shapes (rows) by J^t * J / |J|
1174  fe.CalcCurlShape(ip, curlshape);
1176 
1177  curlshape_J.Mult(tk + d2t[k]*dim, curl_k);
1178  for (int j = 0; j < curl_k.Size(); j++)
1179  {
1180  curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1181  }
1182  }
1183 }
1184 
1186  const double *nk, const Array<int> &d2n, const FiniteElement &fe,
1187  ElementTransformation &Trans, DenseMatrix &curl) const
1188 {
1189  DenseMatrix curl_shape(fe.GetDof(), dim);
1190  Vector curl_k(fe.GetDof());
1191 
1192  curl.SetSize(dof, fe.GetDof());
1193  for (int k = 0; k < dof; k++)
1194  {
1195  fe.CalcCurlShape(Nodes.IntPoint(k), curl_shape);
1196  curl_shape.Mult(nk + d2n[k]*dim, curl_k);
1197  for (int j = 0; j < curl_k.Size(); j++)
1198  {
1199  curl(k,j) = (fabs(curl_k(j)) < 1e-12) ? 0.0 : curl_k(j);
1200  }
1201  }
1202 }
1203 
1205  const double *tk, const Array<int> &d2t,
1206  VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
1207 {
1208  double vk[Geometry::MaxDim];
1209  Vector xk(vk, vc.GetVDim());
1210 
1211  for (int k = 0; k < dof; k++)
1212  {
1213  Trans.SetIntPoint(&Nodes.IntPoint(k));
1214 
1215  vc.Eval(xk, Trans, Nodes.IntPoint(k));
1216  // dof_k = xk^t J tk
1217  dofs(k) = Trans.Jacobian().InnerProduct(tk + d2t[k]*dim, vk);
1218  }
1219 }
1220 
1222  const double *tk, const Array<int> &d2t,
1223  Vector &vc, ElementTransformation &Trans, Vector &dofs) const
1224 {
1225  for (int k = 0; k < dof; k++)
1226  {
1227  Trans.SetIntPoint(&Nodes.IntPoint(k));
1228  // dof_k = xk^t J tk
1229  dofs(k) = Trans.Jacobian().InnerProduct(tk + d2t[k]*dim, &vc[k*dim]);
1230  }
1231 }
1232 
1234  const double *tk, const Array<int> &d2t,
1235  MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
1236 {
1237  // project the rows of the matrix coefficient in an ND space
1238 
1239  const int sdim = T.GetSpaceDim();
1240  MFEM_ASSERT(mc.GetWidth() == sdim, "");
1241  DenseMatrix MQ(mc.GetHeight(), mc.GetWidth());
1242  Vector tk_phys(sdim), dofs_k(MQ.Height());
1243  MFEM_ASSERT(dofs.Size() == dof*MQ.Height(), "");
1244 
1245  for (int k = 0; k < dof; k++)
1246  {
1247  T.SetIntPoint(&Nodes.IntPoint(k));
1248  mc.Eval(MQ, T, Nodes.IntPoint(k));
1249  // tk_phys = J tk
1250  T.Jacobian().Mult(tk + d2t[k]*dim, tk_phys);
1251  MQ.Mult(tk_phys, dofs_k);
1252  for (int r = 0; r < MQ.Height(); r++)
1253  {
1254  dofs(k+dof*r) = dofs_k(r);
1255  }
1256  }
1257 }
1258 
1260  const double *tk, const Array<int> &d2t, const FiniteElement &fe,
1261  ElementTransformation &Trans, DenseMatrix &I) const
1262 {
1263  if (fe.GetRangeType() == SCALAR)
1264  {
1265  int sdim = Trans.GetSpaceDim();
1266  double vk[Geometry::MaxDim];
1267  Vector shape(fe.GetDof());
1268 
1269  I.SetSize(dof, sdim*fe.GetDof());
1270  for (int k = 0; k < dof; k++)
1271  {
1272  const IntegrationPoint &ip = Nodes.IntPoint(k);
1273 
1274  fe.CalcShape(ip, shape);
1275  Trans.SetIntPoint(&ip);
1276  // Transform ND edge tengents from reference to physical space
1277  // vk = J tk
1278  Trans.Jacobian().Mult(tk + d2t[k]*dim, vk);
1279  if (fe.GetMapType() == INTEGRAL)
1280  {
1281  double w = 1.0/Trans.Weight();
1282  for (int d = 0; d < sdim; d++)
1283  {
1284  vk[d] *= w;
1285  }
1286  }
1287 
1288  for (int j = 0; j < shape.Size(); j++)
1289  {
1290  double s = shape(j);
1291  if (fabs(s) < 1e-12)
1292  {
1293  s = 0.0;
1294  }
1295  // Project scalar basis function multiplied by each coordinate
1296  // direction onto the transformed edge tangents
1297  for (int d = 0; d < sdim; d++)
1298  {
1299  I(k, j + d*shape.Size()) = s*vk[d];
1300  }
1301  }
1302  }
1303  }
1304  else
1305  {
1306  int sdim = Trans.GetSpaceDim();
1307  double vk[Geometry::MaxDim];
1308  DenseMatrix vshape(fe.GetDof(), sdim);
1309  Vector vshapetk(fe.GetDof());
1310 
1311  I.SetSize(dof, fe.GetDof());
1312  for (int k = 0; k < dof; k++)
1313  {
1314  const IntegrationPoint &ip = Nodes.IntPoint(k);
1315 
1316  Trans.SetIntPoint(&ip);
1317  // Transform ND edge tangents from reference to physical space
1318  // vk = J tk
1319  Trans.Jacobian().Mult(tk + d2t[k]*dim, vk);
1320  // Compute fe basis functions in physical space
1321  fe.CalcVShape(Trans, vshape);
1322  // Project fe basis functions onto transformed edge tangents
1323  vshape.Mult(vk, vshapetk);
1324  for (int j=0; j<vshapetk.Size(); j++)
1325  {
1326  I(k, j) = vshapetk(j);
1327  }
1328  }
1329  }
1330 }
1331 
1333  const double *tk, const Array<int> &d2t, const FiniteElement &fe,
1334  ElementTransformation &Trans, DenseMatrix &grad) const
1335 {
1336  MFEM_ASSERT(fe.GetMapType() == VALUE, "");
1337 
1338  DenseMatrix dshape(fe.GetDof(), fe.GetDim());
1339  Vector grad_k(fe.GetDof());
1340 
1341  grad.SetSize(dof, fe.GetDof());
1342  for (int k = 0; k < dof; k++)
1343  {
1344  fe.CalcDShape(Nodes.IntPoint(k), dshape);
1345  dshape.Mult(tk + d2t[k]*dim, grad_k);
1346  for (int j = 0; j < grad_k.Size(); j++)
1347  {
1348  grad(k,j) = (fabs(grad_k(j)) < 1e-12) ? 0.0 : grad_k(j);
1349  }
1350  }
1351 }
1352 
1354  const VectorFiniteElement &cfe, ElementTransformation &Trans,
1355  DenseMatrix &I) const
1356 {
1357  Vector v(dim);
1358  IntegrationPoint tr_ip;
1359 
1360  const int fs = dof, cs = cfe.GetDof();
1361  I.SetSize(fs, cs);
1362  DenseMatrix fine_shape(fs, dim), coarse_shape(cs, cfe.GetDim());
1363  DenseMatrix fine_mass(fs), fine_coarse_mass(fs, cs); // initialized with 0
1364  const int ir_order =
1365  std::max(GetOrder(), this->GetOrder()) + this->GetOrder();
1366  const IntegrationRule &ir = IntRules.Get(this->GetGeomType(), ir_order);
1367 
1369  const DenseMatrix &adjJ = Trans.AdjugateJacobian();
1370  for (int i = 0; i < ir.GetNPoints(); i++)
1371  {
1372  const IntegrationPoint &ip = ir.IntPoint(i);
1373  double w = ip.weight;
1374  this->CalcVShape(ip, fine_shape);
1375  Trans.Transform(ip, v);
1376  tr_ip.Set(v.GetData(), dim);
1377  cfe.CalcVShape(tr_ip, coarse_shape);
1378 
1379  AddMult_a_AAt(w, fine_shape, fine_mass);
1380  for (int k=0; k<fs; ++k)
1381  {
1382  for (int j=0; j<cs; ++j)
1383  {
1384  double Mkj = 0.0;
1385  for (int d1=0; d1<dim; ++d1)
1386  {
1387  for (int d2=0; d2<dim; ++d2)
1388  {
1389  Mkj += w*fine_shape(k,d1)*adjJ(d2,d1)*coarse_shape(j,d2);
1390  }
1391  }
1392  fine_coarse_mass(k,j) += (fabs(Mkj) < 1e-12) ? 0.0 : Mkj;
1393  }
1394  }
1395  }
1396  DenseMatrixInverse fine_mass_inv(fine_mass);
1397  fine_mass_inv.Mult(fine_coarse_mass, I);
1398 }
1399 
1401  const VectorFiniteElement &cfe, const double *nk, const Array<int> &d2n,
1402  ElementTransformation &Trans, DenseMatrix &I) const
1403 {
1404  MFEM_ASSERT(map_type == cfe.GetMapType(), "");
1405 
1406  if (!is_nodal) { return LocalL2Projection_RT(cfe, Trans, I); }
1407 
1408  double vk[Geometry::MaxDim];
1409  Vector xk(vk, dim);
1410  IntegrationPoint ip;
1411 #ifdef MFEM_THREAD_SAFE
1412  DenseMatrix vshape(cfe.GetDof(), cfe.GetDim());
1413 #else
1414  DenseMatrix vshape(cfe.vshape.Data(), cfe.GetDof(), cfe.GetDim());
1415 #endif
1416  I.SetSize(dof, vshape.Height());
1417 
1418  // assuming Trans is linear; this should be ok for all refinement types
1420  const DenseMatrix &adjJ = Trans.AdjugateJacobian();
1421  for (int k = 0; k < dof; k++)
1422  {
1423  Trans.Transform(Nodes.IntPoint(k), xk);
1424  ip.Set3(vk);
1425  cfe.CalcVShape(ip, vshape);
1426  // xk = |J| J^{-t} n_k
1427  adjJ.MultTranspose(nk + d2n[k]*dim, vk);
1428  // I_k = vshape_k.adj(J)^t.n_k, k=1,...,dof
1429  for (int j = 0; j < vshape.Height(); j++)
1430  {
1431  double Ikj = 0.;
1432  for (int i = 0; i < dim; i++)
1433  {
1434  Ikj += vshape(j, i) * vk[i];
1435  }
1436  I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1437  }
1438  }
1439 }
1440 
1442  const VectorFiniteElement &cfe,
1443  ElementTransformation &Trans, DenseMatrix &I) const
1444 {
1445  Vector v(dim);
1446  IntegrationPoint tr_ip;
1447 
1448  const int fs = dof, cs = cfe.GetDof();
1449  I.SetSize(fs, cs);
1450  DenseMatrix fine_shape(fs, dim), coarse_shape(cs, cfe.GetDim());
1451  DenseMatrix fine_mass(fs), fine_coarse_mass(fs, cs); // initialized with 0
1452  const int ir_order =
1453  std::max(GetOrder(), this->GetOrder()) + this->GetOrder();
1454  const IntegrationRule &ir = IntRules.Get(this->GetGeomType(), ir_order);
1455 
1457  const DenseMatrix &J = Trans.Jacobian();
1458  for (int i = 0; i < ir.GetNPoints(); i++)
1459  {
1460  const IntegrationPoint &ip = ir.IntPoint(i);
1461  this->CalcVShape(ip, fine_shape);
1462  Trans.Transform(ip, v);
1463  tr_ip.Set(v.GetData(), dim);
1464  cfe.CalcVShape(tr_ip, coarse_shape);
1465 
1466  AddMult_a_AAt(ip.weight, fine_shape, fine_mass);
1467  for (int k=0; k<fs; ++k)
1468  {
1469  for (int j=0; j<cs; ++j)
1470  {
1471  double Mkj = 0.0;
1472  for (int d1=0; d1<dim; ++d1)
1473  {
1474  for (int d2=0; d2<dim; ++d2)
1475  {
1476  Mkj += ip.weight*fine_shape(k,d1)*J(d1,d2)*coarse_shape(j,d2);
1477  }
1478  }
1479  fine_coarse_mass(k,j) += (fabs(Mkj) < 1e-12) ? 0.0 : Mkj;
1480  }
1481  }
1482  }
1483  DenseMatrixInverse fine_mass_inv(fine_mass);
1484  fine_mass_inv.Mult(fine_coarse_mass, I);
1485 }
1486 
1488  const VectorFiniteElement &cfe, const double *tk, const Array<int> &d2t,
1489  ElementTransformation &Trans, DenseMatrix &I) const
1490 {
1491  if (!is_nodal) { return LocalL2Projection_ND(cfe, Trans, I); }
1492 
1493  double vk[Geometry::MaxDim];
1494  Vector xk(vk, dim);
1495  IntegrationPoint ip;
1496 #ifdef MFEM_THREAD_SAFE
1497  DenseMatrix vshape(cfe.GetDof(), cfe.GetDim());
1498 #else
1499  DenseMatrix vshape(cfe.vshape.Data(), cfe.GetDof(), cfe.GetDim());
1500 #endif
1501  I.SetSize(dof, vshape.Height());
1502 
1503  // assuming Trans is linear; this should be ok for all refinement types
1505  const DenseMatrix &J = Trans.Jacobian();
1506  for (int k = 0; k < dof; k++)
1507  {
1508  Trans.Transform(Nodes.IntPoint(k), xk);
1509  ip.Set3(vk);
1510  cfe.CalcVShape(ip, vshape);
1511  // xk = J t_k
1512  J.Mult(tk + d2t[k]*dim, vk);
1513  // I_k = vshape_k.J.t_k, k=1,...,Dof
1514  for (int j = 0; j < vshape.Height(); j++)
1515  {
1516  double Ikj = 0.;
1517  for (int i = 0; i < dim; i++)
1518  {
1519  Ikj += vshape(j, i) * vk[i];
1520  }
1521  I(k, j) = (fabs(Ikj) < 1e-12) ? 0.0 : Ikj;
1522  }
1523  }
1524 }
1525 
1527  const double *nk, const Array<int> &d2n, ElementTransformation &Trans,
1528  DenseMatrix &R) const
1529 {
1530  double pt_data[Geometry::MaxDim];
1531  IntegrationPoint ip;
1532  Vector pt(pt_data, dim);
1533 
1534 #ifdef MFEM_THREAD_SAFE
1536 #endif
1537 
1539  const DenseMatrix &J = Trans.Jacobian();
1540  const double weight = Trans.Weight();
1541  for (int j = 0; j < dof; j++)
1542  {
1543  InvertLinearTrans(Trans, Nodes.IntPoint(j), pt);
1544  ip.Set(pt_data, dim);
1545  if (Geometries.CheckPoint(geom_type, ip)) // do we need an epsilon here?
1546  {
1547  CalcVShape(ip, vshape);
1548  J.MultTranspose(nk+dim*d2n[j], pt_data);
1549  pt /= weight;
1550  for (int k = 0; k < dof; k++)
1551  {
1552  double R_jk = 0.0;
1553  for (int d = 0; d < dim; d++)
1554  {
1555  R_jk += vshape(k,d)*pt_data[d];
1556  }
1557  R(j,k) = R_jk;
1558  }
1559  }
1560  else
1561  {
1562  // Set the whole row to avoid valgrind warnings in R.Threshold().
1563  R.SetRow(j, infinity());
1564  }
1565  }
1566  R.Threshold(1e-12);
1567 }
1568 
1570  const double *tk, const Array<int> &d2t, ElementTransformation &Trans,
1571  DenseMatrix &R) const
1572 {
1573  double pt_data[Geometry::MaxDim];
1574  IntegrationPoint ip;
1575  Vector pt(pt_data, dim);
1576 
1577 #ifdef MFEM_THREAD_SAFE
1579 #endif
1580 
1582  const DenseMatrix &Jinv = Trans.InverseJacobian();
1583  for (int j = 0; j < dof; j++)
1584  {
1585  InvertLinearTrans(Trans, Nodes.IntPoint(j), pt);
1586  ip.Set(pt_data, dim);
1587  if (Geometries.CheckPoint(geom_type, ip)) // do we need an epsilon here?
1588  {
1589  CalcVShape(ip, vshape);
1590  Jinv.Mult(tk+dim*d2t[j], pt_data);
1591  for (int k = 0; k < dof; k++)
1592  {
1593  double R_jk = 0.0;
1594  for (int d = 0; d < dim; d++)
1595  {
1596  R_jk += vshape(k,d)*pt_data[d];
1597  }
1598  R(j,k) = R_jk;
1599  }
1600  }
1601  else
1602  {
1603  // Set the whole row to avoid valgrind warnings in R.Threshold().
1604  R.SetRow(j, infinity());
1605  }
1606  }
1607  R.Threshold(1e-12);
1608 }
1609 
1610 
1611 Poly_1D::Basis::Basis(const int p, const double *nodes, EvalType etype)
1612  : etype(etype), auxiliary_basis(NULL), scale_integrated(false)
1613 {
1614  switch (etype)
1615  {
1616  case ChangeOfBasis:
1617  {
1618  x.SetSize(p + 1);
1619  w.SetSize(p + 1);
1620  DenseMatrix A(p + 1);
1621  for (int i = 0; i <= p; i++)
1622  {
1623  CalcBasis(p, nodes[i], A.GetColumn(i));
1624  }
1625  Ai.Factor(A);
1626  // mfem::out << "Poly_1D::Basis(" << p << ",...) : "; Ai.TestInversion();
1627  break;
1628  }
1629  case Barycentric:
1630  {
1631  x.SetSize(p + 1);
1632  w.SetSize(p + 1);
1633  x = nodes;
1634  w = 1.0;
1635  for (int i = 0; i <= p; i++)
1636  {
1637  for (int j = 0; j < i; j++)
1638  {
1639  double xij = x(i) - x(j);
1640  w(i) *= xij;
1641  w(j) *= -xij;
1642  }
1643  }
1644  for (int i = 0; i <= p; i++)
1645  {
1646  w(i) = 1.0/w(i);
1647  }
1648 
1649 #ifdef MFEM_DEBUG
1650  // Make sure the nodes are increasing
1651  for (int i = 0; i < p; i++)
1652  {
1653  if (x(i) >= x(i+1))
1654  {
1655  mfem_error("Poly_1D::Basis::Basis : nodes are not increasing!");
1656  }
1657  }
1658 #endif
1659  break;
1660  }
1661  case Positive:
1662  x.SetDataAndSize(NULL, p + 1); // use x to store (p + 1)
1663  break;
1664  case Integrated:
1665  auxiliary_basis = new Basis(
1667  u_aux.SetSize(p+2);
1668  d_aux.SetSize(p+2);
1669  d2_aux.SetSize(p+2);
1670  break;
1671  default: break;
1672  }
1673 }
1674 
1675 void Poly_1D::Basis::Eval(const double y, Vector &u) const
1676 {
1677  switch (etype)
1678  {
1679  case ChangeOfBasis:
1680  {
1681  CalcBasis(Ai.Width() - 1, y, x);
1682  Ai.Mult(x, u);
1683  break;
1684  }
1685  case Barycentric:
1686  {
1687  int i, k, p = x.Size() - 1;
1688  double l, lk;
1689 
1690  if (p == 0)
1691  {
1692  u(0) = 1.0;
1693  return;
1694  }
1695 
1696  lk = 1.0;
1697  for (k = 0; k < p; k++)
1698  {
1699  if (y >= (x(k) + x(k+1))/2)
1700  {
1701  lk *= y - x(k);
1702  }
1703  else
1704  {
1705  for (i = k+1; i <= p; i++)
1706  {
1707  lk *= y - x(i);
1708  }
1709  break;
1710  }
1711  }
1712  l = lk * (y - x(k));
1713 
1714  for (i = 0; i < k; i++)
1715  {
1716  u(i) = l * w(i) / (y - x(i));
1717  }
1718  u(k) = lk * w(k);
1719  for (i++; i <= p; i++)
1720  {
1721  u(i) = l * w(i) / (y - x(i));
1722  }
1723  break;
1724  }
1725  case Positive:
1726  CalcBernstein(x.Size() - 1, y, u);
1727  break;
1728  case Integrated:
1729  auxiliary_basis->Eval(y, u_aux, d_aux);
1730  EvalIntegrated(d_aux, u);
1731  break;
1732  default: break;
1733  }
1734 }
1735 
1736 void Poly_1D::Basis::Eval(const double y, Vector &u, Vector &d) const
1737 {
1738  switch (etype)
1739  {
1740  case ChangeOfBasis:
1741  {
1742  CalcBasis(Ai.Width() - 1, y, x, w);
1743  Ai.Mult(x, u);
1744  Ai.Mult(w, d);
1745  break;
1746  }
1747  case Barycentric:
1748  {
1749  int i, k, p = x.Size() - 1;
1750  double l, lp, lk, sk, si;
1751 
1752  if (p == 0)
1753  {
1754  u(0) = 1.0;
1755  d(0) = 0.0;
1756  return;
1757  }
1758 
1759  lk = 1.0;
1760  for (k = 0; k < p; k++)
1761  {
1762  if (y >= (x(k) + x(k+1))/2)
1763  {
1764  lk *= y - x(k);
1765  }
1766  else
1767  {
1768  for (i = k+1; i <= p; i++)
1769  {
1770  lk *= y - x(i);
1771  }
1772  break;
1773  }
1774  }
1775  l = lk * (y - x(k));
1776 
1777  sk = 0.0;
1778  for (i = 0; i < k; i++)
1779  {
1780  si = 1.0/(y - x(i));
1781  sk += si;
1782  u(i) = l * si * w(i);
1783  }
1784  u(k) = lk * w(k);
1785  for (i++; i <= p; i++)
1786  {
1787  si = 1.0/(y - x(i));
1788  sk += si;
1789  u(i) = l * si * w(i);
1790  }
1791  lp = l * sk + lk;
1792 
1793  for (i = 0; i < k; i++)
1794  {
1795  d(i) = (lp * w(i) - u(i))/(y - x(i));
1796  }
1797  d(k) = sk * u(k);
1798  for (i++; i <= p; i++)
1799  {
1800  d(i) = (lp * w(i) - u(i))/(y - x(i));
1801  }
1802  break;
1803  }
1804  case Positive:
1805  CalcBernstein(x.Size() - 1, y, u, d);
1806  break;
1807  case Integrated:
1808  auxiliary_basis->Eval(y, u_aux, d_aux, d2_aux);
1809  EvalIntegrated(d_aux,u);
1810  EvalIntegrated(d2_aux,d);
1811  break;
1812  default: break;
1813  }
1814 }
1815 
1816 void Poly_1D::Basis::Eval(const double y, Vector &u, Vector &d,
1817  Vector &d2) const
1818 {
1819  MFEM_VERIFY(etype == Barycentric,
1820  "Basis::Eval with second order derivatives not implemented for"
1821  " etype = " << etype);
1822  switch (etype)
1823  {
1824  case ChangeOfBasis:
1825  {
1826  CalcBasis(Ai.Width() - 1, y, x, w);
1827  Ai.Mult(x, u);
1828  Ai.Mult(w, d);
1829  // set d2 (not implemented yet)
1830  break;
1831  }
1832  case Barycentric:
1833  {
1834  int i, k, p = x.Size() - 1;
1835  double l, lp, lp2, lk, sk, si, sk2;
1836 
1837  if (p == 0)
1838  {
1839  u(0) = 1.0;
1840  d(0) = 0.0;
1841  d2(0) = 0.0;
1842  return;
1843  }
1844 
1845  lk = 1.0;
1846  for (k = 0; k < p; k++)
1847  {
1848  if (y >= (x(k) + x(k+1))/2)
1849  {
1850  lk *= y - x(k);
1851  }
1852  else
1853  {
1854  for (i = k+1; i <= p; i++)
1855  {
1856  lk *= y - x(i);
1857  }
1858  break;
1859  }
1860  }
1861  l = lk * (y - x(k));
1862 
1863  sk = 0.0;
1864  sk2 = 0.0;
1865  for (i = 0; i < k; i++)
1866  {
1867  si = 1.0/(y - x(i));
1868  sk += si;
1869  sk2 -= si * si;
1870  u(i) = l * si * w(i);
1871  }
1872  u(k) = lk * w(k);
1873  for (i++; i <= p; i++)
1874  {
1875  si = 1.0/(y - x(i));
1876  sk += si;
1877  sk2 -= si * si;
1878  u(i) = l * si * w(i);
1879  }
1880  lp = l * sk + lk;
1881  lp2 = lp * sk + l * sk2 + sk * lk;
1882 
1883  for (i = 0; i < k; i++)
1884  {
1885  d(i) = (lp * w(i) - u(i))/(y - x(i));
1886  d2(i) = (lp2 * w(i) - 2 * d(i))/(y - x(i));
1887  }
1888  d(k) = sk * u(k);
1889  d2(k) = sk2 * u(k) + sk * d(k);
1890  for (i++; i <= p; i++)
1891  {
1892  d(i) = (lp * w(i) - u(i))/(y - x(i));
1893  d2(i) = (lp2 * w(i) - 2 * d(i))/(y - x(i));
1894  }
1895  break;
1896  }
1897  case Positive:
1898  CalcBernstein(x.Size() - 1, y, u, d);
1899  break;
1900  case Integrated:
1901  MFEM_ABORT("Integrated basis must be evaluated with EvalIntegrated");
1902  break;
1903  default: break;
1904  }
1905 }
1906 
1907 void Poly_1D::Basis::EvalIntegrated(const Vector &d_aux_, Vector &u) const
1908 {
1909  MFEM_VERIFY(etype == Integrated,
1910  "EvalIntegrated is only valid for Integrated basis type");
1911  int p = d_aux_.Size() - 1;
1912  // See Gerritsma, M. (2010). "Edge functions for spectral element methods",
1913  // in Lecture Notes in Computational Science and Engineering, 199--207.
1914  u[0] = -d_aux_[0];
1915  for (int j=1; j<p; ++j)
1916  {
1917  u[j] = u[j-1] - d_aux_[j];
1918  }
1919  // If scale_integrated is true, the degrees of freedom represent mean values,
1920  // otherwise they represent subcell integrals. Generally, scale_integrated
1921  // should be true for MapType::VALUE, and false for other map types.
1922  if (scale_integrated)
1923  {
1924  Vector &aux_nodes = auxiliary_basis->x;
1925  for (int j=0; j<aux_nodes.Size()-1; ++j)
1926  {
1927  u[j] *= aux_nodes[j+1] - aux_nodes[j];
1928  }
1929  }
1930 }
1931 
1932 void Poly_1D::Basis::ScaleIntegrated(bool scale_integrated_)
1933 {
1934  scale_integrated = scale_integrated_;
1935 }
1936 
1938 {
1939  delete auxiliary_basis;
1940 }
1941 
1942 const int *Poly_1D::Binom(const int p)
1943 {
1944  if (binom.NumCols() <= p)
1945  {
1946  binom.SetSize(p + 1, p + 1);
1947  for (int i = 0; i <= p; i++)
1948  {
1949  binom(i,0) = binom(i,i) = 1;
1950  for (int j = 1; j < i; j++)
1951  {
1952  binom(i,j) = binom(i-1,j) + binom(i-1,j-1);
1953  }
1954  }
1955  }
1956  return binom[p];
1957 }
1958 
1959 void Poly_1D::ChebyshevPoints(const int p, double *x)
1960 {
1961  for (int i = 0; i <= p; i++)
1962  {
1963  // x[i] = 0.5*(1. + cos(M_PI*(p - i + 0.5)/(p + 1)));
1964  double s = sin(M_PI_2*(i + 0.5)/(p + 1));
1965  x[i] = s*s;
1966  }
1967 }
1968 
1969 void Poly_1D::CalcMono(const int p, const double x, double *u)
1970 {
1971  double xn;
1972  u[0] = xn = 1.;
1973  for (int n = 1; n <= p; n++)
1974  {
1975  u[n] = (xn *= x);
1976  }
1977 }
1978 
1979 void Poly_1D::CalcMono(const int p, const double x, double *u, double *d)
1980 {
1981  double xn;
1982  u[0] = xn = 1.;
1983  d[0] = 0.;
1984  for (int n = 1; n <= p; n++)
1985  {
1986  d[n] = n * xn;
1987  u[n] = (xn *= x);
1988  }
1989 }
1990 
1991 void Poly_1D::CalcBinomTerms(const int p, const double x, const double y,
1992  double *u)
1993 {
1994  if (p == 0)
1995  {
1996  u[0] = 1.;
1997  }
1998  else
1999  {
2000  int i;
2001  const int *b = Binom(p);
2002  double z = x;
2003 
2004  for (i = 1; i < p; i++)
2005  {
2006  u[i] = b[i]*z;
2007  z *= x;
2008  }
2009  u[p] = z;
2010  z = y;
2011  for (i--; i > 0; i--)
2012  {
2013  u[i] *= z;
2014  z *= y;
2015  }
2016  u[0] = z;
2017  }
2018 }
2019 
2020 void Poly_1D::CalcBinomTerms(const int p, const double x, const double y,
2021  double *u, double *d)
2022 {
2023  if (p == 0)
2024  {
2025  u[0] = 1.;
2026  d[0] = 0.;
2027  }
2028  else
2029  {
2030  int i;
2031  const int *b = Binom(p);
2032  const double xpy = x + y, ptx = p*x;
2033  double z = 1.;
2034 
2035  for (i = 1; i < p; i++)
2036  {
2037  d[i] = b[i]*z*(i*xpy - ptx);
2038  z *= x;
2039  u[i] = b[i]*z;
2040  }
2041  d[p] = p*z;
2042  u[p] = z*x;
2043  z = 1.;
2044  for (i--; i > 0; i--)
2045  {
2046  d[i] *= z;
2047  z *= y;
2048  u[i] *= z;
2049  }
2050  d[0] = -p*z;
2051  u[0] = z*y;
2052  }
2053 }
2054 
2055 void Poly_1D::CalcDBinomTerms(const int p, const double x, const double y,
2056  double *d)
2057 {
2058  if (p == 0)
2059  {
2060  d[0] = 0.;
2061  }
2062  else
2063  {
2064  int i;
2065  const int *b = Binom(p);
2066  const double xpy = x + y, ptx = p*x;
2067  double z = 1.;
2068 
2069  for (i = 1; i < p; i++)
2070  {
2071  d[i] = b[i]*z*(i*xpy - ptx);
2072  z *= x;
2073  }
2074  d[p] = p*z;
2075  z = 1.;
2076  for (i--; i > 0; i--)
2077  {
2078  d[i] *= z;
2079  z *= y;
2080  }
2081  d[0] = -p*z;
2082  }
2083 }
2084 
2085 void Poly_1D::CalcLegendre(const int p, const double x, double *u)
2086 {
2087  // use the recursive definition for [-1,1]:
2088  // (n+1)*P_{n+1}(z) = (2*n+1)*z*P_n(z)-n*P_{n-1}(z)
2089  double z;
2090  u[0] = 1.;
2091  if (p == 0) { return; }
2092  u[1] = z = 2.*x - 1.;
2093  for (int n = 1; n < p; n++)
2094  {
2095  u[n+1] = ((2*n + 1)*z*u[n] - n*u[n-1])/(n + 1);
2096  }
2097 }
2098 
2099 void Poly_1D::CalcLegendre(const int p, const double x, double *u, double *d)
2100 {
2101  // use the recursive definition for [-1,1]:
2102  // (n+1)*P_{n+1}(z) = (2*n+1)*z*P_n(z)-n*P_{n-1}(z)
2103  // for the derivative use, z in [-1,1]:
2104  // P'_{n+1}(z) = (2*n+1)*P_n(z)+P'_{n-1}(z)
2105  double z;
2106  u[0] = 1.;
2107  d[0] = 0.;
2108  if (p == 0) { return; }
2109  u[1] = z = 2.*x - 1.;
2110  d[1] = 2.;
2111  for (int n = 1; n < p; n++)
2112  {
2113  u[n+1] = ((2*n + 1)*z*u[n] - n*u[n-1])/(n + 1);
2114  d[n+1] = (4*n + 2)*u[n] + d[n-1];
2115  }
2116 }
2117 
2118 void Poly_1D::CalcChebyshev(const int p, const double x, double *u)
2119 {
2120  // recursive definition, z in [-1,1]
2121  // T_0(z) = 1, T_1(z) = z
2122  // T_{n+1}(z) = 2*z*T_n(z) - T_{n-1}(z)
2123  double z;
2124  u[0] = 1.;
2125  if (p == 0) { return; }
2126  u[1] = z = 2.*x - 1.;
2127  for (int n = 1; n < p; n++)
2128  {
2129  u[n+1] = 2*z*u[n] - u[n-1];
2130  }
2131 }
2132 
2133 void Poly_1D::CalcChebyshev(const int p, const double x, double *u, double *d)
2134 {
2135  // recursive definition, z in [-1,1]
2136  // T_0(z) = 1, T_1(z) = z
2137  // T_{n+1}(z) = 2*z*T_n(z) - T_{n-1}(z)
2138  // T'_n(z) = n*U_{n-1}(z)
2139  // U_0(z) = 1 U_1(z) = 2*z
2140  // U_{n+1}(z) = 2*z*U_n(z) - U_{n-1}(z)
2141  // U_n(z) = z*U_{n-1}(z) + T_n(z) = z*T'_n(z)/n + T_n(z)
2142  // T'_{n+1}(z) = (n + 1)*(z*T'_n(z)/n + T_n(z))
2143  double z;
2144  u[0] = 1.;
2145  d[0] = 0.;
2146  if (p == 0) { return; }
2147  u[1] = z = 2.*x - 1.;
2148  d[1] = 2.;
2149  for (int n = 1; n < p; n++)
2150  {
2151  u[n+1] = 2*z*u[n] - u[n-1];
2152  d[n+1] = (n + 1)*(z*d[n]/n + 2*u[n]);
2153  }
2154 }
2155 
2156 void Poly_1D::CalcChebyshev(const int p, const double x, double *u, double *d,
2157  double *dd)
2158 {
2159  // recursive definition, z in [-1,1]
2160  // T_0(z) = 1, T_1(z) = z
2161  // T_{n+1}(z) = 2*z*T_n(z) - T_{n-1}(z)
2162  // T'_n(z) = n*U_{n-1}(z)
2163  // U_0(z) = 1 U_1(z) = 2*z
2164  // U_{n+1}(z) = 2*z*U_n(z) - U_{n-1}(z)
2165  // U_n(z) = z*U_{n-1}(z) + T_n(z) = z*T'_n(z)/n + T_n(z)
2166  // T'_{n+1}(z) = (n + 1)*(z*T'_n(z)/n + T_n(z))
2167  // T''_{n+1}(z) = (n + 1)*(2*(n + 1)*T'_n(z) + z*T''_n(z)) / n
2168  double z;
2169  u[0] = 1.;
2170  d[0] = 0.;
2171  dd[0]= 0.;
2172  if (p == 0) { return; }
2173  u[1] = z = 2.*x - 1.;
2174  d[1] = 2.;
2175  dd[1] = 0;
2176  for (int n = 1; n < p; n++)
2177  {
2178  u[n+1] = 2*z*u[n] - u[n-1];
2179  d[n+1] = (n + 1)*(z*d[n]/n + 2*u[n]);
2180  dd[n+1] = (n + 1)*(2.*(n + 1)*d[n] + z*dd[n])/n;
2181  }
2182 }
2183 
2184 const double *Poly_1D::GetPoints(const int p, const int btype)
2185 {
2186  BasisType::Check(btype);
2187  const int qtype = BasisType::GetQuadrature1D(btype);
2188 
2189  if (qtype == Quadrature1D::Invalid) { return NULL; }
2190 
2191  if (points_container.find(btype) == points_container.end())
2192  {
2193  points_container[btype] = new Array<double*>(h_mt);
2194  }
2195  Array<double*> &pts = *points_container[btype];
2196  if (pts.Size() <= p)
2197  {
2198  pts.SetSize(p + 1, NULL);
2199  }
2200  if (pts[p] == NULL)
2201  {
2202  pts[p] = new double[p + 1];
2203  quad_func.GivePolyPoints(p+1, pts[p], qtype);
2204  }
2205  return pts[p];
2206 }
2207 
2208 Poly_1D::Basis &Poly_1D::GetBasis(const int p, const int btype)
2209 {
2210  BasisType::Check(btype);
2211 
2212  if ( bases_container.find(btype) == bases_container.end() )
2213  {
2214  // we haven't been asked for basis or points of this type yet
2215  bases_container[btype] = new Array<Basis*>(h_mt);
2216  }
2217  Array<Basis*> &bases = *bases_container[btype];
2218  if (bases.Size() <= p)
2219  {
2220  bases.SetSize(p + 1, NULL);
2221  }
2222  if (bases[p] == NULL)
2223  {
2224  EvalType etype;
2225  if (btype == BasisType::Positive) { etype = Positive; }
2226  else if (btype == BasisType::IntegratedGLL) { etype = Integrated; }
2227  else { etype = Barycentric; }
2228  bases[p] = new Basis(p, GetPoints(p, btype), etype);
2229  }
2230  return *bases[p];
2231 }
2232 
2234 {
2235  for (PointsMap::iterator it = points_container.begin();
2236  it != points_container.end() ; ++it)
2237  {
2238  Array<double*>& pts = *it->second;
2239  for ( int i = 0 ; i < pts.Size() ; ++i )
2240  {
2241  delete [] pts[i];
2242  }
2243  delete it->second;
2244  }
2245 
2246  for (BasisMap::iterator it = bases_container.begin();
2247  it != bases_container.end() ; ++it)
2248  {
2249  Array<Basis*>& bases = *it->second;
2250  for ( int i = 0 ; i < bases.Size() ; ++i )
2251  {
2252  delete bases[i];
2253  }
2254  delete it->second;
2255  }
2256 }
2257 
2258 
2259 TensorBasisElement::TensorBasisElement(const int dims, const int p,
2260  const int btype, const DofMapType dmtype)
2261  : b_type(btype),
2262  basis1d(poly1d.GetBasis(p, b_type))
2263 {
2264  if (dmtype == H1_DOF_MAP || dmtype == Sr_DOF_MAP)
2265  {
2266  switch (dims)
2267  {
2268  case 1:
2269  {
2270  dof_map.SetSize(p + 1);
2271  dof_map[0] = 0;
2272  dof_map[p] = 1;
2273  for (int i = 1; i < p; i++)
2274  {
2275  dof_map[i] = i+1;
2276  }
2277  break;
2278  }
2279  case 2:
2280  {
2281  const int p1 = p + 1;
2282  dof_map.SetSize(p1*p1);
2283 
2284  // vertices
2285  dof_map[0 + 0*p1] = 0;
2286  dof_map[p + 0*p1] = 1;
2287  dof_map[p + p*p1] = 2;
2288  dof_map[0 + p*p1] = 3;
2289 
2290  // edges
2291  int o = 4;
2292  for (int i = 1; i < p; i++)
2293  {
2294  dof_map[i + 0*p1] = o++;
2295  }
2296  for (int i = 1; i < p; i++)
2297  {
2298  dof_map[p + i*p1] = o++;
2299  }
2300  for (int i = 1; i < p; i++)
2301  {
2302  dof_map[(p-i) + p*p1] = o++;
2303  }
2304  for (int i = 1; i < p; i++)
2305  {
2306  dof_map[0 + (p-i)*p1] = o++;
2307  }
2308 
2309  // interior
2310  for (int j = 1; j < p; j++)
2311  {
2312  for (int i = 1; i < p; i++)
2313  {
2314  dof_map[i + j*p1] = o++;
2315  }
2316  }
2317  break;
2318  }
2319  case 3:
2320  {
2321  const int p1 = p + 1;
2322  dof_map.SetSize(p1*p1*p1);
2323 
2324  // vertices
2325  dof_map[0 + (0 + 0*p1)*p1] = 0;
2326  dof_map[p + (0 + 0*p1)*p1] = 1;
2327  dof_map[p + (p + 0*p1)*p1] = 2;
2328  dof_map[0 + (p + 0*p1)*p1] = 3;
2329  dof_map[0 + (0 + p*p1)*p1] = 4;
2330  dof_map[p + (0 + p*p1)*p1] = 5;
2331  dof_map[p + (p + p*p1)*p1] = 6;
2332  dof_map[0 + (p + p*p1)*p1] = 7;
2333 
2334  // edges (see Hexahedron::edges in mesh/hexahedron.cpp).
2335  // edges (see Constants<Geometry::CUBE>::Edges in fem/geom.cpp).
2336  int o = 8;
2337  for (int i = 1; i < p; i++)
2338  {
2339  dof_map[i + (0 + 0*p1)*p1] = o++; // (0,1)
2340  }
2341  for (int i = 1; i < p; i++)
2342  {
2343  dof_map[p + (i + 0*p1)*p1] = o++; // (1,2)
2344  }
2345  for (int i = 1; i < p; i++)
2346  {
2347  dof_map[i + (p + 0*p1)*p1] = o++; // (3,2)
2348  }
2349  for (int i = 1; i < p; i++)
2350  {
2351  dof_map[0 + (i + 0*p1)*p1] = o++; // (0,3)
2352  }
2353  for (int i = 1; i < p; i++)
2354  {
2355  dof_map[i + (0 + p*p1)*p1] = o++; // (4,5)
2356  }
2357  for (int i = 1; i < p; i++)
2358  {
2359  dof_map[p + (i + p*p1)*p1] = o++; // (5,6)
2360  }
2361  for (int i = 1; i < p; i++)
2362  {
2363  dof_map[i + (p + p*p1)*p1] = o++; // (7,6)
2364  }
2365  for (int i = 1; i < p; i++)
2366  {
2367  dof_map[0 + (i + p*p1)*p1] = o++; // (4,7)
2368  }
2369  for (int i = 1; i < p; i++)
2370  {
2371  dof_map[0 + (0 + i*p1)*p1] = o++; // (0,4)
2372  }
2373  for (int i = 1; i < p; i++)
2374  {
2375  dof_map[p + (0 + i*p1)*p1] = o++; // (1,5)
2376  }
2377  for (int i = 1; i < p; i++)
2378  {
2379  dof_map[p + (p + i*p1)*p1] = o++; // (2,6)
2380  }
2381  for (int i = 1; i < p; i++)
2382  {
2383  dof_map[0 + (p + i*p1)*p1] = o++; // (3,7)
2384  }
2385 
2386  // faces (see Mesh::GenerateFaces in mesh/mesh.cpp)
2387  for (int j = 1; j < p; j++)
2388  {
2389  for (int i = 1; i < p; i++)
2390  {
2391  dof_map[i + ((p-j) + 0*p1)*p1] = o++; // (3,2,1,0)
2392  }
2393  }
2394  for (int j = 1; j < p; j++)
2395  {
2396  for (int i = 1; i < p; i++)
2397  {
2398  dof_map[i + (0 + j*p1)*p1] = o++; // (0,1,5,4)
2399  }
2400  }
2401  for (int j = 1; j < p; j++)
2402  {
2403  for (int i = 1; i < p; i++)
2404  {
2405  dof_map[p + (i + j*p1)*p1] = o++; // (1,2,6,5)
2406  }
2407  }
2408  for (int j = 1; j < p; j++)
2409  {
2410  for (int i = 1; i < p; i++)
2411  {
2412  dof_map[(p-i) + (p + j*p1)*p1] = o++; // (2,3,7,6)
2413  }
2414  }
2415  for (int j = 1; j < p; j++)
2416  {
2417  for (int i = 1; i < p; i++)
2418  {
2419  dof_map[0 + ((p-i) + j*p1)*p1] = o++; // (3,0,4,7)
2420  }
2421  }
2422  for (int j = 1; j < p; j++)
2423  {
2424  for (int i = 1; i < p; i++)
2425  {
2426  dof_map[i + (j + p*p1)*p1] = o++; // (4,5,6,7)
2427  }
2428  }
2429 
2430  // interior
2431  for (int k = 1; k < p; k++)
2432  {
2433  for (int j = 1; j < p; j++)
2434  {
2435  for (int i = 1; i < p; i++)
2436  {
2437  dof_map[i + (j + k*p1)*p1] = o++;
2438  }
2439  }
2440  }
2441  break;
2442  }
2443  default:
2444  MFEM_ABORT("invalid dimension: " << dims);
2445  break;
2446  }
2447  }
2448  else if (dmtype == L2_DOF_MAP)
2449  {
2450  // leave dof_map empty, indicating that the dofs are ordered
2451  // lexicographically, i.e. the dof_map is identity
2452  }
2453  else
2454  {
2455  MFEM_ABORT("invalid DofMapType: " << dmtype);
2456  }
2457 }
2458 
2460  const FiniteElement &fe, const IntegrationRule &ir,
2461  DofToQuad::Mode mode, const Poly_1D::Basis &basis, bool closed,
2462  Array<DofToQuad*> &dof2quad_array)
2463 {
2464  MFEM_VERIFY(mode == DofToQuad::TENSOR, "invalid mode requested");
2465 
2466  for (int i = 0; i < dof2quad_array.Size(); i++)
2467  {
2468  const DofToQuad &d2q = *dof2quad_array[i];
2469  if (d2q.IntRule == &ir && d2q.mode == mode) { return d2q; }
2470  }
2471 
2472  DofToQuad *d2q = new DofToQuad;
2473  const int ndof = closed ? fe.GetOrder() + 1 : fe.GetOrder();
2474  const int nqpt = (int)floor(pow(ir.GetNPoints(), 1.0/fe.GetDim()) + 0.5);
2475  d2q->FE = &fe;
2476  d2q->IntRule = &ir;
2477  d2q->mode = mode;
2478  d2q->ndof = ndof;
2479  d2q->nqpt = nqpt;
2480  d2q->B.SetSize(nqpt*ndof);
2481  d2q->Bt.SetSize(ndof*nqpt);
2482  d2q->G.SetSize(nqpt*ndof);
2483  d2q->Gt.SetSize(ndof*nqpt);
2484  Vector val(ndof), grad(ndof);
2485  for (int i = 0; i < nqpt; i++)
2486  {
2487  // The first 'nqpt' points in 'ir' have the same x-coordinates as those
2488  // of the 1D rule.
2489  basis.Eval(ir.IntPoint(i).x, val, grad);
2490  for (int j = 0; j < ndof; j++)
2491  {
2492  d2q->B[i+nqpt*j] = d2q->Bt[j+ndof*i] = val(j);
2493  d2q->G[i+nqpt*j] = d2q->Gt[j+ndof*i] = grad(j);
2494  }
2495  }
2496  dof2quad_array.Append(d2q);
2497  return *d2q;
2498 }
2499 
2501  const int p,
2502  const int btype,
2503  const DofMapType dmtype)
2504  : NodalFiniteElement(dims, GetTensorProductGeometry(dims), Pow(p + 1, dims),
2505  p, dims > 1 ? FunctionSpace::Qk : FunctionSpace::Pk),
2506  TensorBasisElement(dims, p, btype, dmtype)
2507 {
2509 }
2510 
2512 {
2514  // If we are using the "integrated" basis, the basis functions should be
2515  // scaled for MapType::VALUE, and not scaled for MapType::INTEGRAL. This
2516  // ensures spectral equivalence of the mass matrix with its low-order-refined
2517  // counterpart (cf. LORDiscretization)
2518  if (basis1d.IsIntegratedType())
2519  {
2521  }
2522 }
2523 
2525  Array<int> &face_map) const
2526 {
2527  internal::GetTensorFaceMap(dim, order, face_id, face_map);
2528 }
2529 
2531  const int d,
2532  const int p,
2533  const int cbtype,
2534  const int obtype,
2535  const int M,
2536  const DofMapType dmtype)
2537  : VectorFiniteElement(dims, GetTensorProductGeometry(dims), d,
2538  p, M, FunctionSpace::Qk),
2539  TensorBasisElement(dims, p, VerifyNodal(VerifyClosed(cbtype)), dmtype),
2540  obasis1d(poly1d.GetBasis(p - 1, VerifyOpen(obtype)))
2541 {
2542  MFEM_VERIFY(dims > 1, "Constructor for VectorTensorFiniteElement with both "
2543  "open and closed bases is not valid for 1D elements.");
2544 }
2545 
2547  const int d,
2548  const int p,
2549  const int obtype,
2550  const int M,
2551  const DofMapType dmtype)
2552  : VectorFiniteElement(dims, GetTensorProductGeometry(dims), d,
2553  p, M, FunctionSpace::Pk),
2554  TensorBasisElement(dims, p, VerifyOpen(obtype), dmtype),
2555  obasis1d(poly1d.GetBasis(p, VerifyOpen(obtype)))
2556 {
2557  MFEM_VERIFY(dims == 1, "Constructor for VectorTensorFiniteElement without "
2558  "closed basis is only valid for 1D elements.");
2559 }
2560 
2562 {
2563  for (int i = 0; i < dof2quad_array_open.Size(); i++)
2564  {
2565  delete dof2quad_array_open[i];
2566  }
2567 }
2568 
2569 }
Abstract class for all finite elements.
Definition: fe_base.hpp:233
int GetHeight() const
Get the height of the matrix.
void Set(const double *p, const int dim)
Definition: intrules.hpp:43
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
Definition: densemat.cpp:2761
virtual void CalcPhysLinLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Definition: fe_base.cpp:244
const DenseMatrix & AdjugateJacobian()
Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoin...
Definition: eltrans.hpp:135
void trans(const Vector &u, Vector &x)
Definition: ex27.cpp:412
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:253
bool IsIntegratedType() const
Returns true if the basis is "integrated", false otherwise.
Definition: fe_base.hpp:1023
const class FiniteElement * FE
The FiniteElement that created and owns this object.
Definition: fe_base.hpp:141
Array< int > lex_ordering
Definition: fe_base.hpp:711
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
static int Check(int b_type)
If the input does not represents a valid BasisType, abort with an error; otherwise return the input...
Definition: fe_base.hpp:45
EvalType
One-dimensional basis evaluation type.
Definition: fe_base.hpp:967
void ScaleIntegrated(bool scale_integrated_)
Set whether the "integrated" basis should be scaled by the subcell sizes. Has no effect for non-integ...
Definition: fe_base.cpp:1932
NodalTensorFiniteElement(const int dims, const int p, const int btype, const DofMapType dmtype)
Definition: fe_base.cpp:2500
virtual void ProjectFromNodes(Vector &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector of values at the finite element nodes and a transformation, compute its projection (ap...
Definition: fe_base.cpp:138
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:980
Basis(const int p, const double *nodes, EvalType etype=Barycentric)
Create a nodal or positive (Bernstein) basis of degree p.
Definition: fe_base.cpp:1611
void LocalL2Projection_RT(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1353
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition: fe_base.hpp:161
Base class for vector Coefficients that optionally depend on time and space.
void SetRow(int r, const double *row)
Definition: densemat.cpp:2177
void InvRightScaling(const Vector &s)
InvRightScaling: this = this * diag(1./s);.
Definition: densemat.cpp:394
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in reference space at the given p...
Definition: fe_base.cpp:65
void ProjectMatrixCoefficient_RT(const double *nk, const Array< int > &d2n, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Project the rows of the matrix coefficient in an RT space.
Definition: fe_base.cpp:1017
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 ...
void ProjectGrad_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe_base.cpp:1120
static void GivePolyPoints(const int np, double *pts, const int type)
Definition: intrules.cpp:703
void Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
Definition: densemat.cpp:3943
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const override
Compute the discrete gradient matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:830
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &Hessian) const
Evaluate the Hessians of all shape functions of a scalar finite element in reference space at the giv...
Definition: fe_base.cpp:101
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_base.cpp:107
FiniteElement(int D, Geometry::Type G, int Do, int O, int F=FunctionSpace::Pk)
Construct FiniteElement with given.
Definition: fe_base.cpp:23
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
virtual void Eval(DenseMatrix &K, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the matrix coefficient in the element described by T at the point ip, storing the result in ...
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
int dim
Dimension of reference space.
Definition: fe_base.hpp:236
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
TensorBasisElement(const int dims, const int p, const int btype, const DofMapType dmtype)
Definition: fe_base.cpp:2259
virtual void CalcPhysLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Evaluate the Laplacian of all shape functions of a scalar finite element in reference space at the gi...
Definition: fe_base.cpp:203
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Definition: fe_base.hpp:173
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Array< double > Gt
Transpose of G.
Definition: fe_base.hpp:212
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe_base.hpp:169
void Eval(const double x, Vector &u) const
Evaluate the basis functions at point x in [0,1].
Definition: fe_base.cpp:1675
static int GetQuadrature1D(int b_type)
Get the corresponding Quadrature1D constant, when that makes sense; otherwise return Quadrature1D::In...
Definition: fe_base.hpp:61
VectorFiniteElement(int D, Geometry::Type G, int Do, int O, int M, int F=FunctionSpace::Pk)
Definition: fe_base.cpp:891
const DenseMatrix & InverseJacobian()
Return the inverse of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:145
DenseMatrix vshape
Definition: fe_base.hpp:248
void CalcVShape_RT(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe_base.cpp:957
static const DofToQuad & GetTensorDofToQuad(const FiniteElement &fe, const IntegrationRule &ir, DofToQuad::Mode mode, const Poly_1D::Basis &basis, bool closed, Array< DofToQuad *> &dof2quad_array)
Definition: fe_base.cpp:2459
void Factor()
Factor the current DenseMatrix, *a.
Definition: densemat.cpp:3897
STL namespace.
Geometry::Type geom_type
Geometry::Type of the reference element.
Definition: fe_base.hpp:239
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:98
Intermediate class for finite elements whose basis functions return vector values.
Definition: fe_base.hpp:788
static void CalcLegendre(const int p, const double x, double *u)
Definition: fe_base.cpp:2085
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
Definition: geom.hpp:71
void ProjectCurl_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe_base.cpp:1147
Use change of basis, O(p^2) Evals.
Definition: fe_base.hpp:969
int GetWidth() const
Get the width of the matrix.
Geometry Geometries
Definition: fe.cpp:49
void ProjectMatrixCoefficient_ND(const double *tk, const Array< int > &d2t, MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Project the rows of the matrix coefficient in an ND space.
Definition: fe_base.cpp:1233
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:480
int cdim
Dimension of curl for vector-valued basis functions.
Definition: fe_base.hpp:238
int vdim
Vector dimension of vector-valued basis functions.
Definition: fe_base.hpp:237
static void CalcBasis(const int p, const double x, double *u)
Evaluate the values of a hierarchical 1D basis at point x hierarchical = k-th basis function is degre...
Definition: fe_base.hpp:1087
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:1991
const double * GetPoints(const int p, const int btype)
Get the coordinates of the points of the given BasisType, btype.
Definition: fe_base.cpp:2184
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_base.cpp:126
Class for standard nodal finite elements.
Definition: fe_base.hpp:708
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
VectorTensorFiniteElement(const int dims, const int d, const int p, const int cbtype, const int obtype, const int M, const DofMapType dmtype)
Definition: fe_base.cpp:2530
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition: fe_base.cpp:365
Class for finite elements with basis functions that return scalar values.
Definition: fe_base.hpp:649
void SetSize(int m, int n)
Definition: array.hpp:375
void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const override
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
Definition: fe_base.cpp:749
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:256
virtual void SetMapType(int M)
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition: fe_base.hpp:675
void LocalRestriction_ND(const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &R) const
Definition: fe_base.cpp:1569
virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const
Get the dofs associated with the given face. *dofs is set to an internal array of the local dofc on t...
Definition: fe_base.cpp:96
double * GetData() const
Returns the matrix data array.
Definition: densemat.hpp:115
void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const override
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
Definition: fe_base.cpp:678
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
Poly_1D::Basis & basis1d
Definition: fe_base.hpp:1188
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:214
Fast evaluation of Bernstein polynomials.
Definition: fe_base.hpp:971
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
static void ChebyshevPoints(const int p, double *x)
Compute the points for the Chebyshev polynomials of order p and place them in the already allocated x...
Definition: fe_base.cpp:1959
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
Definition: densemat.cpp:3191
double b
Definition: lissajous.cpp:42
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
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_base.cpp:150
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:335
Mode mode
Describes the contents of the B, Bt, G, and Gt arrays, see Mode.
Definition: fe_base.hpp:165
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
void ProjectCurl_2D(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe_base.cpp:638
Bernstein polynomials.
Definition: fe_base.hpp:33
void LocalL2Projection_ND(const VectorFiniteElement &cfe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1441
void Get(double *p, const int dim) const
Definition: intrules.hpp:57
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:678
virtual void CalcPhysCurlShape(ElementTransformation &Trans, DenseMatrix &curl_shape) const
Evaluate the curl of all shape functions of a vector finite element in physical space at the point de...
Definition: fe_base.cpp:71
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Definition: densemat.cpp:3213
int orders[Geometry::MaxDim]
Anisotropic orders.
Definition: fe_base.hpp:245
void GetColumn(int c, Vector &col) const
Definition: densemat.cpp:1330
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Compute the discrete gradient matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:161
virtual ~FiniteElement()
Deconstruct the FiniteElement.
Definition: fe_base.cpp:500
double FNorm2() const
Compute the square of the Frobenius norm of the matrix.
Definition: densemat.hpp:269
Use barycentric Lagrangian interpolation, O(p) Evals.
Definition: fe_base.hpp:970
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
Definition: fe_base.cpp:119
int GetVDim()
Returns dimension of the vector.
void CalcPhysDivShape(ElementTransformation &Trans, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in physical space at the po...
Definition: fe_base.cpp:58
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:172
Basis & GetBasis(const int p, const int btype)
Get a Poly_1D::Basis object of the given degree and BasisType, btype.
Definition: fe_base.cpp:2208
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Compute the discrete curl matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:168
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
Definition: fe_base.cpp:192
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Compute the discrete divergence matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:175
void EvalIntegrated(const Vector &d, Vector &i) const
Evaluate the "integrated" basis type using pre-computed closed basis derivatives. ...
Definition: fe_base.cpp:1907
static const int MaxDim
Definition: geom.hpp:43
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:206
void Threshold(double eps)
Replace small entries, abs(a_ij) <= eps, with zero.
Definition: densemat.cpp:2207
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:2586
void Set3(const double x1, const double x2, const double x3)
Definition: intrules.hpp:76
IntegrationRule Nodes
Definition: fe_base.hpp:246
void LocalInterpolation_RT(const VectorFiniteElement &cfe, const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1400
Array< double > Bt
Transpose of B.
Definition: fe_base.hpp:190
virtual void ProjectMatrixCoefficient(MatrixCoefficient &mc, ElementTransformation &T, Vector &dofs) const
Given a matrix coefficient and a transformation, compute an approximation ("projection") in the local...
Definition: fe_base.cpp:144
Integrated indicator functions (cf. Gerritsma)
Definition: fe_base.hpp:972
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
void Project_RT(const double *nk, const Array< int > &d2n, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Project a vector coefficient onto the RT basis functions.
Definition: fe_base.cpp:980
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
void LocalInterpolation_ND(const VectorFiniteElement &cfe, const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe_base.cpp:1487
Base class for Matrix Coefficients that optionally depend on time and space.
void ScalarLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I, const ScalarFiniteElement &fine_fe) const
Get matrix I "Interpolation" defined through local L2-projection in the space defined by the fine_fe...
Definition: fe_base.cpp:548
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
Definition: densemat.cpp:2699
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:1942
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:312
Poly_1D poly1d
Definition: fe.cpp:28
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe_base.hpp:136
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:154
void SetMapType(const int map_type_) override
Set the FiniteElement::MapType of the element to either VALUE or INTEGRAL. Also sets the FiniteElemen...
Definition: fe_base.cpp:2511
Integrated GLL indicator functions.
Definition: fe_base.hpp:39
Array< double > B
Basis functions evaluated at quadrature points.
Definition: fe_base.hpp:184
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Definition: fe_base.cpp:182
Class for integration point with weight.
Definition: intrules.hpp:31
Mode
Type of data stored in the arrays B, Bt, G, and Gt.
Definition: fe_base.hpp:149
Full multidimensional representation which does not use tensor product structure. The ordering of the...
Definition: fe_base.hpp:154
void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const override
Compute the discrete divergence matrix from the given FiniteElement onto &#39;this&#39; FiniteElement. The ElementTransformation is included to support cases when the matrix depends on it.
Definition: fe_base.cpp:859
void ProjectCurl_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe_base.cpp:1185
int dof
Number of degrees of freedom.
Definition: fe_base.hpp:243
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
Definition: fe_base.hpp:314
Array< DofToQuad * > dof2quad_array
Container for all DofToQuad objects created by the FiniteElement.
Definition: fe_base.hpp:253
int NumCols() const
Definition: array.hpp:378
Implements CalcDivShape methods.
Definition: fe_base.hpp:296
Array< double > G
Gradients/divergences/curls of basis functions evaluated at quadrature points.
Definition: fe_base.hpp:205
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe_base.cpp:40
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_base.cpp:2524
Class for evaluating 1D nodal, positive (Bernstein), or integrated (Gerritsma) bases.
Definition: fe_base.hpp:978
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition: fe_base.hpp:349
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition: vector.hpp:44
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
void InvertLinearTrans(ElementTransformation &trans, const IntegrationPoint &pt, Vector &x)
Definition: fe_base.cpp:662
Vector data type.
Definition: vector.hpp:58
virtual void GetFaceMap(const int face_id, Array< int > &face_map) const
Return the mapping from lexicographic face DOFs to lexicographic element DOFs for the given local fac...
Definition: fe_base.cpp:494
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:1152
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
void CalcVShape_ND(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe_base.cpp:969
void LocalRestriction_RT(const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &R) const
Definition: fe_base.cpp:1526
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition: fe_base.hpp:340
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
Definition: geom.cpp:405
Describes the function space on each element.
Definition: fe_base.hpp:216
void pts(int iphi, int t, double x[])
RefCoord s[3]
const DenseMatrix & Hessian()
Return the Hessian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:125
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:710
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
Implements CalcDShape methods.
Definition: fe_base.hpp:295
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
virtual void CalcPhysHessian(ElementTransformation &Trans, DenseMatrix &Hessian) const
Evaluate the Hessian of all shape functions of a scalar finite element in reference space at the give...
Definition: fe_base.cpp:288
virtual void GetLocalRestriction(ElementTransformation &Trans, DenseMatrix &R) const
Return a local restriction matrix R (Dof x Dof) mapping fine dofs to coarse dofs. ...
Definition: fe_base.cpp:113
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
void Project_ND(const double *tk, const Array< int > &d2t, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Project a vector coefficient onto the ND basis functions.
Definition: fe_base.cpp:1204
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
void ScalarLocalL2Restriction(ElementTransformation &Trans, DenseMatrix &R, const ScalarFiniteElement &coarse_fe) const
Get restriction matrix R defined through local L2-projection in the space defined by the coarse_fe...
Definition: fe_base.cpp:589
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
Definition: densemat.cpp:3075
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:2055
const IntegrationRule * IntRule
IntegrationRule that defines the quadrature points at which the basis functions of the FE are evaluat...
Definition: fe_base.hpp:146
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Evaluate the divergence of all shape functions of a vector finite element in reference space at the g...
Definition: fe_base.cpp:52
void ProjectGrad_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe_base.cpp:1332
int order
Order/degree of the shape functions.
Definition: fe_base.hpp:243
No derivatives implemented.
Definition: fe_base.hpp:294
Implements CalcCurlShape methods.
Definition: fe_base.hpp:297
void NodalLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I, const ScalarFiniteElement &fine_fe) const
Get the matrix I that defines nodal interpolation between this element and the refined element fine_f...
Definition: fe_base.cpp:509