MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lininteg.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 
13 #include "fem.hpp"
14 #include <cmath>
15 
16 namespace mfem
17 {
18 
20  const Array<int> &markers,
21  Vector &b)
22 {
23  MFEM_ABORT("Not supported.");
24 }
25 
27  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
28 {
29  mfem_error("LinearFormIntegrator::AssembleRHSElementVect(...)");
30 }
31 
33  const FiniteElement &el1, const FiniteElement &el2,
35 {
36  mfem_error("LinearFormIntegrator::AssembleRHSElementVect(...)");
37 }
38 
41  Vector &elvect)
42 {
43  int dof = el.GetDof();
44 
45  shape.SetSize(dof); // vector of size dof
46  elvect.SetSize(dof);
47  elvect = 0.0;
48 
49  const IntegrationRule *ir = IntRule;
50  if (ir == NULL)
51  {
52  // ir = &IntRules.Get(el.GetGeomType(),
53  // oa * el.GetOrder() + ob + Tr.OrderW());
54  ir = &IntRules.Get(el.GetGeomType(), oa * el.GetOrder() + ob);
55  }
56 
57  for (int i = 0; i < ir->GetNPoints(); i++)
58  {
59  const IntegrationPoint &ip = ir->IntPoint(i);
60 
61  Tr.SetIntPoint (&ip);
62  double val = Tr.Weight() * Q.Eval(Tr, ip);
63 
64  el.CalcPhysShape(Tr, shape);
65 
66  add(elvect, ip.weight * val, shape, elvect);
67  }
68 }
69 
71  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
72 {
73  MFEM_ASSERT(delta != NULL, "coefficient must be DeltaCoefficient");
74  elvect.SetSize(fe.GetDof());
75  fe.CalcPhysShape(Trans, elvect);
76  elvect *= delta->EvalDelta(Trans, Trans.GetIntPoint());
77 }
78 
80  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
81 {
82  int dof = el.GetDof();
83  int spaceDim = Tr.GetSpaceDim();
84 
85  dshape.SetSize(dof, spaceDim);
86 
87  elvect.SetSize(dof);
88  elvect = 0.0;
89 
90  const IntegrationRule *ir = IntRule;
91  if (ir == NULL)
92  {
93  int intorder = 2 * el.GetOrder();
94  ir = &IntRules.Get(el.GetGeomType(), intorder);
95  }
96 
97  for (int i = 0; i < ir->GetNPoints(); i++)
98  {
99  const IntegrationPoint &ip = ir->IntPoint(i);
100 
101  Tr.SetIntPoint(&ip);
102  el.CalcPhysDShape(Tr, dshape);
103 
104  Q.Eval(Qvec, Tr, ip);
105  Qvec *= ip.weight * Tr.Weight();
106 
107  dshape.AddMult(Qvec, elvect);
108  }
109 }
110 
112  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
113 {
114  MFEM_ASSERT(vec_delta != NULL,"coefficient must be VectorDeltaCoefficient");
115  int dof = fe.GetDof();
116  int spaceDim = Trans.GetSpaceDim();
117 
118  dshape.SetSize(dof, spaceDim);
119  fe.CalcPhysDShape(Trans, dshape);
120 
121  vec_delta->EvalDelta(Qvec, Trans, Trans.GetIntPoint());
122 
123  elvect.SetSize(dof);
124  dshape.Mult(Qvec, elvect);
125 }
126 
128  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
129 {
130  int dof = el.GetDof();
131 
132  shape.SetSize(dof); // vector of size dof
133  elvect.SetSize(dof);
134  elvect = 0.0;
135 
136  const IntegrationRule *ir = IntRule;
137  if (ir == NULL)
138  {
139  int intorder = oa * el.GetOrder() + ob; // <----------
140  ir = &IntRules.Get(el.GetGeomType(), intorder);
141  }
142 
143  for (int i = 0; i < ir->GetNPoints(); i++)
144  {
145  const IntegrationPoint &ip = ir->IntPoint(i);
146 
147  Tr.SetIntPoint (&ip);
148  double val = Tr.Weight() * Q.Eval(Tr, ip);
149 
150  el.CalcShape(ip, shape);
151 
152  add(elvect, ip.weight * val, shape, elvect);
153  }
154 }
155 
157  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
158 {
159  int dof = el.GetDof();
160 
161  shape.SetSize(dof); // vector of size dof
162  elvect.SetSize(dof);
163  elvect = 0.0;
164 
165  const IntegrationRule *ir = IntRule;
166  if (ir == NULL)
167  {
168  int intorder = oa * el.GetOrder() + ob; // <------ user control
169  ir = &IntRules.Get(Tr.FaceGeom, intorder); // of integration order
170  }
171 
172  for (int i = 0; i < ir->GetNPoints(); i++)
173  {
174  const IntegrationPoint &ip = ir->IntPoint(i);
175 
176  // Set the integration point in the face and the neighboring element
177  Tr.SetAllIntPoints(&ip);
178 
179  // Access the neighboring element's integration point
180  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
181 
182  double val = Tr.Face->Weight() * ip.weight * Q.Eval(*Tr.Face, ip);
183 
184  el.CalcShape(eip, shape);
185 
186  add(elvect, val, shape, elvect);
187  }
188 }
189 
191  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
192 {
193  int dim = el.GetDim()+1;
194  int dof = el.GetDof();
195  Vector nor(dim), Qvec;
196 
197  shape.SetSize(dof);
198  elvect.SetSize(dof);
199  elvect = 0.0;
200 
201  const IntegrationRule *ir = IntRule;
202  if (ir == NULL)
203  {
204  int intorder = oa * el.GetOrder() + ob; // <----------
205  ir = &IntRules.Get(el.GetGeomType(), intorder);
206  }
207 
208  for (int i = 0; i < ir->GetNPoints(); i++)
209  {
210  const IntegrationPoint &ip = ir->IntPoint(i);
211 
212  Tr.SetIntPoint(&ip);
213  if (dim > 1)
214  {
215  CalcOrtho(Tr.Jacobian(), nor);
216  }
217  else
218  {
219  nor[0] = 1.0;
220  }
221  Q.Eval(Qvec, Tr, ip);
222 
223  el.CalcShape(ip, shape);
224 
225  elvect.Add(ip.weight*(Qvec*nor), shape);
226  }
227 }
228 
230  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
231 {
232  int dim = el.GetDim()+1;
233  int dof = el.GetDof();
234  Vector tangent(dim), Qvec;
235 
236  shape.SetSize(dof);
237  elvect.SetSize(dof);
238  elvect = 0.0;
239 
240  if (dim != 2)
241  {
242  mfem_error("These methods make sense only in 2D problems.");
243  }
244 
245  const IntegrationRule *ir = IntRule;
246  if (ir == NULL)
247  {
248  int intorder = oa * el.GetOrder() + ob; // <----------
249  ir = &IntRules.Get(el.GetGeomType(), intorder);
250  }
251 
252  for (int i = 0; i < ir->GetNPoints(); i++)
253  {
254  const IntegrationPoint &ip = ir->IntPoint(i);
255 
256  Tr.SetIntPoint(&ip);
257  const DenseMatrix &Jac = Tr.Jacobian();
258  tangent(0) = Jac(0,0);
259  tangent(1) = Jac(1,0);
260 
261  Q.Eval(Qvec, Tr, ip);
262 
263  el.CalcShape(ip, shape);
264 
265  add(elvect, ip.weight*(Qvec*tangent), shape, elvect);
266  }
267 }
268 
270  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
271 {
272  int vdim = Q.GetVDim();
273  int dof = el.GetDof();
274 
275  double val,cf;
276 
277  shape.SetSize(dof); // vector of size dof
278 
279  elvect.SetSize(dof * vdim);
280  elvect = 0.0;
281 
282  const IntegrationRule *ir = IntRule;
283  if (ir == NULL)
284  {
285  int intorder = 2*el.GetOrder();
286  ir = &IntRules.Get(el.GetGeomType(), intorder);
287  }
288 
289  for (int i = 0; i < ir->GetNPoints(); i++)
290  {
291  const IntegrationPoint &ip = ir->IntPoint(i);
292 
293  Tr.SetIntPoint (&ip);
294  val = Tr.Weight();
295 
296  el.CalcShape(ip, shape);
297  Q.Eval (Qvec, Tr, ip);
298 
299  for (int k = 0; k < vdim; k++)
300  {
301  cf = val * Qvec(k);
302 
303  for (int s = 0; s < dof; s++)
304  {
305  elvect(dof*k+s) += ip.weight * cf * shape(s);
306  }
307  }
308  }
309 }
310 
312  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
313 {
314  MFEM_ASSERT(vec_delta != NULL, "coefficient must be VectorDeltaCoefficient");
315  int vdim = Q.GetVDim();
316  int dof = fe.GetDof();
317 
318  shape.SetSize(dof);
319  fe.CalcPhysShape(Trans, shape);
320 
321  vec_delta->EvalDelta(Qvec, Trans, Trans.GetIntPoint());
322 
323  elvect.SetSize(dof*vdim);
324  DenseMatrix elvec_as_mat(elvect.GetData(), dof, vdim);
325  MultVWt(shape, Qvec, elvec_as_mat);
326 }
327 
329  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
330 {
331  const int dim = el.GetDim();
332  const int dof = el.GetDof();
333  const int vdim = Q.GetVDim();
334  const int sdim = Tr.GetSpaceDim();
335 
336  dshape.SetSize(dof,sdim);
337 
338  elvect.SetSize(dof*(vdim/sdim));
339  elvect = 0.0;
340 
341  const IntegrationRule *ir = IntRule;
342  if (ir == NULL)
343  {
344  int intorder = 2 * el.GetOrder();
345  ir = &IntRules.Get(el.GetGeomType(), intorder);
346  }
347 
348  Vector pelvect(dof);
349  Vector part_x(dim);
350 
351  for (int q = 0; q < ir->GetNPoints(); q++)
352  {
353  const IntegrationPoint &ip = ir->IntPoint(q);
354 
355  Tr.SetIntPoint(&ip);
356  el.CalcPhysDShape(Tr, dshape);
357 
358  Q.Eval(Qvec, Tr, ip);
359  Qvec *= ip.weight * Tr.Weight();
360 
361  for (int k = 0; k < vdim/sdim; k++)
362  {
363  for (int d=0; d < sdim; ++d) { part_x(d) = Qvec(k*sdim+d); }
364  dshape.Mult(part_x, pelvect);
365  for (int s = 0; s < dof; ++s) { elvect(s+k*dof) += pelvect(s); }
366  }
367  }
368 }
369 
372 {
373  MFEM_ABORT("Not implemented!");
374 }
375 
377  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
378 {
379  int vdim = Q.GetVDim();
380  int dof = el.GetDof();
381 
382  shape.SetSize(dof);
383  vec.SetSize(vdim);
384 
385  elvect.SetSize(dof * vdim);
386  elvect = 0.0;
387 
388  const IntegrationRule *ir = IntRule;
389  if (ir == NULL)
390  {
391  int intorder = 2*el.GetOrder();
392  ir = &IntRules.Get(el.GetGeomType(), intorder);
393  }
394 
395  for (int i = 0; i < ir->GetNPoints(); i++)
396  {
397  const IntegrationPoint &ip = ir->IntPoint(i);
398 
399  Q.Eval(vec, Tr, ip);
400  Tr.SetIntPoint (&ip);
401  vec *= Tr.Weight() * ip.weight;
402  el.CalcShape(ip, shape);
403  for (int k = 0; k < vdim; k++)
404  for (int s = 0; s < dof; s++)
405  {
406  elvect(dof*k+s) += vec(k) * shape(s);
407  }
408  }
409 }
410 
412  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
413 {
414  int vdim = Q.GetVDim();
415  int dof = el.GetDof();
416 
417  shape.SetSize(dof);
418  vec.SetSize(vdim);
419 
420  elvect.SetSize(dof * vdim);
421  elvect = 0.0;
422 
423  const IntegrationRule *ir = IntRule;
424  if (ir == NULL)
425  {
426  int intorder = 2*el.GetOrder();
427  ir = &IntRules.Get(Tr.GetGeometryType(), intorder);
428  }
429 
430  for (int i = 0; i < ir->GetNPoints(); i++)
431  {
432  const IntegrationPoint &ip = ir->IntPoint(i);
433 
434  // Set the integration point in the face and the neighboring element
435  Tr.SetAllIntPoints(&ip);
436 
437  // Access the neighboring element's integration point
438  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
439 
440  // Use Tr transformation in case Q depends on boundary attribute
441  Q.Eval(vec, Tr, ip);
442  vec *= Tr.Weight() * ip.weight;
443  el.CalcShape(eip, shape);
444  for (int k = 0; k < vdim; k++)
445  {
446  for (int s = 0; s < dof; s++)
447  {
448  elvect(dof*k+s) += vec(k) * shape(s);
449  }
450  }
451  }
452 }
453 
455  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
456 {
457  int dof = el.GetDof();
458  int spaceDim = Tr.GetSpaceDim();
459  int vdim = std::max(spaceDim, el.GetVDim());
460 
461  vshape.SetSize(dof,vdim);
462  vec.SetSize(vdim);
463 
464  elvect.SetSize(dof);
465  elvect = 0.0;
466 
467  const IntegrationRule *ir = IntRule;
468  if (ir == NULL)
469  {
470  // int intorder = 2*el.GetOrder() - 1; // ok for O(h^{k+1}) conv. in L2
471  int intorder = 2*el.GetOrder();
472  ir = &IntRules.Get(el.GetGeomType(), intorder);
473  }
474 
475  for (int i = 0; i < ir->GetNPoints(); i++)
476  {
477  const IntegrationPoint &ip = ir->IntPoint(i);
478 
479  Tr.SetIntPoint (&ip);
480  el.CalcVShape(Tr, vshape);
481 
482  QF.Eval (vec, Tr, ip);
483  vec *= ip.weight * Tr.Weight();
484  vshape.AddMult (vec, elvect);
485  }
486 }
487 
489  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
490 {
491  MFEM_ASSERT(vec_delta != NULL, "coefficient must be VectorDeltaCoefficient");
492  int dof = fe.GetDof();
493  int spaceDim = Trans.GetSpaceDim();
494 
495  vshape.SetSize(dof, spaceDim);
496  fe.CalcPhysVShape(Trans, vshape);
497 
498  vec_delta->EvalDelta(vec, Trans, Trans.GetIntPoint());
499 
500  elvect.SetSize(dof);
501  vshape.Mult(vec, elvect);
502 }
503 
505  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
506 {
507  int dof = el.GetDof();
508  int spaceDim = Tr.GetSpaceDim();
509  int n=(spaceDim == 3)? spaceDim : 1;
510  curlshape.SetSize(dof,n);
511  vec.SetSize(n);
512 
513  elvect.SetSize(dof);
514  elvect = 0.0;
515 
516  const IntegrationRule *ir = IntRule;
517  if (ir == NULL)
518  {
519  int intorder = 2*el.GetOrder();
520  ir = &IntRules.Get(el.GetGeomType(), intorder);
521  }
522 
523  for (int i = 0; i < ir->GetNPoints(); i++)
524  {
525  const IntegrationPoint &ip = ir->IntPoint(i);
526 
527  Tr.SetIntPoint (&ip);
528  el.CalcPhysCurlShape(Tr, curlshape);
529  QF->Eval(vec, Tr, ip);
530 
531  vec *= ip.weight * Tr.Weight();
532  curlshape.AddMult (vec, elvect);
533  }
534 }
535 
537  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
538 {
539  int spaceDim = Trans.GetSpaceDim();
540  MFEM_ASSERT(vec_delta != NULL,
541  "coefficient must be VectorDeltaCoefficient");
542  int dof = fe.GetDof();
543  int n=(spaceDim == 3)? spaceDim : 1;
544  vec.SetSize(n);
545  curlshape.SetSize(dof, n);
546  elvect.SetSize(dof);
547  fe.CalcPhysCurlShape(Trans, curlshape);
548 
549  vec_delta->EvalDelta(vec, Trans, Trans.GetIntPoint());
550  curlshape.Mult(vec, elvect);
551 }
552 
554  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
555 {
556  int dof = el.GetDof();
557 
558  divshape.SetSize(dof); // vector of size dof
559  elvect.SetSize(dof);
560  elvect = 0.0;
561 
562  const IntegrationRule *ir = IntRule;
563  if (ir == NULL)
564  {
565  int intorder = 2 * el.GetOrder();
566  ir = &IntRules.Get(el.GetGeomType(), intorder);
567  }
568 
569  for (int i = 0; i < ir->GetNPoints(); i++)
570  {
571  const IntegrationPoint &ip = ir->IntPoint(i);
572 
573  Tr.SetIntPoint (&ip);
574  double val = Tr.Weight() * Q.Eval(Tr, ip);
575  el.CalcPhysDivShape(Tr, divshape);
576 
577  add(elvect, ip.weight * val, divshape, elvect);
578  }
579 }
580 
582  const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
583 {
584  MFEM_ASSERT(delta != NULL, "coefficient must be DeltaCoefficient");
585  elvect.SetSize(fe.GetDof());
586  fe.CalcPhysDivShape(Trans, elvect);
587  elvect *= delta->EvalDelta(Trans, Trans.GetIntPoint());
588 }
589 
591  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
592 {
593  int dim = el.GetDim()+1;
594  int dof = el.GetDof();
595 
596  shape.SetSize (dof);
597  nor.SetSize (dim);
598  elvect.SetSize (dim*dof);
599 
600  const IntegrationRule *ir = IntRule;
601  if (ir == NULL)
602  {
603  ir = &IntRules.Get(el.GetGeomType(), el.GetOrder() + 1);
604  }
605 
606  elvect = 0.0;
607  for (int i = 0; i < ir->GetNPoints(); i++)
608  {
609  const IntegrationPoint &ip = ir->IntPoint(i);
610  Tr.SetIntPoint (&ip);
611  CalcOrtho(Tr.Jacobian(), nor);
612  el.CalcShape (ip, shape);
613  nor *= Sign * ip.weight * F -> Eval (Tr, ip);
614  for (int j = 0; j < dof; j++)
615  for (int k = 0; k < dim; k++)
616  {
617  elvect(dof*k+j) += nor(k) * shape(j);
618  }
619  }
620 }
621 
622 
624  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
625 {
626  int dof = el.GetDof();
627 
628  shape.SetSize(dof);
629  elvect.SetSize(dof);
630  elvect = 0.0;
631 
632  const IntegrationRule *ir = IntRule;
633  if (ir == NULL)
634  {
635  int intorder = oa * el.GetOrder() + ob; // <----------
636  ir = &IntRules.Get(el.GetGeomType(), intorder);
637  }
638 
639  for (int i = 0; i < ir->GetNPoints(); i++)
640  {
641  const IntegrationPoint &ip = ir->IntPoint(i);
642  el.CalcShape(ip, shape);
643 
644  double val = ip.weight;
645  if (F)
646  {
647  Tr.SetIntPoint (&ip);
648  val *= F->Eval(Tr, ip);
649  }
650 
651  elvect.Add(val, shape);
652  }
653 }
654 
656  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
657 {
658  int dof = el.GetDof();
659  int dim = el.GetDim();
660  int vdim = el.GetVDim();
661  DenseMatrix vshape(dof, vdim);
662  Vector f_loc(3);
663  Vector f_hat(2);
664 
665  MFEM_VERIFY(vdim == 2, "VectorFEBoundaryTangentLFIntegrator "
666  "must be called with vector basis functions of dimension 2.");
667 
668  elvect.SetSize(dof);
669  elvect = 0.0;
670 
671  const IntegrationRule *ir = IntRule;
672  if (ir == NULL)
673  {
674  int intorder = oa * el.GetOrder() + ob; // <----------
675  ir = &IntRules.Get(el.GetGeomType(), intorder);
676  }
677 
678  for (int i = 0; i < ir->GetNPoints(); i++)
679  {
680  const IntegrationPoint &ip = ir->IntPoint(i);
681 
682  el.CalcVShape(ip, vshape);
683 
684  Tr.SetIntPoint(&ip);
685  f.Eval(f_loc, Tr, ip);
686 
687  if (dim == 2)
688  {
689  Tr.Jacobian().MultTranspose(f_loc, f_hat);
690  }
691  else if (dim == 1)
692  {
693  const DenseMatrix & J = Tr.Jacobian();
694  f_hat(0) = J(0,0) * f_loc(0) + J(1,0) * f_loc(1);
695  f_hat(1) = f_loc(2);
696  }
697  else
698  {
699  f_hat(0) = f_loc(1);
700  f_hat(1) = f_loc(2);
701  }
702 
703  Swap<double>(f_hat(0), f_hat(1));
704  f_hat(0) = -f_hat(0);
705  f_hat *= ip.weight;
706 
707  vshape.AddMult(f_hat, elvect);
708  }
709 }
710 
712  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
713 {
714  mfem_error("BoundaryFlowIntegrator::AssembleRHSElementVect\n"
715  " is not implemented as boundary integrator!\n"
716  " Use LinearForm::AddBdrFaceIntegrator instead of\n"
717  " LinearForm::AddBoundaryIntegrator.");
718 }
719 
721  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
722 {
723  int dim, ndof, order;
724  double un, w, vu_data[3], nor_data[3];
725 
726  dim = el.GetDim();
727  ndof = el.GetDof();
728  Vector vu(vu_data, dim), nor(nor_data, dim);
729 
730  const IntegrationRule *ir = IntRule;
731  if (ir == NULL)
732  {
733  // Assuming order(u)==order(mesh)
734  order = Tr.Elem1->OrderW() + 2*el.GetOrder();
735  if (el.Space() == FunctionSpace::Pk)
736  {
737  order++;
738  }
739  ir = &IntRules.Get(Tr.GetGeometryType(), order);
740  }
741 
742  shape.SetSize(ndof);
743  elvect.SetSize(ndof);
744  elvect = 0.0;
745 
746  for (int p = 0; p < ir->GetNPoints(); p++)
747  {
748  const IntegrationPoint &ip = ir->IntPoint(p);
749 
750  // Set the integration point in the face and the neighboring element
751  Tr.SetAllIntPoints(&ip);
752 
753  // Access the neighboring element's integration point
754  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
755  el.CalcShape(eip, shape);
756 
757  // Use Tr.Elem1 transformation for u so that it matches the coefficient
758  // used with the ConvectionIntegrator and/or the DGTraceIntegrator.
759  u->Eval(vu, *Tr.Elem1, eip);
760 
761  if (dim == 1)
762  {
763  nor(0) = 2*eip.x - 1.0;
764  }
765  else
766  {
767  CalcOrtho(Tr.Jacobian(), nor);
768  }
769 
770  un = vu * nor;
771  w = 0.5*alpha*un - beta*fabs(un);
772  w *= ip.weight*f->Eval(Tr, ip);
773  elvect.Add(w, shape);
774  }
775 }
776 
778  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
779 {
780  mfem_error("DGDirichletLFIntegrator::AssembleRHSElementVect");
781 }
782 
784  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
785 {
786  int dim, ndof;
787  bool kappa_is_nonzero = (kappa != 0.);
788  double w;
789 
790  dim = el.GetDim();
791  ndof = el.GetDof();
792 
793  nor.SetSize(dim);
794  nh.SetSize(dim);
795  ni.SetSize(dim);
796  adjJ.SetSize(dim);
797  if (MQ)
798  {
799  mq.SetSize(dim);
800  }
801 
802  shape.SetSize(ndof);
803  dshape.SetSize(ndof, dim);
804  dshape_dn.SetSize(ndof);
805 
806  elvect.SetSize(ndof);
807  elvect = 0.0;
808 
809  const IntegrationRule *ir = IntRule;
810  if (ir == NULL)
811  {
812  // a simple choice for the integration order; is this OK?
813  int order = 2*el.GetOrder();
814  ir = &IntRules.Get(Tr.GetGeometryType(), order);
815  }
816 
817  for (int p = 0; p < ir->GetNPoints(); p++)
818  {
819  const IntegrationPoint &ip = ir->IntPoint(p);
820 
821  // Set the integration point in the face and the neighboring element
822  Tr.SetAllIntPoints(&ip);
823 
824  // Access the neighboring element's integration point
825  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
826 
827  if (dim == 1)
828  {
829  nor(0) = 2*eip.x - 1.0;
830  }
831  else
832  {
833  CalcOrtho(Tr.Jacobian(), nor);
834  }
835 
836  el.CalcShape(eip, shape);
837  el.CalcDShape(eip, dshape);
838 
839  // compute uD through the face transformation
840  w = ip.weight * uD->Eval(Tr, ip) / Tr.Elem1->Weight();
841  if (!MQ)
842  {
843  if (Q)
844  {
845  w *= Q->Eval(*Tr.Elem1, eip);
846  }
847  ni.Set(w, nor);
848  }
849  else
850  {
851  nh.Set(w, nor);
852  MQ->Eval(mq, *Tr.Elem1, eip);
853  mq.MultTranspose(nh, ni);
854  }
856  adjJ.Mult(ni, nh);
857 
859  elvect.Add(sigma, dshape_dn);
860 
861  if (kappa_is_nonzero)
862  {
863  elvect.Add(kappa*(ni*nor), shape);
864  }
865  }
866 }
867 
869  const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
870 {
871  mfem_error("DGElasticityDirichletLFIntegrator::AssembleRHSElementVect");
872 }
873 
875  const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect)
876 {
877  MFEM_ASSERT(Tr.Elem2No < 0, "interior boundary is not supported");
878 
879 #ifdef MFEM_THREAD_SAFE
880  Vector shape;
884  Vector nor;
887  Vector u_dir;
888 #endif
889 
890  const int dim = el.GetDim();
891  const int ndofs = el.GetDof();
892  const int nvdofs = dim*ndofs;
893 
894  elvect.SetSize(nvdofs);
895  elvect = 0.0;
896 
897  adjJ.SetSize(dim);
898  shape.SetSize(ndofs);
899  dshape.SetSize(ndofs, dim);
900  dshape_ps.SetSize(ndofs, dim);
901  nor.SetSize(dim);
902  dshape_dn.SetSize(ndofs);
903  dshape_du.SetSize(ndofs);
904  u_dir.SetSize(dim);
905 
906  const IntegrationRule *ir = IntRule;
907  if (ir == NULL)
908  {
909  const int order = 2*el.GetOrder(); // <-----
910  ir = &IntRules.Get(Tr.GetGeometryType(), order);
911  }
912 
913  for (int pi = 0; pi < ir->GetNPoints(); ++pi)
914  {
915  const IntegrationPoint &ip = ir->IntPoint(pi);
916 
917  // Set the integration point in the face and the neighboring element
918  Tr.SetAllIntPoints(&ip);
919 
920  // Access the neighboring element's integration point
921  const IntegrationPoint &eip = Tr.GetElement1IntPoint();
922 
923  // Evaluate the Dirichlet b.c. using the face transformation.
924  uD.Eval(u_dir, Tr, ip);
925 
926  el.CalcShape(eip, shape);
927  el.CalcDShape(eip, dshape);
928 
930  Mult(dshape, adjJ, dshape_ps);
931 
932  if (dim == 1)
933  {
934  nor(0) = 2*eip.x - 1.0;
935  }
936  else
937  {
938  CalcOrtho(Tr.Jacobian(), nor);
939  }
940 
941  double wL, wM, jcoef;
942  {
943  const double w = ip.weight / Tr.Elem1->Weight();
944  wL = w * lambda->Eval(*Tr.Elem1, eip);
945  wM = w * mu->Eval(*Tr.Elem1, eip);
946  jcoef = kappa * (wL + 2.0*wM) * (nor*nor);
947  dshape_ps.Mult(nor, dshape_dn);
948  dshape_ps.Mult(u_dir, dshape_du);
949  }
950 
951  // alpha < uD, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
952  // + kappa < h^{-1} (lambda + 2 mu) uD, v >
953 
954  // i = idof + ndofs * im
955  // v_phi(i,d) = delta(im,d) phi(idof)
956  // div(v_phi(i)) = dphi(idof,im)
957  // (grad(v_phi(i)))(k,l) = delta(im,k) dphi(idof,l)
958  //
959  // term 1:
960  // alpha < uD, lambda div(v_phi(i)) n >
961  // alpha lambda div(v_phi(i)) (uD.n) =
962  // alpha lambda dphi(idof,im) (uD.n) --> quadrature -->
963  // ip.weight/det(J1) alpha lambda (uD.nor) dshape_ps(idof,im) =
964  // alpha * wL * (u_dir*nor) * dshape_ps(idof,im)
965  // term 2:
966  // < alpha uD, mu grad(v_phi(i)).n > =
967  // alpha mu uD^T grad(v_phi(i)) n =
968  // alpha mu uD(k) delta(im,k) dphi(idof,l) n(l) =
969  // alpha mu uD(im) dphi(idof,l) n(l) --> quadrature -->
970  // ip.weight/det(J1) alpha mu uD(im) dshape_ps(idof,l) nor(l) =
971  // alpha * wM * u_dir(im) * dshape_dn(idof)
972  // term 3:
973  // < alpha uD, mu (grad(v_phi(i)))^T n > =
974  // alpha mu n^T grad(v_phi(i)) uD =
975  // alpha mu n(k) delta(im,k) dphi(idof,l) uD(l) =
976  // alpha mu n(im) dphi(idof,l) uD(l) --> quadrature -->
977  // ip.weight/det(J1) alpha mu nor(im) dshape_ps(idof,l) uD(l) =
978  // alpha * wM * nor(im) * dshape_du(idof)
979  // term j:
980  // < kappa h^{-1} (lambda + 2 mu) uD, v_phi(i) > =
981  // kappa/h (lambda + 2 mu) uD(k) v_phi(i,k) =
982  // kappa/h (lambda + 2 mu) uD(k) delta(im,k) phi(idof) =
983  // kappa/h (lambda + 2 mu) uD(im) phi(idof) --> quadrature -->
984  // [ 1/h = |nor|/det(J1) ]
985  // ip.weight/det(J1) |nor|^2 kappa (lambda + 2 mu) uD(im) phi(idof) =
986  // jcoef * u_dir(im) * shape(idof)
987 
988  wM *= alpha;
989  const double t1 = alpha * wL * (u_dir*nor);
990  for (int im = 0, i = 0; im < dim; ++im)
991  {
992  const double t2 = wM * u_dir(im);
993  const double t3 = wM * nor(im);
994  const double tj = jcoef * u_dir(im);
995  for (int idof = 0; idof < ndofs; ++idof, ++i)
996  {
997  elvect(i) += (t1*dshape_ps(idof,im) + t2*dshape_dn(idof) +
998  t3*dshape_du(idof) + tj*shape(idof));
999  }
1000  }
1001  }
1002 }
1003 
1004 
1005 
1007 (const FiniteElement &el,
1009  Vector &elvect)
1010 {
1011  int n = el.GetDof();
1012  elvect.SetSize(n);
1013  for (int i = 0; i < n; i++)
1014  {
1015  elvect(i) = dist(generator);
1016  }
1017 
1018  int iel = Tr.ElementNo;
1019 
1020  if (!save_factors || !L[iel])
1021  {
1022  DenseMatrix *M, m;
1023  if (save_factors)
1024  {
1025  L[iel]=new DenseMatrix;
1026  M = L[iel];
1027  }
1028  else
1029  {
1030  M = &m;
1031  }
1032  massinteg.AssembleElementMatrix(el, Tr, *M);
1033  CholeskyFactors chol(M->Data());
1034  chol.Factor(M->Height());
1035  chol.LMult(n,1,elvect);
1036  }
1037  else
1038  {
1039  CholeskyFactors chol(L[iel]->Data());
1040  chol.LMult(n,1,elvect);
1041  }
1042 }
1043 
1044 
1046  const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
1047 {
1048  const IntegrationRule *ir =
1050 
1051  const int nqp = ir->GetNPoints();
1052  const int vdim = vqfc.GetVDim();
1053  const int ndofs = fe.GetDof();
1054  Vector shape(ndofs);
1055  Vector temp(vdim);
1056  elvect.SetSize(vdim * ndofs);
1057  elvect = 0.0;
1058  for (int q = 0; q < nqp; q++)
1059  {
1060  const IntegrationPoint &ip = ir->IntPoint(q);
1061  Tr.SetIntPoint(&ip);
1062  const double w = Tr.Weight() * ip.weight;
1063  vqfc.Eval(temp, Tr, ip);
1064  fe.CalcShape(ip, shape);
1065  for (int ind = 0; ind < vdim; ind++)
1066  {
1067  for (int nd = 0; nd < ndofs; nd++)
1068  {
1069  elvect(nd + ind * ndofs) += w * shape(nd) * temp(ind);
1070  }
1071  }
1072  }
1073 }
1074 
1075 
1078  Vector &elvect)
1079 {
1080  const IntegrationRule *ir =
1082 
1083  const int nqp = ir->GetNPoints();
1084  const int ndofs = fe.GetDof();
1085  Vector shape(ndofs);
1086  elvect.SetSize(ndofs);
1087  elvect = 0.0;
1088  for (int q = 0; q < nqp; q++)
1089  {
1090  const IntegrationPoint &ip = ir->IntPoint(q);
1091  Tr.SetIntPoint (&ip);
1092  const double w = Tr.Weight() * ip.weight;
1093  double temp = qfc.Eval(Tr, ip);
1094  fe.CalcShape(ip, shape);
1095  shape *= (w * temp);
1096  elvect += shape;
1097  }
1098 }
1099 
1100 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:111
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:623
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:581
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:1007
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:39
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
Definition: densemat.cpp:3108
virtual int OrderW() const =0
Return the order of the determinant of the Jacobian (weight) of the transformation.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:162
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:553
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:57
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 SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:711
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:868
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:79
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 ...
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:70
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2440
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
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:269
int Space() const
Returns the type of FunctionSpace on the element.
Definition: fe_base.hpp:337
void LMult(int m, int n, double *X) const
Definition: densemat.cpp:3627
const IntegrationRule * IntRule
Definition: lininteg.hpp:27
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:98
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:209
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2654
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:127
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:313
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:70
Polynomials of order k.
Definition: fe_base.hpp:223
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:201
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:370
const QuadratureFunction & GetQuadFunction() const
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
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
double b
Definition: lissajous.cpp:42
void SetAllIntPoints(const IntegrationPoint *face_ip)
Set the integration point in the Face and the two neighboring elements, if present.
Definition: eltrans.hpp:566
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 AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:376
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...
VectorDeltaCoefficient * vec_delta
Definition: lininteg.hpp:66
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:536
MatrixCoefficient * MQ
Definition: lininteg.hpp:528
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
Definition: lininteg.cpp:19
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:191
int GetVDim()
Returns dimension of the vector.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
void AddMult(const Vector &x, Vector &y) const
y += A.x
Definition: densemat.cpp:224
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:229
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:777
void CalcPhysVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Equivalent to the CalcVShape() method with the same arguments.
Definition: fe_base.hpp:412
int GetVDim() const
Returns the vector dimension for vector-valued finite elements.
Definition: fe_base.hpp:314
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)=0
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
Definition: qspace.hpp:72
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
virtual bool Factor(int m, double TOL=0.0)
Compute the Cholesky factorization of the current matrix.
Definition: densemat.cpp:3575
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:181
DeltaCoefficient * delta
Definition: lininteg.hpp:65
const IntegrationPoint & GetElement1IntPoint()
Get a const reference to the integration point in neighboring element 1 corresponding to the currentl...
Definition: eltrans.hpp:577
virtual void EvalDelta(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)
Return the specified direction vector multiplied by the value returned by DeltaCoefficient::EvalDelta...
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:590
virtual double EvalDelta(ElementTransformation &T, const IntegrationPoint &ip)
The value of the function assuming we are evaluating at the delta center.
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:267
const QuadratureFunction & GetQuadFunction() const
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:655
virtual int GetSpaceDim() const =0
Get the dimension of the target (physical) space.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:252
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
Definition: qfunction.hpp:70
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:488
int dim
Definition: ex24.cpp:53
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:1045
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:190
ElementTransformation * Elem1
Definition: eltrans.hpp:522
ElementTransformation * Face
Definition: eltrans.hpp:523
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:504
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:1076
Vector data type.
Definition: vector.hpp:60
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:173
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:454
RefCoord s[3]
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:311
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:39
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition: lininteg.cpp:328
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...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379