MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilininteg.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_BILININTEG
13 #define MFEM_BILININTEG
14 
15 #include "../config/config.hpp"
16 #include "nonlininteg.hpp"
17 
18 namespace mfem
19 {
20 
21 /// Abstract base class BilinearFormIntegrator
23 {
24 protected:
26 
28  { IntRule = ir; }
29 
30 public:
31  /// Given a particular Finite Element computes the element matrix elmat.
32  virtual void AssembleElementMatrix(const FiniteElement &el,
34  DenseMatrix &elmat);
35 
36  /** Compute the local matrix representation of a bilinear form
37  a(u,v) defined on different trial (given by u) and test
38  (given by v) spaces. The rows in the local matrix correspond
39  to the test dofs and the columns -- to the trial dofs. */
40  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
41  const FiniteElement &test_fe,
43  DenseMatrix &elmat);
44 
45  virtual void AssembleFaceMatrix(const FiniteElement &el1,
46  const FiniteElement &el2,
48  DenseMatrix &elmat);
49 
50  /** Abstract method used for assembling TraceFaceIntegrators in a
51  MixedBilinearForm. */
52  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
53  const FiniteElement &test_fe1,
54  const FiniteElement &test_fe2,
56  DenseMatrix &elmat);
57 
58  /// Perform the local action of the BilinearFormIntegrator
59  virtual void AssembleElementVector(const FiniteElement &el,
61  const Vector &elfun, Vector &elvect);
62 
63  virtual void AssembleElementGrad(const FiniteElement &el,
65  const Vector &elfun, DenseMatrix &elmat)
66  { AssembleElementMatrix(el, Tr, elmat); }
67 
68  virtual void ComputeElementFlux(const FiniteElement &el,
70  Vector &u,
71  const FiniteElement &fluxelem,
72  Vector &flux, int with_coef = 1) { }
73 
74  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
76  Vector &flux, Vector *d_energy = NULL)
77  { return 0.0; }
78 
79  void SetIntRule(const IntegrationRule *ir) { IntRule = ir; }
80 
82 };
83 
85 {
86 private:
87  int own_bfi;
89 
90  DenseMatrix bfi_elmat;
91 
92 public:
93  TransposeIntegrator (BilinearFormIntegrator *_bfi, int _own_bfi = 1)
94  { bfi = _bfi; own_bfi = _own_bfi; }
95 
96  virtual void AssembleElementMatrix(const FiniteElement &el,
98  DenseMatrix &elmat);
99 
100  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
101  const FiniteElement &test_fe,
103  DenseMatrix &elmat);
104 
106  virtual void AssembleFaceMatrix(const FiniteElement &el1,
107  const FiniteElement &el2,
109  DenseMatrix &elmat);
110 
111  virtual ~TransposeIntegrator() { if (own_bfi) { delete bfi; } }
112 };
113 
115 {
116 private:
117  int own_bfi;
119 
120 public:
121  LumpedIntegrator (BilinearFormIntegrator *_bfi, int _own_bfi = 1)
122  { bfi = _bfi; own_bfi = _own_bfi; }
123 
124  virtual void AssembleElementMatrix(const FiniteElement &el,
126  DenseMatrix &elmat);
127 
128  virtual ~LumpedIntegrator() { if (own_bfi) { delete bfi; } }
129 };
130 
131 /// Integrator that inverts the matrix assembled by another integrator.
133 {
134 private:
135  int own_integrator;
136  BilinearFormIntegrator *integrator;
137 
138 public:
139  InverseIntegrator(BilinearFormIntegrator *integ, int own_integ = 1)
140  { integrator = integ; own_integrator = own_integ; }
141 
142  virtual void AssembleElementMatrix(const FiniteElement &el,
144  DenseMatrix &elmat);
145 
146  virtual ~InverseIntegrator() { if (own_integrator) { delete integrator; } }
147 };
148 
149 /// Integrator defining a sum of multiple Integrators.
151 {
152 private:
153  int own_integrators;
154  DenseMatrix elem_mat;
155  Array<BilinearFormIntegrator*> integrators;
156 
157 public:
158  SumIntegrator(int own_integs = 1) { own_integrators = own_integs; }
159 
161  { integrators.Append(integ); }
162 
163  virtual void AssembleElementMatrix(const FiniteElement &el,
165  DenseMatrix &elmat);
166 
167  virtual ~SumIntegrator();
168 };
169 
170 /** An abstract class for integrating the product of two scalar basis functions
171  with an optional scalar coefficient. */
173 {
174 public:
175 
176  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
177  const FiniteElement &test_fe,
179  DenseMatrix &elmat);
180 
181  /// Support for use in BilinearForm. Can be used only when appropriate.
182  virtual void AssembleElementMatrix(const FiniteElement &fe,
183  ElementTransformation &Trans,
184  DenseMatrix &elmat)
185  { AssembleElementMatrix2(fe, fe, Trans, elmat); }
186 
187 protected:
188  /// This parameter can be set by derived methods to enable single shape
189  /// evaluation in case CalcTestShape() and CalcTrialShape() return the same
190  /// result if given the same FiniteElement. The default is false.
192 
193  MixedScalarIntegrator() : same_calc_shape(false), Q(NULL) {}
195 
196  inline virtual bool VerifyFiniteElementTypes(
197  const FiniteElement & trial_fe,
198  const FiniteElement & test_fe) const
199  {
200  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
202  }
203 
204  inline virtual const char * FiniteElementTypeFailureMessage() const
205  {
206  return "MixedScalarIntegrator: "
207  "Trial and test spaces must both be scalar fields.";
208  }
209 
210  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
211  const FiniteElement & test_fe,
212  ElementTransformation &Trans)
213  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
214 
215 
216  inline virtual void CalcTestShape(const FiniteElement & test_fe,
217  ElementTransformation &Trans,
218  Vector & shape)
219  { test_fe.CalcPhysShape(Trans, shape); }
220 
221  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
222  ElementTransformation &Trans,
223  Vector & shape)
224  { trial_fe.CalcPhysShape(Trans, shape); }
225 
226 private:
227 
228  Coefficient *Q;
229 
230 #ifndef MFEM_THREAD_SAFE
231  Vector test_shape;
232  Vector trial_shape;
233 #endif
234 
235 };
236 
237 /** An abstract class for integrating the inner product of two vector basis
238  functions with an optional scalar, vector, or matrix coefficient. */
240 {
241 public:
242 
243  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
244  const FiniteElement &test_fe,
246  DenseMatrix &elmat);
247 
248  /// Support for use in BilinearForm. Can be used only when appropriate.
249  virtual void AssembleElementMatrix(const FiniteElement &fe,
250  ElementTransformation &Trans,
251  DenseMatrix &elmat)
252  { AssembleElementMatrix2(fe, fe, Trans, elmat); }
253 
254 protected:
255  /// This parameter can be set by derived methods to enable single shape
256  /// evaluation in case CalcTestShape() and CalcTrialShape() return the same
257  /// result if given the same FiniteElement. The default is false.
259 
261  : same_calc_shape(false), Q(NULL), VQ(NULL), DQ(NULL), MQ(NULL) {}
263  : same_calc_shape(false), Q(&q), VQ(NULL), DQ(NULL), MQ(NULL) {}
265  : same_calc_shape(false), Q(NULL), VQ(diag?NULL:&dq), DQ(diag?&dq:NULL),
266  MQ(NULL) {}
268  : same_calc_shape(false), Q(NULL), VQ(NULL), DQ(NULL), MQ(&mq) {}
269 
270  inline virtual bool VerifyFiniteElementTypes(
271  const FiniteElement & trial_fe,
272  const FiniteElement & test_fe) const
273  {
274  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
276  }
277 
278  inline virtual const char * FiniteElementTypeFailureMessage() const
279  {
280  return "MixedVectorIntegrator: "
281  "Trial and test spaces must both be vector fields";
282  }
283 
284  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
285  const FiniteElement & test_fe,
286  ElementTransformation &Trans)
287  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
288 
289 
290  inline virtual void CalcTestShape(const FiniteElement & test_fe,
291  ElementTransformation &Trans,
292  DenseMatrix & shape)
293  { test_fe.CalcVShape(Trans, shape); }
294 
295  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
296  ElementTransformation &Trans,
297  DenseMatrix & shape)
298  { trial_fe.CalcVShape(Trans, shape); }
299 
300 private:
301 
302  Coefficient *Q;
303  VectorCoefficient *VQ;
304  VectorCoefficient *DQ;
305  MatrixCoefficient *MQ;
306 
307 #ifndef MFEM_THREAD_SAFE
308  Vector V;
309  Vector D;
310  DenseMatrix M;
311  DenseMatrix test_shape;
312  DenseMatrix trial_shape;
313  DenseMatrix test_shape_tmp;
314 #endif
315 
316 };
317 
318 /** An abstract class for integrating the product of a scalar basis function and
319  the inner product of a vector basis function with a vector coefficient. In
320  2D the inner product can be replaced with a cross product. */
322 {
323 public:
324 
325  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
326  const FiniteElement &test_fe,
328  DenseMatrix &elmat);
329 
330 protected:
331 
332  MixedScalarVectorIntegrator(VectorCoefficient &vq, bool _transpose = false,
333  bool _cross_2d = false)
334  : VQ(&vq), transpose(_transpose) , cross_2d(_cross_2d) {}
335 
336  inline virtual bool VerifyFiniteElementTypes(
337  const FiniteElement & trial_fe,
338  const FiniteElement & test_fe) const
339  {
340  return ((transpose &&
342  test_fe.GetRangeType() == mfem::FiniteElement::SCALAR ) ||
343  (!transpose &&
346  );
347  }
348 
349  inline virtual const char * FiniteElementTypeFailureMessage() const
350  {
351  if ( transpose )
352  {
353  return "MixedScalarVectorIntegrator: "
354  "Trial space must be a vector field "
355  "and the test space must be a scalar field";
356  }
357  else
358  {
359  return "MixedScalarVectorIntegrator: "
360  "Trial space must be a scalar field "
361  "and the test space must be a vector field";
362  }
363  }
364 
365  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
366  const FiniteElement & test_fe,
367  ElementTransformation &Trans)
368  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
369 
370 
371  inline virtual void CalcVShape(const FiniteElement & vector_fe,
372  ElementTransformation &Trans,
373  DenseMatrix & shape)
374  { vector_fe.CalcVShape(Trans, shape); }
375 
376  inline virtual void CalcShape(const FiniteElement & scalar_fe,
377  ElementTransformation &Trans,
378  Vector & shape)
379  { scalar_fe.CalcPhysShape(Trans, shape); }
380 
381 private:
382 
383  VectorCoefficient *VQ;
384  bool transpose;
385  bool cross_2d; // In 2D use a cross product rather than a dot product
386 
387 #ifndef MFEM_THREAD_SAFE
388  Vector V;
389  DenseMatrix vshape;
390  Vector shape;
391  Vector vshape_tmp;
392 #endif
393 
394 };
395 
396 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 1D, 2D,
397  or 3D and where Q is an optional scalar coefficient, u and v are each in H1
398  or L2. */
400 {
401 public:
404  : MixedScalarIntegrator(q) { same_calc_shape = true; }
405 };
406 
407 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 2D, or
408  3D and where Q is a vector coefficient, u is in H1 or L2 and v is in H(Curl)
409  or H(Div). */
411 {
412 public:
415 };
416 
417 /** Class for integrating the bilinear form a(u,v) := (Q D u, v) in 1D where Q
418  is an optional scalar coefficient, u is in H1, and v is in L2. */
420 {
421 public:
424  : MixedScalarIntegrator(q) {}
425 
426 protected:
427  inline virtual bool VerifyFiniteElementTypes(
428  const FiniteElement & trial_fe,
429  const FiniteElement & test_fe) const
430  {
431  return (trial_fe.GetDim() == 1 && test_fe.GetDim() == 1 &&
432  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
434  }
435 
436  inline virtual const char * FiniteElementTypeFailureMessage() const
437  {
438  return "MixedScalarDerivativeIntegrator: "
439  "Trial and test spaces must both be scalar fields in 1D "
440  "and the trial space must implement CaldDShape.";
441  }
442 
443  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
445  Vector & shape)
446  {
447  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
448  trial_fe.CalcPhysDShape(Trans, dshape);
449  }
450 };
451 
452 /** Class for integrating the bilinear form a(u,v) := -(Q u, D v) in 1D where Q
453  is an optional scalar coefficient, u is in L2, and v is in H1. */
455 {
456 public:
459  : MixedScalarIntegrator(q) {}
460 
461 protected:
462  inline virtual bool VerifyFiniteElementTypes(
463  const FiniteElement & trial_fe,
464  const FiniteElement & test_fe) const
465  {
466  return (trial_fe.GetDim() == 1 && test_fe.GetDim() == 1 &&
469  }
470 
471  inline virtual const char * FiniteElementTypeFailureMessage() const
472  {
473  return "MixedScalarWeakDerivativeIntegrator: "
474  "Trial and test spaces must both be scalar fields in 1D "
475  "and the test space must implement CalcDShape with "
476  "map type \"VALUE\".";
477  }
478 
479  inline virtual void CalcTestShape(const FiniteElement & test_fe,
481  Vector & shape)
482  {
483  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
484  test_fe.CalcPhysDShape(Trans, dshape);
485  shape *= -1.0;
486  }
487 };
488 
489 /** Class for integrating the bilinear form a(u,v) := (Q div u, v) in either 2D
490  or 3D where Q is an optional scalar coefficient, u is in H(Div), and v is a
491  scalar field. */
493 {
494 public:
497  : MixedScalarIntegrator(q) {}
498 
499 protected:
500  inline virtual bool VerifyFiniteElementTypes(
501  const FiniteElement & trial_fe,
502  const FiniteElement & test_fe) const
503  {
504  return (trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
506  }
507 
508  inline virtual const char * FiniteElementTypeFailureMessage() const
509  {
510  return "MixedScalarDivergenceIntegrator: "
511  "Trial must be H(Div) and the test space must be a "
512  "scalar field";
513  }
514 
515  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
516  const FiniteElement & test_fe,
518  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
519 
520  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
522  Vector & shape)
523  { trial_fe.CalcPhysDivShape(Trans, shape); }
524 };
525 
526 /** Class for integrating the bilinear form a(u,v) := (V div u, v) in either 2D
527  or 3D where V is a vector coefficient, u is in H(Div), and v is a vector
528  field. */
530 {
531 public:
534 
535 protected:
536  inline virtual bool VerifyFiniteElementTypes(
537  const FiniteElement & trial_fe,
538  const FiniteElement & test_fe) const
539  {
540  return (trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
542  }
543 
544  inline virtual const char * FiniteElementTypeFailureMessage() const
545  {
546  return "MixedVectorDivergenceIntegrator: "
547  "Trial must be H(Div) and the test space must be a "
548  "vector field";
549  }
550 
551  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
552  const FiniteElement & test_fe,
554  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
555 
556  inline virtual void CalcShape(const FiniteElement & scalar_fe,
558  Vector & shape)
559  { scalar_fe.CalcPhysDivShape(Trans, shape); }
560 };
561 
562 /** Class for integrating the bilinear form a(u,v) := -(Q u, div v) in either 2D
563  or 3D where Q is an optional scalar coefficient, u is in L2 or H1, and v is
564  in H(Div). */
566 {
567 public:
570  : MixedScalarIntegrator(q) {}
571 
572 protected:
573  inline virtual bool VerifyFiniteElementTypes(
574  const FiniteElement & trial_fe,
575  const FiniteElement & test_fe) const
576  {
577  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
578  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
579  }
580 
581  inline virtual const char * FiniteElementTypeFailureMessage() const
582  {
583  return "MixedScalarWeakGradientIntegrator: "
584  "Trial space must be a scalar field "
585  "and the test space must be H(Div)";
586  }
587 
588  virtual void CalcTestShape(const FiniteElement & test_fe,
590  Vector & shape)
591  {
592  test_fe.CalcPhysDivShape(Trans, shape);
593  shape *= -1.0;
594  }
595 };
596 
597 /** Class for integrating the bilinear form a(u,v) := (Q curl u, v) in 2D where
598  Q is an optional scalar coefficient, u is in H(Curl), and v is in L2 or
599  H1. */
601 {
602 public:
605  : MixedScalarIntegrator(q) {}
606 
607 protected:
608  inline virtual bool VerifyFiniteElementTypes(
609  const FiniteElement & trial_fe,
610  const FiniteElement & test_fe) const
611  {
612  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
613  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
615  }
616 
617  inline virtual const char * FiniteElementTypeFailureMessage() const
618  {
619  return "MixedScalarCurlIntegrator: "
620  "Trial must be H(Curl) and the test space must be a "
621  "scalar field";
622  }
623 
624  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
625  const FiniteElement & test_fe,
627  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
628 
629  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
631  Vector & shape)
632  {
633  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
634  trial_fe.CalcPhysCurlShape(Trans, dshape);
635  }
636 };
637 
638 /** Class for integrating the bilinear form a(u,v) := (Q u, curl v) in 2D where
639  Q is an optional scalar coefficient, u is in L2 or H1, and v is in
640  H(Curl). */
642 {
643 public:
646  : MixedScalarIntegrator(q) {}
647 
648 protected:
649  inline virtual bool VerifyFiniteElementTypes(
650  const FiniteElement & trial_fe,
651  const FiniteElement & test_fe) const
652  {
653  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
656  }
657 
658  inline virtual const char * FiniteElementTypeFailureMessage() const
659  {
660  return "MixedScalarWeakCurlIntegrator: "
661  "Trial space must be a scalar field "
662  "and the test space must be H(Curl)";
663  }
664 
665  inline virtual void CalcTestShape(const FiniteElement & test_fe,
667  Vector & shape)
668  {
669  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
670  test_fe.CalcPhysCurlShape(Trans, dshape);
671  }
672 };
673 
674 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 2D or
675  3D and where Q is an optional coefficient (of type scalar, matrix, or
676  diagonal matrix) u and v are each in H(Curl) or H(Div). */
678 {
679 public:
682  : MixedVectorIntegrator(q) { same_calc_shape = true; }
684  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
686  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
687 };
688 
689 /** Class for integrating the bilinear form a(u,v) := (V x u, v) in 3D and where
690  V is a vector coefficient u and v are each in H(Curl) or H(Div). */
692 {
693 public:
695  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
696 };
697 
698 /** Class for integrating the bilinear form a(u,v) := (V . u, v) in 2D or 3D and
699  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in H1 or
700  L2. */
702 {
703 public:
705  : MixedScalarVectorIntegrator(vq, true) {}
706 
707  inline virtual bool VerifyFiniteElementTypes(
708  const FiniteElement & trial_fe,
709  const FiniteElement & test_fe) const
710  {
711  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
713  }
714 
715  inline virtual const char * FiniteElementTypeFailureMessage() const
716  {
717  return "MixedDotProductIntegrator: "
718  "Trial space must be a vector field "
719  "and the test space must be a scalar field";
720  }
721 };
722 
723 /** Class for integrating the bilinear form a(u,v) := (-V . u, Div v) in 2D or
724  3D and where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
725  RT. */
727 {
728 public:
730  : MixedScalarVectorIntegrator(vq, true) {}
731 
732  inline virtual bool VerifyFiniteElementTypes(
733  const FiniteElement & trial_fe,
734  const FiniteElement & test_fe) const
735  {
736  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
738  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
739  }
740 
741  inline virtual const char * FiniteElementTypeFailureMessage() const
742  {
743  return "MixedWeakGradDotIntegrator: "
744  "Trial space must be a vector field "
745  "and the test space must be a vector field with a divergence";
746  }
747 
748  inline virtual void CalcShape(const FiniteElement & scalar_fe,
750  Vector & shape)
751  { scalar_fe.CalcPhysDivShape(Trans, shape); shape *= -1.0; }
752 };
753 
754 /** Class for integrating the bilinear form a(u,v) := (V x u, Grad v) in 3D and
755  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in H1. */
757 {
758 public:
760  : MixedVectorIntegrator(vq, false) {}
761 
762  inline virtual bool VerifyFiniteElementTypes(
763  const FiniteElement & trial_fe,
764  const FiniteElement & test_fe) const
765  {
766  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
770  }
771 
772  inline virtual const char * FiniteElementTypeFailureMessage() const
773  {
774  return "MixedWeakDivCrossIntegrator: "
775  "Trial space must be a vector field in 3D "
776  "and the test space must be a scalar field with a gradient";
777  }
778 
779  inline virtual void CalcTestShape(const FiniteElement & test_fe,
781  DenseMatrix & shape)
782  { test_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
783 };
784 
785 /** Class for integrating the bilinear form a(u,v) := (Q Grad u, Grad v) in 3D
786  or in 2D and where Q is a scalar or matrix coefficient u and v are both in
787  H1. */
789 {
790 public:
793  : MixedVectorIntegrator(q) { same_calc_shape = true; }
795  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
797  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
798 
799  inline virtual bool VerifyFiniteElementTypes(
800  const FiniteElement & trial_fe,
801  const FiniteElement & test_fe) const
802  {
803  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
804  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
807  }
808 
809  inline virtual const char * FiniteElementTypeFailureMessage() const
810  {
811  return "MixedGradGradIntegrator: "
812  "Trial and test spaces must both be scalar fields "
813  "with a gradient operator.";
814  }
815 
816  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
817  const FiniteElement & test_fe,
819  {
820  // Same as DiffusionIntegrator
821  return test_fe.Space() == FunctionSpace::Pk ?
822  trial_fe.GetOrder() + test_fe.GetOrder() - 2 :
823  trial_fe.GetOrder() + test_fe.GetOrder() + test_fe.GetDim() - 1;
824  }
825 
826  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
828  DenseMatrix & shape)
829  { trial_fe.CalcPhysDShape(Trans, shape); }
830 
831  inline virtual void CalcTestShape(const FiniteElement & test_fe,
833  DenseMatrix & shape)
834  { test_fe.CalcPhysDShape(Trans, shape); }
835 };
836 
837 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, Grad v) in 3D
838  or in 2D and where V is a vector coefficient u and v are both in H1. */
840 {
841 public:
843  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
844 
845  inline virtual bool VerifyFiniteElementTypes(
846  const FiniteElement & trial_fe,
847  const FiniteElement & test_fe) const
848  {
849  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
850  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
853  }
854 
855  inline virtual const char * FiniteElementTypeFailureMessage() const
856  {
857  return "MixedCrossGradGradIntegrator: "
858  "Trial and test spaces must both be scalar fields "
859  "with a gradient operator.";
860  }
861 
862  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
864  DenseMatrix & shape)
865  { trial_fe.CalcPhysDShape(Trans, shape); }
866 
867  inline virtual void CalcTestShape(const FiniteElement & test_fe,
869  DenseMatrix & shape)
870  { test_fe.CalcPhysDShape(Trans, shape); }
871 };
872 
873 /** Class for integrating the bilinear form a(u,v) := (Q Curl u, Curl v) in 3D
874  and where Q is a scalar or matrix coefficient u and v are both in
875  H(Curl). */
877 {
878 public:
881  : MixedVectorIntegrator(q) { same_calc_shape = true; }
883  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
885  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
886 
887  inline virtual bool VerifyFiniteElementTypes(
888  const FiniteElement & trial_fe,
889  const FiniteElement & test_fe) const
890  {
891  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
893  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
896  }
897 
898  inline virtual const char * FiniteElementTypeFailureMessage() const
899  {
900  return "MixedCurlCurlIntegrator"
901  "Trial and test spaces must both be vector fields in 3D "
902  "with a curl.";
903  }
904 
905  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
907  DenseMatrix & shape)
908  { trial_fe.CalcPhysCurlShape(Trans, shape); }
909 
910  inline virtual void CalcTestShape(const FiniteElement & test_fe,
912  DenseMatrix & shape)
913  { test_fe.CalcPhysCurlShape(Trans, shape); }
914 };
915 
916 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, Curl v) in 3D
917  and where V is a vector coefficient u and v are both in H(Curl). */
919 {
920 public:
922  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
923 
924  inline virtual bool VerifyFiniteElementTypes(
925  const FiniteElement & trial_fe,
926  const FiniteElement & test_fe) const
927  {
928  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
930  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
933  }
934 
935  inline virtual const char * FiniteElementTypeFailureMessage() const
936  {
937  return "MixedCrossCurlCurlIntegrator: "
938  "Trial and test spaces must both be vector fields in 3D "
939  "with a curl.";
940  }
941 
942  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
944  DenseMatrix & shape)
945  { trial_fe.CalcPhysCurlShape(Trans, shape); }
946 
947  inline virtual void CalcTestShape(const FiniteElement & test_fe,
949  DenseMatrix & shape)
950  { test_fe.CalcPhysCurlShape(Trans, shape); }
951 };
952 
953 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, Grad v) in 3D
954  and where V is a vector coefficient u is in H(Curl) and v is in H1. */
956 {
957 public:
959  : MixedVectorIntegrator(vq, false) {}
960 
961  inline virtual bool VerifyFiniteElementTypes(
962  const FiniteElement & trial_fe,
963  const FiniteElement & test_fe) const
964  {
965  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
967  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
970  }
971 
972  inline virtual const char * FiniteElementTypeFailureMessage() const
973  {
974  return "MixedCrossCurlGradIntegrator"
975  "Trial space must be a vector field in 3D with a curl"
976  "and the test space must be a scalar field with a gradient";
977  }
978 
979  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
981  DenseMatrix & shape)
982  { trial_fe.CalcPhysCurlShape(Trans, shape); }
983 
984  inline virtual void CalcTestShape(const FiniteElement & test_fe,
986  DenseMatrix & shape)
987  { test_fe.CalcPhysDShape(Trans, shape); }
988 };
989 
990 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, Curl v) in 3D
991  and where V is a scalar coefficient u is in H1 and v is in H(Curl). */
993 {
994 public:
996  : MixedVectorIntegrator(vq, false) {}
997 
998  inline virtual bool VerifyFiniteElementTypes(
999  const FiniteElement & trial_fe,
1000  const FiniteElement & test_fe) const
1001  {
1002  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1003  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1004  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1006  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1007  }
1008 
1009  inline virtual const char * FiniteElementTypeFailureMessage() const
1010  {
1011  return "MixedCrossGradCurlIntegrator"
1012  "Trial space must be a scalar field in 3D with a gradient"
1013  "and the test space must be a vector field with a curl";
1014  }
1015 
1016  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1018  DenseMatrix & shape)
1019  { trial_fe.CalcPhysDShape(Trans, shape); }
1020 
1021  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1023  DenseMatrix & shape)
1024  { test_fe.CalcPhysCurlShape(Trans, shape); }
1025 };
1026 
1027 /** Class for integrating the bilinear form a(u,v) := (V x u, Curl v) in 3D and
1028  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
1029  H(Curl). */
1031 {
1032 public:
1034  : MixedVectorIntegrator(vq, false) {}
1035 
1036  inline virtual bool VerifyFiniteElementTypes(
1037  const FiniteElement & trial_fe,
1038  const FiniteElement & test_fe) const
1039  {
1040  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1041  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1043  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1044  }
1045 
1046  inline virtual const char * FiniteElementTypeFailureMessage() const
1047  {
1048  return "MixedWeakCurlCrossIntegrator: "
1049  "Trial space must be a vector field in 3D "
1050  "and the test space must be a vector field with a curl";
1051  }
1052 
1053  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1055  DenseMatrix & shape)
1056  { test_fe.CalcPhysCurlShape(Trans, shape); }
1057 };
1058 
1059 /** Class for integrating the bilinear form a(u,v) := (V x u, Curl v) in 2D and
1060  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
1061  H(Curl). */
1063 {
1064 public:
1066  : MixedScalarVectorIntegrator(vq, true, true) {}
1067 
1068  inline virtual bool VerifyFiniteElementTypes(
1069  const FiniteElement & trial_fe,
1070  const FiniteElement & test_fe) const
1071  {
1072  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1073  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1075  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1076  }
1077 
1078  inline virtual const char * FiniteElementTypeFailureMessage() const
1079  {
1080  return "MixedScalarWeakCurlCrossIntegrator: "
1081  "Trial space must be a vector field in 2D "
1082  "and the test space must be a vector field with a curl";
1083  }
1084 
1085  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1087  Vector & shape)
1088  {
1089  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
1090  scalar_fe.CalcPhysCurlShape(Trans, dshape);
1091  }
1092 };
1093 
1094 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, v) in 3D or
1095  in 2D and where V is a vector coefficient u is in H1 and v is in H(Curl) or
1096  H(Div). */
1098 {
1099 public:
1101  : MixedVectorIntegrator(vq, false) {}
1102 
1103  inline virtual bool VerifyFiniteElementTypes(
1104  const FiniteElement & trial_fe,
1105  const FiniteElement & test_fe) const
1106  {
1107  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1108  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1109  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1111  }
1112 
1113  inline virtual const char * FiniteElementTypeFailureMessage() const
1114  {
1115  return "MixedCrossGradIntegrator: "
1116  "Trial space must be a scalar field with a gradient operator"
1117  " and the test space must be a vector field both in 3D.";
1118  }
1119 
1120  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1122  DenseMatrix & shape)
1123  { trial_fe.CalcPhysDShape(Trans, shape); }
1124 
1125  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1127  DenseMatrix & shape)
1128  { test_fe.CalcVShape(Trans, shape); }
1129 };
1130 
1131 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, v) in 3D and
1132  where V is a vector coefficient u is in H(Curl) and v is in H(Curl) or
1133  H(Div). */
1135 {
1136 public:
1138  : MixedVectorIntegrator(vq, false) {}
1139 
1140  inline virtual bool VerifyFiniteElementTypes(
1141  const FiniteElement & trial_fe,
1142  const FiniteElement & test_fe) const
1143  {
1144  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1145  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1146  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1148  }
1149 
1150  inline virtual const char * FiniteElementTypeFailureMessage() const
1151  {
1152  return "MixedCrossCurlIntegrator: "
1153  "Trial space must be a vector field in 3D with a curl "
1154  "and the test space must be a vector field";
1155  }
1156 
1157  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1159  DenseMatrix & shape)
1160  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1161 };
1162 
1163 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, v) in 2D and
1164  where V is a vector coefficient u is in H(Curl) and v is in H(Curl) or
1165  H(Div). */
1167 {
1168 public:
1170  : MixedScalarVectorIntegrator(vq, false, true) {}
1171 
1172  inline virtual bool VerifyFiniteElementTypes(
1173  const FiniteElement & trial_fe,
1174  const FiniteElement & test_fe) const
1175  {
1176  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1177  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1178  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1180  }
1181 
1182  inline virtual const char * FiniteElementTypeFailureMessage() const
1183  {
1184  return "MixedCrossCurlIntegrator: "
1185  "Trial space must be a vector field in 2D with a curl "
1186  "and the test space must be a vector field";
1187  }
1188 
1189  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1191  Vector & shape)
1192  {
1193  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
1194  scalar_fe.CalcPhysCurlShape(Trans, dshape); shape *= -1.0;
1195  }
1196 };
1197 
1198 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, v) in 2D and
1199  where V is a vector coefficient u is in H1 and v is in H1 or L2. */
1201 {
1202 public:
1204  : MixedScalarVectorIntegrator(vq, true, true) {}
1205 
1206  inline virtual bool VerifyFiniteElementTypes(
1207  const FiniteElement & trial_fe,
1208  const FiniteElement & test_fe) const
1209  {
1210  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1211  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1212  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1214  }
1215 
1216  inline virtual const char * FiniteElementTypeFailureMessage() const
1217  {
1218  return "MixedScalarCrossGradIntegrator: "
1219  "Trial space must be a scalar field in 2D with a gradient "
1220  "and the test space must be a scalar field";
1221  }
1222 
1223  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1225  DenseMatrix & shape)
1226  { vector_fe.CalcPhysDShape(Trans, shape); }
1227 };
1228 
1229 /** Class for integrating the bilinear form a(u,v) := (V x u, v) in 2D and where
1230  V is a vector coefficient u is in ND or RT and v is in H1 or L2. */
1232 {
1233 public:
1235  : MixedScalarVectorIntegrator(vq, true, true) {}
1236 
1237  inline virtual bool VerifyFiniteElementTypes(
1238  const FiniteElement & trial_fe,
1239  const FiniteElement & test_fe) const
1240  {
1241  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1242  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1244  }
1245 
1246  inline virtual const char * FiniteElementTypeFailureMessage() const
1247  {
1248  return "MixedScalarCrossProductIntegrator: "
1249  "Trial space must be a vector field in 2D "
1250  "and the test space must be a scalar field";
1251  }
1252 };
1253 
1254 /** Class for integrating the bilinear form a(u,v) := (V x z u, v) in 2D and
1255  where V is a vector coefficient u is in H1 or L2 and v is in ND or RT. */
1257 {
1258 public:
1260  : MixedScalarVectorIntegrator(vq, false, true) {}
1261 
1262  inline virtual bool VerifyFiniteElementTypes(
1263  const FiniteElement & trial_fe,
1264  const FiniteElement & test_fe) const
1265  {
1266  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1267  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1269  }
1270 
1271  inline virtual const char * FiniteElementTypeFailureMessage() const
1272  {
1273  return "MixedScalarWeakCrossProductIntegrator: "
1274  "Trial space must be a scalar field in 2D "
1275  "and the test space must be a vector field";
1276  }
1277 
1278  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1280  Vector & shape)
1281  { scalar_fe.CalcPhysShape(Trans, shape); shape *= -1.0; }
1282 };
1283 
1284 /** Class for integrating the bilinear form a(u,v) := (V . Grad u, v) in 2D or
1285  3D and where V is a vector coefficient, u is in H1 and v is in H1 or L2. */
1287 {
1288 public:
1290  : MixedScalarVectorIntegrator(vq, true) {}
1291 
1292  inline virtual bool VerifyFiniteElementTypes(
1293  const FiniteElement & trial_fe,
1294  const FiniteElement & test_fe) const
1295  {
1296  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1297  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1299  }
1300 
1301  inline virtual const char * FiniteElementTypeFailureMessage() const
1302  {
1303  return "MixedDirectionalDerivativeIntegrator: "
1304  "Trial space must be a scalar field with a gradient "
1305  "and the test space must be a scalar field";
1306  }
1307 
1308  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1310  DenseMatrix & shape)
1311  { vector_fe.CalcPhysDShape(Trans, shape); }
1312 };
1313 
1314 /** Class for integrating the bilinear form a(u,v) := (-V . Grad u, Div v) in 2D
1315  or 3D and where V is a vector coefficient, u is in H1 and v is in RT. */
1317 {
1318 public:
1320  : MixedScalarVectorIntegrator(vq, true) {}
1321 
1322  inline virtual bool VerifyFiniteElementTypes(
1323  const FiniteElement & trial_fe,
1324  const FiniteElement & test_fe) const
1325  {
1326  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1327  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1329  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
1330  }
1331 
1332  inline virtual const char * FiniteElementTypeFailureMessage() const
1333  {
1334  return "MixedGradDivIntegrator: "
1335  "Trial space must be a scalar field with a gradient"
1336  "and the test space must be a vector field with a divergence";
1337  }
1338 
1339  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1341  DenseMatrix & shape)
1342  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1343 
1344  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1346  Vector & shape)
1347  { scalar_fe.CalcPhysDivShape(Trans, shape); }
1348 };
1349 
1350 /** Class for integrating the bilinear form a(u,v) := (-V Div u, Grad v) in 2D
1351  or 3D and where V is a vector coefficient, u is in RT and v is in H1. */
1353 {
1354 public:
1356  : MixedScalarVectorIntegrator(vq, false) {}
1357 
1358  inline virtual bool VerifyFiniteElementTypes(
1359  const FiniteElement & trial_fe,
1360  const FiniteElement & test_fe) const
1361  {
1362  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1363  trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
1366  );
1367  }
1368 
1369  inline virtual const char * FiniteElementTypeFailureMessage() const
1370  {
1371  return "MixedDivGradIntegrator: "
1372  "Trial space must be a vector field with a divergence"
1373  "and the test space must be a scalar field with a gradient";
1374  }
1375 
1376  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1378  DenseMatrix & shape)
1379  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1380 
1381  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1383  Vector & shape)
1384  { scalar_fe.CalcPhysDivShape(Trans, shape); }
1385 };
1386 
1387 /** Class for integrating the bilinear form a(u,v) := (-V u, Grad v) in 2D or 3D
1388  and where V is a vector coefficient, u is in H1 and v is in H1. */
1390 {
1391 public:
1393  : MixedScalarVectorIntegrator(vq, false) {}
1394 
1395  inline virtual bool VerifyFiniteElementTypes(
1396  const FiniteElement & trial_fe,
1397  const FiniteElement & test_fe) const
1398  {
1399  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1401  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1402  }
1403 
1404  inline virtual const char * FiniteElementTypeFailureMessage() const
1405  {
1406  return "MixedScalarWeakDivergenceIntegrator: "
1407  "Trial space must be a scalar field "
1408  "and the test space must be a scalar field with a gradient";
1409  }
1410 
1411  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1413  DenseMatrix & shape)
1414  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1415 };
1416 
1417 /** Class for integrating the bilinear form a(u,v) := (Q grad u, v) in either 2D
1418  or 3D and where Q is an optional coefficient (of type scalar, matrix, or
1419  diagonal matrix) u is in H1 and v is in H(Curl) or H(Div). */
1421 {
1422 public:
1425  : MixedVectorIntegrator(q) {}
1427  : MixedVectorIntegrator(dq, true) {}
1429  : MixedVectorIntegrator(mq) {}
1430 
1431 protected:
1432  inline virtual bool VerifyFiniteElementTypes(
1433  const FiniteElement & trial_fe,
1434  const FiniteElement & test_fe) const
1435  {
1436  return (trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1438  }
1439 
1440  inline virtual const char * FiniteElementTypeFailureMessage() const
1441  {
1442  return "MixedVectorGradientIntegrator: "
1443  "Trial spaces must be H1 and the test space must be a "
1444  "vector field in 2D or 3D";
1445  }
1446 
1447  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1449  DenseMatrix & shape)
1450  {
1451  trial_fe.CalcPhysDShape(Trans, shape);
1452  }
1453 
1454 private:
1455  DenseMatrix Jinv;
1456 };
1457 
1458 /** Class for integrating the bilinear form a(u,v) := (Q curl u, v) in 3D and
1459  where Q is an optional coefficient (of type scalar, matrix, or diagonal
1460  matrix) u is in H(Curl) and v is in H(Div) or H(Curl). */
1462 {
1463 public:
1466  : MixedVectorIntegrator(q) {}
1468  : MixedVectorIntegrator(dq, true) {}
1470  : MixedVectorIntegrator(mq) {}
1471 
1472 protected:
1473  inline virtual bool VerifyFiniteElementTypes(
1474  const FiniteElement & trial_fe,
1475  const FiniteElement & test_fe) const
1476  {
1477  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1478  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1480  }
1481 
1482  inline virtual const char * FiniteElementTypeFailureMessage() const
1483  {
1484  return "MixedVectorCurlIntegrator: "
1485  "Trial space must be H(Curl) and the test space must be a "
1486  "vector field in 3D";
1487  }
1488 
1489  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1491  DenseMatrix & shape)
1492  {
1493  trial_fe.CalcPhysCurlShape(Trans, shape);
1494  }
1495 };
1496 
1497 /** Class for integrating the bilinear form a(u,v) := (Q u, curl v) in 3D and
1498  where Q is an optional coefficient (of type scalar, matrix, or diagonal
1499  matrix) u is in H(Div) or H(Curl) and v is in H(Curl). */
1501 {
1502 public:
1505  : MixedVectorIntegrator(q) {}
1507  : MixedVectorIntegrator(dq, true) {}
1509  : MixedVectorIntegrator(mq) {}
1510 
1511 protected:
1512  inline virtual bool VerifyFiniteElementTypes(
1513  const FiniteElement & trial_fe,
1514  const FiniteElement & test_fe) const
1515  {
1516  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1517  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1518  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1519  }
1520 
1521  inline virtual const char * FiniteElementTypeFailureMessage() const
1522  {
1523  return "MixedVectorWeakCurlIntegrator: "
1524  "Trial space must be vector field in 3D and the "
1525  "test space must be H(Curl)";
1526  }
1527 
1528  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1530  DenseMatrix & shape)
1531  {
1532  test_fe.CalcPhysCurlShape(Trans, shape);
1533  }
1534 };
1535 
1536 /** Class for integrating the bilinear form a(u,v) := - (Q u, grad v) in either
1537  2D or 3D and where Q is an optional coefficient (of type scalar, matrix, or
1538  diagonal matrix) u is in H(Div) or H(Curl) and v is in H1. */
1540 {
1541 public:
1544  : MixedVectorIntegrator(q) {}
1546  : MixedVectorIntegrator(dq, true) {}
1548  : MixedVectorIntegrator(mq) {}
1549 
1550 protected:
1551  inline virtual bool VerifyFiniteElementTypes(
1552  const FiniteElement & trial_fe,
1553  const FiniteElement & test_fe) const
1554  {
1555  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1556  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1557  }
1558 
1559  inline virtual const char * FiniteElementTypeFailureMessage() const
1560  {
1561  return "MixedVectorWeakDivergenceIntegrator: "
1562  "Trial space must be vector field and the "
1563  "test space must be H1";
1564  }
1565 
1566  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1568  DenseMatrix & shape)
1569  {
1570  test_fe.CalcPhysDShape(Trans, shape);
1571  shape *= -1.0;
1572  }
1573 };
1574 
1575 /** Class for integrating the bilinear form a(u,v) := (Q grad u, grad v) where Q
1576  can be a scalar or a matrix coefficient. */
1578 {
1579 private:
1580  Vector vec, pointflux, shape;
1581 #ifndef MFEM_THREAD_SAFE
1582  DenseMatrix dshape, dshapedxt, invdfdx, mq;
1583  DenseMatrix te_dshape, te_dshapedxt;
1584 #endif
1585  Coefficient *Q;
1586  MatrixCoefficient *MQ;
1587 
1588 public:
1589  /// Construct a diffusion integrator with coefficient Q = 1
1590  DiffusionIntegrator() { Q = NULL; MQ = NULL; }
1591 
1592  /// Construct a diffusion integrator with a scalar coefficient q
1593  DiffusionIntegrator (Coefficient &q) : Q(&q) { MQ = NULL; }
1594 
1595  /// Construct a diffusion integrator with a matrix coefficient q
1596  DiffusionIntegrator (MatrixCoefficient &q) : MQ(&q) { Q = NULL; }
1597 
1598  /** Given a particular Finite Element
1599  computes the element stiffness matrix elmat. */
1600  virtual void AssembleElementMatrix(const FiniteElement &el,
1602  DenseMatrix &elmat);
1603  /** Given a trial and test Finite Element computes the element stiffness
1604  matrix elmat. */
1605  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1606  const FiniteElement &test_fe,
1608  DenseMatrix &elmat);
1609 
1610  /// Perform the local action of the BilinearFormIntegrator
1611  virtual void AssembleElementVector(const FiniteElement &el,
1613  const Vector &elfun, Vector &elvect);
1614 
1615  virtual void ComputeElementFlux(const FiniteElement &el,
1617  Vector &u, const FiniteElement &fluxelem,
1618  Vector &flux, int with_coef = 1);
1619 
1620  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
1622  Vector &flux, Vector *d_energy = NULL);
1623 };
1624 
1625 /** Class for local mass matrix assembling a(u,v) := (Q u, v) */
1627 {
1628 protected:
1631 
1632 public:
1634  : BilinearFormIntegrator(ir) { Q = NULL; }
1635  /// Construct a mass integrator with coefficient q
1637  : BilinearFormIntegrator(ir), Q(&q) { }
1638 
1639  /** Given a particular Finite Element
1640  computes the element mass matrix elmat. */
1641  virtual void AssembleElementMatrix(const FiniteElement &el,
1643  DenseMatrix &elmat);
1644  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1645  const FiniteElement &test_fe,
1647  DenseMatrix &elmat);
1648 };
1649 
1651 {
1652 public:
1654 
1656 
1657  virtual void AssembleFaceMatrix(const FiniteElement &el1,
1658  const FiniteElement &el2,
1660  DenseMatrix &elmat);
1661 };
1662 
1663 /// alpha (q . grad u, v)
1665 {
1666 private:
1667  DenseMatrix dshape, adjJ, Q_ir;
1668  Vector shape, vec2, BdFidxT;
1669  VectorCoefficient &Q;
1670  double alpha;
1671 
1672 public:
1674  : Q(q) { alpha = a; }
1675  virtual void AssembleElementMatrix(const FiniteElement &,
1677  DenseMatrix &);
1678 };
1679 
1680 /// alpha (q . grad u, v) using the "group" FE discretization
1682 {
1683 private:
1684  DenseMatrix dshape, adjJ, Q_nodal, grad;
1685  Vector shape;
1686  VectorCoefficient &Q;
1687  double alpha;
1688 
1689 public:
1691  : Q(q) { alpha = a; }
1692  virtual void AssembleElementMatrix(const FiniteElement &,
1694  DenseMatrix &);
1695 };
1696 
1697 /** Class for integrating the bilinear form a(u,v) := (Q u, v),
1698  where u=(u1,...,un) and v=(v1,...,vn); ui and vi are defined
1699  by scalar FE through standard transformation. */
1701 {
1702 private:
1703  Vector shape, te_shape, vec;
1704  DenseMatrix partelmat;
1705  DenseMatrix mcoeff;
1706  Coefficient *Q;
1707  VectorCoefficient *VQ;
1708  MatrixCoefficient *MQ;
1709 
1710  int Q_order;
1711 
1712 public:
1713  /// Construct an integrator with coefficient 1.0
1715  { Q = NULL; VQ = NULL; MQ = NULL; Q_order = 0; }
1716  /** Construct an integrator with scalar coefficient q.
1717  If possible, save memory by using a scalar integrator since
1718  the resulting matrix is block diagonal with the same diagonal
1719  block repeated. */
1721  : Q(&q) { VQ = NULL; MQ = NULL; Q_order = qo; }
1723  : BilinearFormIntegrator(ir), Q(&q)
1724  { VQ = NULL; MQ = NULL; Q_order = 0; }
1725  /// Construct an integrator with diagonal coefficient q
1727  : VQ(&q) { Q = NULL; MQ = NULL; Q_order = qo; }
1728  /// Construct an integrator with matrix coefficient q
1730  : MQ(&q) { Q = NULL; VQ = NULL; Q_order = qo; }
1731 
1732  virtual void AssembleElementMatrix(const FiniteElement &el,
1734  DenseMatrix &elmat);
1735  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1736  const FiniteElement &test_fe,
1738  DenseMatrix &elmat);
1739 };
1740 
1741 
1742 /** Class for integrating (div u, p) where u is a vector field given
1743  by VectorFiniteElement through Piola transformation (for RT
1744  elements); p is scalar function given by FiniteElement through
1745  standard transformation. Here, u is the trial function and p is
1746  the test function.
1747  Note: the element matrix returned by AssembleElementMatrix2
1748  does NOT depend on the ElementTransformation Trans. */
1750 {
1751 private:
1752  Coefficient *Q;
1753 #ifndef MFEM_THREAD_SAFE
1754  Vector divshape, shape;
1755 #endif
1756 public:
1759  virtual void AssembleElementMatrix(const FiniteElement &el,
1761  DenseMatrix &elmat) { }
1762  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1763  const FiniteElement &test_fe,
1765  DenseMatrix &elmat);
1766 };
1767 
1768 
1769 /** Integrator for `(-Q u, grad v)` for Nedelec (`u`) and H1 (`v`) elements.
1770  This is equivalent to a weak divergence of the Nedelec basis functions. */
1772 {
1773 private:
1774  Coefficient *Q;
1775 #ifndef MFEM_THREAD_SAFE
1776  DenseMatrix dshape;
1777  DenseMatrix dshapedxt;
1778  DenseMatrix vshape;
1779  DenseMatrix invdfdx;
1780 #endif
1781 public:
1784  virtual void AssembleElementMatrix(const FiniteElement &el,
1786  DenseMatrix &elmat) { }
1787  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1788  const FiniteElement &test_fe,
1790  DenseMatrix &elmat);
1791 };
1792 
1793 /** Integrator for (curl u, v) for Nedelec and RT elements. If the trial and
1794  test spaces are switched, assembles the form (u, curl v). */
1796 {
1797 private:
1798  Coefficient *Q;
1799 #ifndef MFEM_THREAD_SAFE
1800  DenseMatrix curlshapeTrial;
1801  DenseMatrix vshapeTest;
1802  DenseMatrix curlshapeTrial_dFT;
1803 #endif
1804 public:
1805  VectorFECurlIntegrator() { Q = NULL; }
1807  virtual void AssembleElementMatrix(const FiniteElement &el,
1809  DenseMatrix &elmat) { }
1810  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1811  const FiniteElement &test_fe,
1813  DenseMatrix &elmat);
1814 };
1815 
1816 
1817 /// Class for integrating (Q D_i(u), v); u and v are scalars
1819 {
1820 private:
1821  Coefficient & Q;
1822  int xi;
1823  DenseMatrix dshape, dshapedxt, invdfdx;
1824  Vector shape, dshapedxi;
1825 public:
1826  DerivativeIntegrator(Coefficient &q, int i) : Q(q), xi(i) { }
1827  virtual void AssembleElementMatrix(const FiniteElement &el,
1829  DenseMatrix &elmat)
1830  { AssembleElementMatrix2(el,el,Trans,elmat); }
1831  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1832  const FiniteElement &test_fe,
1834  DenseMatrix &elmat);
1835 };
1836 
1837 /// Integrator for (curl u, curl v) for Nedelec elements
1839 {
1840 private:
1841  Vector vec, pointflux;
1842 #ifndef MFEM_THREAD_SAFE
1843  DenseMatrix curlshape, curlshape_dFt;
1844  DenseMatrix vshape, projcurl;
1845 #endif
1846  Coefficient *Q;
1847 
1848 public:
1849  CurlCurlIntegrator() { Q = NULL; }
1850  /// Construct a bilinear form integrator for Nedelec elements
1852 
1853  /* Given a particular Finite Element, compute the
1854  element curl-curl matrix elmat */
1855  virtual void AssembleElementMatrix(const FiniteElement &el,
1857  DenseMatrix &elmat);
1858 
1859  virtual void ComputeElementFlux(const FiniteElement &el,
1861  Vector &u, const FiniteElement &fluxelem,
1862  Vector &flux, int with_coef);
1863 
1864  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
1866  Vector &flux, Vector *d_energy = NULL);
1867 };
1868 
1869 /** Integrator for (curl u, curl v) for FE spaces defined by 'dim' copies of a
1870  scalar FE space. */
1872 {
1873 private:
1874 #ifndef MFEM_THREAD_SAFE
1875  DenseMatrix dshape_hat, dshape, curlshape, Jadj, grad_hat, grad;
1876 #endif
1877  Coefficient *Q;
1878 
1879 public:
1881 
1883 
1884  /// Assemble an element matrix
1885  virtual void AssembleElementMatrix(const FiniteElement &el,
1887  DenseMatrix &elmat);
1888  /// Compute element energy: (1/2) (curl u, curl u)_E
1889  virtual double GetElementEnergy(const FiniteElement &el,
1891  const Vector &elfun);
1892 };
1893 
1894 /// Integrator for (Q u, v) for VectorFiniteElements
1896 {
1897 private:
1898  Coefficient *Q;
1899  VectorCoefficient *VQ;
1900  MatrixCoefficient *MQ;
1901  void Init(Coefficient *q, VectorCoefficient *vq, MatrixCoefficient *mq)
1902  { Q = q; VQ = vq; MQ = mq; }
1903 
1904 #ifndef MFEM_THREAD_SAFE
1905  Vector shape;
1906  Vector D;
1907  DenseMatrix K;
1908  DenseMatrix test_vshape;
1909  DenseMatrix trial_vshape;
1910 #endif
1911 
1912 public:
1913  VectorFEMassIntegrator() { Init(NULL, NULL, NULL); }
1914  VectorFEMassIntegrator(Coefficient *_q) { Init(_q, NULL, NULL); }
1915  VectorFEMassIntegrator(Coefficient &q) { Init(&q, NULL, NULL); }
1916  VectorFEMassIntegrator(VectorCoefficient *_vq) { Init(NULL, _vq, NULL); }
1917  VectorFEMassIntegrator(VectorCoefficient &vq) { Init(NULL, &vq, NULL); }
1918  VectorFEMassIntegrator(MatrixCoefficient *_mq) { Init(NULL, NULL, _mq); }
1919  VectorFEMassIntegrator(MatrixCoefficient &mq) { Init(NULL, NULL, &mq); }
1920 
1921  virtual void AssembleElementMatrix(const FiniteElement &el,
1923  DenseMatrix &elmat);
1924  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1925  const FiniteElement &test_fe,
1927  DenseMatrix &elmat);
1928 };
1929 
1930 /** Integrator for (Q div u, p) where u=(v1,...,vn) and all vi are in the same
1931  scalar FE space; p is also in a (different) scalar FE space. */
1933 {
1934 private:
1935  Coefficient *Q;
1936 
1937  Vector shape;
1938  Vector divshape;
1939  DenseMatrix dshape;
1940  DenseMatrix gshape;
1941  DenseMatrix Jadj;
1942 
1943 public:
1947 
1948  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1949  const FiniteElement &test_fe,
1951  DenseMatrix &elmat);
1952 };
1953 
1954 /// (Q div u, div v) for RT elements
1956 {
1957 private:
1958  Coefficient *Q;
1959 
1960 #ifndef MFEM_THREAD_SAFE
1961  Vector divshape;
1962 #endif
1963 
1964 public:
1965  DivDivIntegrator() { Q = NULL; }
1967 
1968  virtual void AssembleElementMatrix(const FiniteElement &el,
1970  DenseMatrix &elmat);
1971 };
1972 
1973 /** Integrator for
1974  (Q grad u, grad v) = sum_i (Q grad u_i, grad v_i) e_i e_i^T
1975  for FE spaces defined by 'dim' copies of a scalar FE space. Where e_i
1976  is the unit vector in the i-th direction. The resulting local element
1977  matrix is a block-diagonal matrix consisting of 'dim' copies of a scalar
1978  diffusion matrix in each diagonal block. */
1980 {
1981 private:
1982  Coefficient *Q;
1983 
1984  DenseMatrix Jinv;
1985  DenseMatrix dshape;
1986  DenseMatrix gshape;
1987  DenseMatrix pelmat;
1988 
1989 public:
1992 
1993  virtual void AssembleElementMatrix(const FiniteElement &el,
1995  DenseMatrix &elmat);
1996  virtual void AssembleElementVector(const FiniteElement &el,
1998  const Vector &elfun, Vector &elvect);
1999 };
2000 
2001 /** Integrator for the linear elasticity form:
2002  a(u,v) = (lambda div(u), div(v)) + (2 mu e(u), e(v)),
2003  where e(v) = (1/2) (grad(v) + grad(v)^T).
2004  This is a 'Vector' integrator, i.e. defined for FE spaces
2005  using multiple copies of a scalar FE space. */
2007 {
2008 private:
2009  double q_lambda, q_mu;
2010  Coefficient *lambda, *mu;
2011 
2012 #ifndef MFEM_THREAD_SAFE
2013  DenseMatrix dshape, Jinv, gshape, pelmat;
2014  Vector divshape;
2015 #endif
2016 
2017 public:
2019  { lambda = &l; mu = &m; }
2020  /** With this constructor lambda = q_l * m and mu = q_m * m;
2021  if dim * q_l + 2 * q_m = 0 then trace(sigma) = 0. */
2022  ElasticityIntegrator(Coefficient &m, double q_l, double q_m)
2023  { lambda = NULL; mu = &m; q_lambda = q_l; q_mu = q_m; }
2024 
2025  virtual void AssembleElementMatrix(const FiniteElement &,
2027  DenseMatrix &);
2028 };
2029 
2030 /** Integrator for the DG form:
2031  alpha < rho_u (u.n) {v},[w] > + beta < rho_u |u.n| [v],[w] >,
2032  where v and w are the trial and test variables, respectively, and rho/u are
2033  given scalar/vector coefficients. The vector coefficient, u, is assumed to
2034  be continuous across the faces and when given the scalar coefficient, rho,
2035  is assumed to be discontinuous. The integrator uses the upwind value of rho,
2036  rho_u, which is value from the side into which the vector coefficient, u,
2037  points. */
2039 {
2040 private:
2041  Coefficient *rho;
2042  VectorCoefficient *u;
2043  double alpha, beta;
2044 
2045  Vector shape1, shape2;
2046 
2047 public:
2048  /// Construct integrator with rho = 1.
2049  DGTraceIntegrator(VectorCoefficient &_u, double a, double b)
2050  { rho = NULL; u = &_u; alpha = a; beta = b; }
2051 
2053  double a, double b)
2054  { rho = &_rho; u = &_u; alpha = a; beta = b; }
2055 
2057  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2058  const FiniteElement &el2,
2060  DenseMatrix &elmat);
2061 };
2062 
2063 /** Integrator for the DG form:
2064 
2065  - < {(Q grad(u)).n}, [v] > + sigma < [u], {(Q grad(v)).n} >
2066  + kappa < {h^{-1} Q} [u], [v] >,
2067 
2068  where Q is a scalar or matrix diffusion coefficient and u, v are the trial
2069  and test spaces, respectively. The parameters sigma and kappa determine the
2070  DG method to be used (when this integrator is added to the "broken"
2071  DiffusionIntegrator):
2072  * sigma = -1, kappa >= kappa0: symm. interior penalty (IP or SIPG) method,
2073  * sigma = +1, kappa > 0: non-symmetric interior penalty (NIPG) method,
2074  * sigma = +1, kappa = 0: the method of Baumann and Oden. */
2076 {
2077 protected:
2080  double sigma, kappa;
2081 
2082  // these are not thread-safe!
2085 
2086 public:
2087  DGDiffusionIntegrator(const double s, const double k)
2088  : Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
2089  DGDiffusionIntegrator(Coefficient &q, const double s, const double k)
2090  : Q(&q), MQ(NULL), sigma(s), kappa(k) { }
2091  DGDiffusionIntegrator(MatrixCoefficient &q, const double s, const double k)
2092  : Q(NULL), MQ(&q), sigma(s), kappa(k) { }
2094  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2095  const FiniteElement &el2,
2097  DenseMatrix &elmat);
2098 };
2099 
2100 /** Integrator for the DG elasticity form, for the formulations see:
2101  - PhD Thesis of Jonas De Basabe, High-Order Finite %Element Methods for
2102  Seismic Wave Propagation, UT Austin, 2009, p. 23, and references therein
2103  - Peter Hansbo and Mats G. Larson, Discontinuous Galerkin and the
2104  Crouzeix-Raviart %Element: Application to Elasticity, PREPRINT 2000-09,
2105  p.3
2106 
2107  \f[
2108  - \left< \{ \tau(u) \}, [v] \right> + \alpha \left< \{ \tau(v) \}, [u]
2109  \right> + \kappa \left< h^{-1} \{ \lambda + 2 \mu \} [u], [v] \right>
2110  \f]
2111 
2112  where \f$ \left<u, v\right> = \int_{F} u \cdot v \f$, and \f$ F \f$ is a
2113  face which is either a boundary face \f$ F_b \f$ of an element \f$ K \f$ or
2114  an interior face \f$ F_i \f$ separating elements \f$ K_1 \f$ and \f$ K_2 \f$.
2115 
2116  In the bilinear form above \f$ \tau(u) \f$ is traction, and it's also
2117  \f$ \tau(u) = \sigma(u) \cdot \vec{n} \f$, where \f$ \sigma(u) \f$ is
2118  stress, and \f$ \vec{n} \f$ is the unit normal vector w.r.t. to \f$ F \f$.
2119 
2120  In other words, we have
2121  \f[
2122  - \left< \{ \sigma(u) \cdot \vec{n} \}, [v] \right> + \alpha \left< \{
2123  \sigma(v) \cdot \vec{n} \}, [u] \right> + \kappa \left< h^{-1} \{
2124  \lambda + 2 \mu \} [u], [v] \right>
2125  \f]
2126 
2127  For isotropic media
2128  \f[
2129  \begin{split}
2130  \sigma(u) &= \lambda \nabla \cdot u I + 2 \mu \varepsilon(u) \\
2131  &= \lambda \nabla \cdot u I + 2 \mu \frac{1}{2} (\nabla u + \nabla
2132  u^T) \\
2133  &= \lambda \nabla \cdot u I + \mu (\nabla u + \nabla u^T)
2134  \end{split}
2135  \f]
2136 
2137  where \f$ I \f$ is identity matrix, \f$ \lambda \f$ and \f$ \mu \f$ are Lame
2138  coefficients (see ElasticityIntegrator), \f$ u, v \f$ are the trial and test
2139  functions, respectively.
2140 
2141  The parameters \f$ \alpha \f$ and \f$ \kappa \f$ determine the DG method to
2142  use (when this integrator is added to the "broken" ElasticityIntegrator):
2143 
2144  - IIPG, \f$\alpha = 0\f$,
2145  C. Dawson, S. Sun, M. Wheeler, Compatible algorithms for coupled flow and
2146  transport, Comp. Meth. Appl. Mech. Eng., 193(23-26), 2565-2580, 2004.
2147 
2148  - SIPG, \f$\alpha = -1\f$,
2149  M. Grote, A. Schneebeli, D. Schotzau, Discontinuous Galerkin Finite
2150  %Element Method for the Wave Equation, SINUM, 44(6), 2408-2431, 2006.
2151 
2152  - NIPG, \f$\alpha = 1\f$,
2153  B. Riviere, M. Wheeler, V. Girault, A Priori Error Estimates for Finite
2154  %Element Methods Based on Discontinuous Approximation Spaces for Elliptic
2155  Problems, SINUM, 39(3), 902-931, 2001.
2156 
2157  This is a '%Vector' integrator, i.e. defined for FE spaces using multiple
2158  copies of a scalar FE space.
2159  */
2161 {
2162 public:
2163  DGElasticityIntegrator(double alpha_, double kappa_)
2164  : lambda(NULL), mu(NULL), alpha(alpha_), kappa(kappa_) {}
2165 
2167  double alpha_, double kappa_)
2168  : lambda(&lambda_), mu(&mu_), alpha(alpha_), kappa(kappa_) {}
2169 
2171  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2172  const FiniteElement &el2,
2174  DenseMatrix &elmat);
2175 
2176 protected:
2178  double alpha, kappa;
2179 
2180 #ifndef MFEM_THREAD_SAFE
2181  // values of all scalar basis functions for one component of u (which is a
2182  // vector) at the integration point in the reference space
2184  // values of derivatives of all scalar basis functions for one component
2185  // of u (which is a vector) at the integration point in the reference space
2187  // Adjugate of the Jacobian of the transformation: adjJ = det(J) J^{-1}
2189  // gradient of shape functions in the real (physical, not reference)
2190  // coordinates, scaled by det(J):
2191  // dshape_ps(jdof,jm) = sum_{t} adjJ(t,jm)*dshape(jdof,t)
2193  Vector nor; // nor = |weight(J_face)| n
2194  Vector nL1, nL2; // nL1 = (lambda1 * ip.weight / detJ1) nor
2195  Vector nM1, nM2; // nM1 = (mu1 * ip.weight / detJ1) nor
2196  Vector dshape1_dnM, dshape2_dnM; // dshape1_dnM = dshape1_ps . nM1
2197  // 'jmat' corresponds to the term: kappa <h^{-1} {lambda + 2 mu} [u], [v]>
2199 #endif
2200 
2201  static void AssembleBlock(
2202  const int dim, const int row_ndofs, const int col_ndofs,
2203  const int row_offset, const int col_offset,
2204  const double jmatcoef, const Vector &col_nL, const Vector &col_nM,
2205  const Vector &row_shape, const Vector &col_shape,
2206  const Vector &col_dshape_dnM, const DenseMatrix &col_dshape,
2207  DenseMatrix &elmat, DenseMatrix &jmat);
2208 };
2209 
2210 /** Integrator for the DPG form: < v, [w] > over all faces (the interface) where
2211  the trial variable v is defined on the interface and the test variable w is
2212  defined inside the elements, generally in a DG space. */
2214 {
2215 private:
2216  Vector face_shape, shape1, shape2;
2217 
2218 public:
2221  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
2222  const FiniteElement &test_fe1,
2223  const FiniteElement &test_fe2,
2225  DenseMatrix &elmat);
2226 };
2227 
2228 /** Integrator for the form: < v, [w.n] > over all faces (the interface) where
2229  the trial variable v is defined on the interface and the test variable w is
2230  in an H(div)-conforming space. */
2232 {
2233 private:
2234  Vector face_shape, normal, shape1_n, shape2_n;
2235  DenseMatrix shape1, shape2;
2236 
2237 public:
2240  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
2241  const FiniteElement &test_fe1,
2242  const FiniteElement &test_fe2,
2244  DenseMatrix &elmat);
2245 };
2246 
2247 /** Abstract class to serve as a base for local interpolators to be used in the
2248  DiscreteLinearOperator class. */
2250 
2251 
2252 /** Class for constructing the gradient as a DiscreteLinearOperator from an
2253  H1-conforming space to an H(curl)-conforming space. The range space can be
2254  vector L2 space as well. */
2256 {
2257 public:
2258  virtual void AssembleElementMatrix2(const FiniteElement &h1_fe,
2259  const FiniteElement &nd_fe,
2261  DenseMatrix &elmat)
2262  { nd_fe.ProjectGrad(h1_fe, Trans, elmat); }
2263 };
2264 
2265 
2266 /** Class for constructing the identity map as a DiscreteLinearOperator. This
2267  is the discrete embedding matrix when the domain space is a subspace of
2268  the range space. Otherwise, a dof projection matrix is constructed. */
2270 {
2271 public:
2272  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2273  const FiniteElement &ran_fe,
2275  DenseMatrix &elmat)
2276  { ran_fe.Project(dom_fe, Trans, elmat); }
2277 };
2278 
2279 
2280 /** Class for constructing the (local) discrete curl matrix which can be used
2281  as an integrator in a DiscreteLinearOperator object to assemble the global
2282  discrete curl matrix. */
2284 {
2285 public:
2286  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2287  const FiniteElement &ran_fe,
2289  DenseMatrix &elmat)
2290  { ran_fe.ProjectCurl(dom_fe, Trans, elmat); }
2291 };
2292 
2293 
2294 /** Class for constructing the (local) discrete divergence matrix which can
2295  be used as an integrator in a DiscreteLinearOperator object to assemble
2296  the global discrete divergence matrix.
2297 
2298  Note: Since the dofs in the L2_FECollection are nodal values, the local
2299  discrete divergence matrix (with an RT-type domain space) will depend on
2300  the transformation. On the other hand, the local matrix returned by
2301  VectorFEDivergenceIntegrator is independent of the transformation. */
2303 {
2304 public:
2305  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2306  const FiniteElement &ran_fe,
2308  DenseMatrix &elmat)
2309  { ran_fe.ProjectDiv(dom_fe, Trans, elmat); }
2310 };
2311 
2312 
2313 /** A trace face interpolator class for interpolating the normal component of
2314  the domain space, e.g. vector H1, into the range space, e.g. the trace of
2315  RT which uses FiniteElement::INTEGRAL map type. */
2317 {
2318 public:
2319  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2320  const FiniteElement &ran_fe,
2322  DenseMatrix &elmat);
2323 };
2324 
2325 }
2326 
2327 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:46
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:462
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
Definition: bilininteg.hpp:182
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:855
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
Definition: bilininteg.hpp:249
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:998
VectorFEMassIntegrator(Coefficient &q)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:65
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:961
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:115
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:617
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:588
MixedCrossCurlGradIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:958
DGTraceIntegrator(VectorCoefficient &_u, double a, double b)
Construct integrator with rho = 1.
MixedGradDivIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
DGDiffusionIntegrator(Coefficient &q, const double s, const double k)
Integrator defining a sum of multiple Integrators.
Definition: bilininteg.hpp:150
SumIntegrator(int own_integs=1)
Definition: bilininteg.hpp:158
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
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.cpp:38
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
DiffusionIntegrator(Coefficient &q)
Construct a diffusion integrator with a scalar coefficient q.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:270
MixedScalarCrossProductIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:862
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Definition: fe.cpp:156
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:979
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:365
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:115
MixedVectorWeakDivergenceIntegrator(VectorCoefficient &dq)
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
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.cpp:59
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:748
ElasticityIntegrator(Coefficient &l, Coefficient &m)
VectorFEMassIntegrator(VectorCoefficient *_vq)
VectorDiffusionIntegrator(Coefficient &q)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:216
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:443
Integrator for (curl u, curl v) for Nedelec elements.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorCurlCurlIntegrator(Coefficient &q)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:479
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
MixedScalarMassIntegrator(Coefficient &q)
Definition: bilininteg.hpp:403
VectorMassIntegrator(Coefficient &q, const IntegrationRule *ir)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:831
alpha (q . grad u, v) using the &quot;group&quot; FE discretization
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:944
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:124
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:573
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
DivDivIntegrator(Coefficient &q)
int Size() const
Returns the size of the vector.
Definition: vector.hpp:106
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorWeakCurlIntegrator(MatrixCoefficient &mq)
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Polynomials of order k.
Definition: fe.hpp:34
virtual const char * FiniteElementTypeFailureMessage() const
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:910
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:82
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:91
BilinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: bilininteg.hpp:27
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:127
void AddIntegrator(BilinearFormIntegrator *integ)
Definition: bilininteg.hpp:160
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:809
MixedVectorCurlIntegrator(VectorCoefficient &dq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:221
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Definition: bilininteg.hpp:74
MixedCrossProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:694
MixedScalarVectorIntegrator(VectorCoefficient &vq, bool _transpose=false, bool _cross_2d=false)
Definition: bilininteg.hpp:332
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.cpp:140
double * GetData() const
Definition: vector.hpp:114
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorGradientIntegrator(Coefficient &q)
(Q div u, div v) for RT elements
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:732
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:349
MixedScalarCrossGradIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:826
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
Definition: bilininteg.cpp:595
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:805
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:799
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.cpp:73
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:129
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Implements CalcDivShape methods.
Definition: fe.hpp:101
MixedCurlCurlIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:884
MatrixCoefficient * MQ
VectorFEMassIntegrator(VectorCoefficient &vq)
VectorMassIntegrator()
Construct an integrator with coefficient 1.0.
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
VectorMassIntegrator(VectorCoefficient &q, int qo=0)
Construct an integrator with diagonal coefficient q.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:196
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:376
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:297
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:73
MixedGradGradIntegrator(VectorCoefficient &dq)
Definition: bilininteg.hpp:794
Integrator that inverts the matrix assembled by another integrator.
Definition: bilininteg.hpp:132
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:148
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
ElasticityIntegrator(Coefficient &m, double q_l, double q_m)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:721
MixedWeakDivCrossIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:759
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
int dim
Definition: ex3.cpp:47
VectorDivergenceIntegrator(Coefficient &q)
MixedCrossCurlCurlIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:921
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:98
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:984
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedVectorWeakCurlIntegrator(VectorCoefficient &dq)
MixedScalarWeakCurlIntegrator(Coefficient &q)
Definition: bilininteg.hpp:645
MixedCurlCurlIntegrator(Coefficient &q)
Definition: bilininteg.hpp:880
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:500
Implements CalcDShape methods.
Definition: fe.hpp:100
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:105
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: bilininteg.hpp:63
VectorFEDivergenceIntegrator(Coefficient &q)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:905
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:762
void SetIntRule(const IntegrationRule *ir)
Definition: bilininteg.hpp:79
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
VectorFECurlIntegrator(Coefficient &q)
virtual void AssembleElementMatrix2(const FiniteElement &h1_fe, const FiniteElement &nd_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
Definition: bilininteg.hpp:68
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:278
DiffusionIntegrator()
Construct a diffusion integrator with coefficient Q = 1.
MixedVectorMassIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:685
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
MixedScalarDivergenceIntegrator(Coefficient &q)
Definition: bilininteg.hpp:496
MixedVectorIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:267
Implements CalcCurlShape methods.
Definition: fe.hpp:102
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:290
MixedVectorMassIntegrator(Coefficient &q)
Definition: bilininteg.hpp:681
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.cpp:174
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedDirectionalDerivativeIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
MixedCrossGradIntegrator(VectorCoefficient &vq)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:551
virtual const char * FiniteElementTypeFailureMessage() const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:649
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:942
MixedVectorCurlIntegrator(MatrixCoefficient &mq)
MixedCurlCurlIntegrator(VectorCoefficient &dq)
Definition: bilininteg.hpp:882
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:665
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:508
BoundaryMassIntegrator(Coefficient &q)
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:890
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:887
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
Definition: bilininteg.cpp:56
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
Definition: bilininteg.cpp:517
MixedWeakGradDotIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:729
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:715
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:556
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:295
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorDivergenceIntegrator(Coefficient *_q)
MixedVectorIntegrator(Coefficient &q)
Definition: bilininteg.hpp:262
DiffusionIntegrator(MatrixCoefficient &q)
Construct a diffusion integrator with a matrix coefficient q.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:658
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedGradGradIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:796
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:846
VectorFEMassIntegrator(MatrixCoefficient &mq)
MixedVectorMassIntegrator(VectorCoefficient &dq)
Definition: bilininteg.hpp:683
MixedScalarWeakCurlCrossIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:608
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:23
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:536
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
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.cpp:164
virtual const char * FiniteElementTypeFailureMessage() const
MixedScalarWeakDivergenceIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:437
DGDiffusionIntegrator(const double s, const double k)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:779
LumpedIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
Definition: bilininteg.hpp:121
MixedCrossGradCurlIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:995
MixedVectorDivergenceIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:532
MassIntegrator(Coefficient &q, const IntegrationRule *ir=NULL)
Construct a mass integrator with coefficient q.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:204
MixedCrossGradGradIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:842
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:515
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:924
MixedVectorIntegrator(VectorCoefficient &dq, bool diag=true)
Definition: bilininteg.hpp:264
MixedScalarCrossCurlIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Assemble an element matrix.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:39
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:581
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:741
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:772
VectorMassIntegrator(Coefficient &q, int qo=0)
VectorFEMassIntegrator(Coefficient *_q)
DGElasticityIntegrator(double alpha_, double kappa_)
MixedScalarDerivativeIntegrator(Coefficient &q)
Definition: bilininteg.hpp:423
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorFEMassIntegrator(MatrixCoefficient *_mq)
MixedDotProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:704
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorCurlIntegrator(Coefficient &q)
MixedWeakCurlCrossIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:471
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:764
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:188
VectorMassIntegrator(MatrixCoefficient &q, int qo=0)
Construct an integrator with matrix coefficient q.
MixedScalarWeakCrossProductIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:284
GroupConvectionIntegrator(VectorCoefficient &q, double a=1.0)
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:371
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:867
MixedVectorGradientIntegrator(VectorCoefficient &dq)
static void AssembleBlock(const int dim, const int row_ndofs, const int col_ndofs, const int row_offset, const int col_offset, const double jmatcoef, const Vector &col_nL, const Vector &col_nM, const Vector &row_shape, const Vector &col_shape, const Vector &col_dshape_dnM, const DenseMatrix &col_dshape, DenseMatrix &elmat, DenseMatrix &jmat)
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual ~SumIntegrator()
Definition: bilininteg.cpp:118
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
MixedVectorWeakDivergenceIntegrator(MatrixCoefficient &mq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:436
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:898
virtual const char * FiniteElementTypeFailureMessage() const
DGDiffusionIntegrator(MatrixCoefficient &q, const double s, const double k)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:629
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
Compute element energy: (1/2) (curl u, curl u)_E.
MixedCrossCurlIntegrator(VectorCoefficient &vq)
Vector data type.
Definition: vector.hpp:36
int GetDerivType() const
Definition: fe.hpp:135
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:544
MixedVectorProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:413
virtual const char * FiniteElementTypeFailureMessage() const
ConvectionIntegrator(VectorCoefficient &q, double a=1.0)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
MixedScalarCurlIntegrator(Coefficient &q)
Definition: bilininteg.hpp:604
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:31
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:947
DerivativeIntegrator(Coefficient &q, int i)
DGTraceIntegrator(Coefficient &_rho, VectorCoefficient &_u, double a, double b)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:707
CurlCurlIntegrator(Coefficient &q)
Construct a bilinear form integrator for Nedelec elements.
virtual const char * FiniteElementTypeFailureMessage() const
MixedGradGradIntegrator(Coefficient &q)
Definition: bilininteg.hpp:792
MixedVectorGradientIntegrator(MatrixCoefficient &mq)
const IntegrationRule * IntRule
Definition: bilininteg.hpp:25
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:816
virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
int GetRangeType() const
Definition: fe.hpp:129
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:520
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:365
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:972
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:935
MassIntegrator(const IntegrationRule *ir=NULL)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:210
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Definition: bilininteg.cpp:653
Integrator for (Q u, v) for VectorFiniteElements.
MixedVectorWeakCurlIntegrator(Coefficient &q)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:845
TransposeIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
Definition: bilininteg.hpp:93
MixedDivGradIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:427
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:624
MixedScalarIntegrator(Coefficient &q)
Definition: bilininteg.hpp:194
DGElasticityIntegrator(Coefficient &lambda_, Coefficient &mu_, double alpha_, double kappa_)
Class for integrating (Q D_i(u), v); u and v are scalars.
InverseIntegrator(BilinearFormIntegrator *integ, int own_integ=1)
Definition: bilininteg.hpp:139
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:336
alpha (q . grad u, v)