MFEM  v4.1.0
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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_BILININTEG
13 #define MFEM_BILININTEG
14 
15 #include "../config/config.hpp"
16 #include "nonlininteg.hpp"
17 #include "fespace.hpp"
18 #include "libceed/ceed.hpp"
19 
20 namespace mfem
21 {
22 
23 /// Abstract base class BilinearFormIntegrator
25 {
26 protected:
28  : NonlinearFormIntegrator(ir) { }
29 
30 public:
31  // TODO: add support for other assembly levels (in addition to PA) and their
32  // actions.
33 
34  // TODO: for mixed meshes the quadrature rules to be used by methods like
35  // AssemblePA() can be given as a QuadratureSpace, e.g. using a new method:
36  // SetQuadratureSpace().
37 
38  // TODO: the methods for the various assembly levels make sense even in the
39  // base class NonlinearFormIntegrator, except that not all assembly levels
40  // make sense for the action of the nonlinear operator (but they all make
41  // sense for its Jacobian).
42 
44 
45  /// Method defining partial assembly.
46  /** The result of the partial assembly is stored internally so that it can be
47  used later in the methods AddMultPA() and AddMultTransposePA(). */
48  virtual void AssemblePA(const FiniteElementSpace &fes);
49  /** Used with BilinearFormIntegrators that have different spaces. */
50  virtual void AssemblePA(const FiniteElementSpace &trial_fes,
51  const FiniteElementSpace &test_fes);
52 
53  virtual void AssemblePAInteriorFaces(const FiniteElementSpace &fes);
54 
55  virtual void AssemblePABoundaryFaces(const FiniteElementSpace &fes);
56 
57  /// Assemble diagonal and add it to Vector @a diag.
58  virtual void AssembleDiagonalPA(Vector &diag);
59 
60  /// Method for partially assembled action.
61  /** Perform the action of integrator on the input @a x and add the result to
62  the output @a y. Both @a x and @a y are E-vectors, i.e. they represent
63  the element-wise discontinuous version of the FE space.
64 
65  This method can be called only after the method AssemblePA() has been
66  called. */
67  virtual void AddMultPA(const Vector &x, Vector &y) const;
68 
69  /// Method for partially assembled transposed action.
70  /** Perform the transpose action of integrator on the input @a x and add the
71  result to the output @a y. Both @a x and @a y are E-vectors, i.e. they
72  represent the element-wise discontinuous version of the FE space.
73 
74  This method can be called only after the method AssemblePA() has been
75  called. */
76  virtual void AddMultTransposePA(const Vector &x, Vector &y) const;
77 
78  /// Given a particular Finite Element computes the element matrix elmat.
79  virtual void AssembleElementMatrix(const FiniteElement &el,
81  DenseMatrix &elmat);
82 
83  /** Compute the local matrix representation of a bilinear form
84  a(u,v) defined on different trial (given by u) and test
85  (given by v) spaces. The rows in the local matrix correspond
86  to the test dofs and the columns -- to the trial dofs. */
87  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
88  const FiniteElement &test_fe,
90  DenseMatrix &elmat);
91 
92  virtual void AssembleFaceMatrix(const FiniteElement &el1,
93  const FiniteElement &el2,
95  DenseMatrix &elmat);
96 
97  /** Abstract method used for assembling TraceFaceIntegrators in a
98  MixedBilinearForm. */
99  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
100  const FiniteElement &test_fe1,
101  const FiniteElement &test_fe2,
103  DenseMatrix &elmat);
104 
105  /// Perform the local action of the BilinearFormIntegrator
106  virtual void AssembleElementVector(const FiniteElement &el,
108  const Vector &elfun, Vector &elvect);
109 
110  virtual void AssembleElementGrad(const FiniteElement &el,
112  const Vector &elfun, DenseMatrix &elmat)
113  { AssembleElementMatrix(el, Tr, elmat); }
114 
115  virtual void AssembleFaceGrad(const FiniteElement &el1,
116  const FiniteElement &el2,
118  const Vector &elfun, DenseMatrix &elmat)
119  { AssembleFaceMatrix(el1, el2, Tr, elmat); }
120 
121  /** @brief Virtual method required for Zienkiewicz-Zhu type error estimators.
122 
123  The purpose of the method is to compute a local "flux" finite element
124  function given a local finite element solution. The "flux" function has
125  to be computed in terms of its coefficients (represented by the Vector
126  @a flux) which multiply the basis functions defined by the FiniteElement
127  @a fluxelem. Typically, the "flux" function will have more than one
128  component and consequently @a flux should be store the coefficients of
129  all components: first all coefficient for component 0, then all
130  coefficients for component 1, etc. What the "flux" function represents
131  depends on the specific integrator. For example, in the case of
132  DiffusionIntegrator, the flux is the gradient of the solution multiplied
133  by the diffusion coefficient.
134 
135  @param[in] el FiniteElement of the solution.
136  @param[in] Trans The ElementTransformation describing the physical
137  position of the mesh element.
138  @param[in] u Solution coefficients representing the expansion of the
139  solution function in the basis of @a el.
140  @param[in] fluxelem FiniteElement of the "flux".
141  @param[out] flux "Flux" coefficients representing the expansion of the
142  "flux" function in the basis of @a fluxelem. The size
143  of @a flux as a Vector has to be set by this method,
144  e.g. using Vector::SetSize().
145  @param[in] with_coef If zero (the default value is 1) the implementation
146  of the method may choose not to scale the "flux"
147  function by any coefficients describing the
148  integrator.
149  */
150  virtual void ComputeElementFlux(const FiniteElement &el,
152  Vector &u,
153  const FiniteElement &fluxelem,
154  Vector &flux, bool with_coef = true) { }
155 
156  /** @brief Virtual method required for Zienkiewicz-Zhu type error estimators.
157 
158  The purpose of this method is to compute a local number that measures the
159  energy of a given "flux" function (see ComputeElementFlux() for a
160  description of the "flux" function). Typically, the energy of a "flux"
161  function should be equal to a_local(u,u), if the "flux" is defined from
162  a solution u; here a_local(.,.) denotes the element-local bilinear
163  form represented by the integrator.
164 
165  @param[in] fluxelem FiniteElement of the "flux".
166  @param[in] Trans The ElementTransformation describing the physical
167  position of the mesh element.
168  @param[in] flux "Flux" coefficients representing the expansion of the
169  "flux" function in the basis of @a fluxelem.
170  @param[out] d_energy If not NULL, the given Vector should be set to
171  represent directional energy split that can be used
172  for anisotropic error estimation.
173  @returns The computed energy.
174  */
175  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
177  Vector &flux, Vector *d_energy = NULL)
178  { return 0.0; }
179 
181 };
182 
184 {
185 private:
186  int own_bfi;
188 
189  DenseMatrix bfi_elmat;
190 
191 public:
192  TransposeIntegrator (BilinearFormIntegrator *_bfi, int _own_bfi = 1)
193  { bfi = _bfi; own_bfi = _own_bfi; }
194 
195  virtual void AssembleElementMatrix(const FiniteElement &el,
197  DenseMatrix &elmat);
198 
199  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
200  const FiniteElement &test_fe,
202  DenseMatrix &elmat);
203 
205  virtual void AssembleFaceMatrix(const FiniteElement &el1,
206  const FiniteElement &el2,
208  DenseMatrix &elmat);
209 
211 
212  virtual void AssemblePA(const FiniteElementSpace& fes)
213  {
214  bfi->AssemblePA(fes);
215  }
216 
218  {
219  bfi->AssemblePAInteriorFaces(fes);
220  }
221 
223  {
224  bfi->AssemblePABoundaryFaces(fes);
225  }
226 
227  virtual void AddMultTransposePA(const Vector &x, Vector &y) const
228  {
229  bfi->AddMultPA(x, y);
230  }
231 
232  virtual void AddMultPA(const Vector& x, Vector& y) const
233  {
234  bfi->AddMultTransposePA(x, y);
235  }
236 
237  virtual ~TransposeIntegrator() { if (own_bfi) { delete bfi; } }
238 };
239 
241 {
242 private:
243  int own_bfi;
245 
246 public:
247  LumpedIntegrator (BilinearFormIntegrator *_bfi, int _own_bfi = 1)
248  { bfi = _bfi; own_bfi = _own_bfi; }
249 
250  virtual void AssembleElementMatrix(const FiniteElement &el,
252  DenseMatrix &elmat);
253 
254  virtual ~LumpedIntegrator() { if (own_bfi) { delete bfi; } }
255 };
256 
257 /// Integrator that inverts the matrix assembled by another integrator.
259 {
260 private:
261  int own_integrator;
262  BilinearFormIntegrator *integrator;
263 
264 public:
265  InverseIntegrator(BilinearFormIntegrator *integ, int own_integ = 1)
266  { integrator = integ; own_integrator = own_integ; }
267 
268  virtual void AssembleElementMatrix(const FiniteElement &el,
270  DenseMatrix &elmat);
271 
272  virtual ~InverseIntegrator() { if (own_integrator) { delete integrator; } }
273 };
274 
275 /// Integrator defining a sum of multiple Integrators.
277 {
278 private:
279  int own_integrators;
280  DenseMatrix elem_mat;
281  Array<BilinearFormIntegrator*> integrators;
282 
283 public:
284  SumIntegrator(int own_integs = 1) { own_integrators = own_integs; }
285 
287  { integrators.Append(integ); }
288 
289  virtual void AssembleElementMatrix(const FiniteElement &el,
291  DenseMatrix &elmat);
292 
293  virtual ~SumIntegrator();
294 };
295 
296 /** An abstract class for integrating the product of two scalar basis functions
297  with an optional scalar coefficient. */
299 {
300 public:
301 
302  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
303  const FiniteElement &test_fe,
305  DenseMatrix &elmat);
306 
307  /// Support for use in BilinearForm. Can be used only when appropriate.
308  virtual void AssembleElementMatrix(const FiniteElement &fe,
309  ElementTransformation &Trans,
310  DenseMatrix &elmat)
311  { AssembleElementMatrix2(fe, fe, Trans, elmat); }
312 
313 protected:
314  /// This parameter can be set by derived methods to enable single shape
315  /// evaluation in case CalcTestShape() and CalcTrialShape() return the same
316  /// result if given the same FiniteElement. The default is false.
318 
321 
322  inline virtual bool VerifyFiniteElementTypes(
323  const FiniteElement & trial_fe,
324  const FiniteElement & test_fe) const
325  {
326  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
328  }
329 
330  inline virtual const char * FiniteElementTypeFailureMessage() const
331  {
332  return "MixedScalarIntegrator: "
333  "Trial and test spaces must both be scalar fields.";
334  }
335 
336  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
337  const FiniteElement & test_fe,
338  ElementTransformation &Trans)
339  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
340 
341 
342  inline virtual void CalcTestShape(const FiniteElement & test_fe,
343  ElementTransformation &Trans,
344  Vector & shape)
345  { test_fe.CalcPhysShape(Trans, shape); }
346 
347  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
348  ElementTransformation &Trans,
349  Vector & shape)
350  { trial_fe.CalcPhysShape(Trans, shape); }
351 
353 
354 private:
355 
356 #ifndef MFEM_THREAD_SAFE
357  Vector test_shape;
358  Vector trial_shape;
359 #endif
360 
361 };
362 
363 /** An abstract class for integrating the inner product of two vector basis
364  functions with an optional scalar, vector, or matrix coefficient. */
366 {
367 public:
368 
369  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
370  const FiniteElement &test_fe,
372  DenseMatrix &elmat);
373 
374  /// Support for use in BilinearForm. Can be used only when appropriate.
375  virtual void AssembleElementMatrix(const FiniteElement &fe,
376  ElementTransformation &Trans,
377  DenseMatrix &elmat)
378  { AssembleElementMatrix2(fe, fe, Trans, elmat); }
379 
380 protected:
381  /// This parameter can be set by derived methods to enable single shape
382  /// evaluation in case CalcTestShape() and CalcTrialShape() return the same
383  /// result if given the same FiniteElement. The default is false.
385 
387  : same_calc_shape(false), Q(NULL), VQ(NULL), DQ(NULL), MQ(NULL) {}
389  : same_calc_shape(false), Q(&q), VQ(NULL), DQ(NULL), MQ(NULL) {}
391  : same_calc_shape(false), Q(NULL), VQ(diag?NULL:&dq), DQ(diag?&dq:NULL),
392  MQ(NULL) {}
394  : same_calc_shape(false), Q(NULL), VQ(NULL), DQ(NULL), MQ(&mq) {}
395 
396  inline virtual bool VerifyFiniteElementTypes(
397  const FiniteElement & trial_fe,
398  const FiniteElement & test_fe) const
399  {
400  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
402  }
403 
404  inline virtual const char * FiniteElementTypeFailureMessage() const
405  {
406  return "MixedVectorIntegrator: "
407  "Trial and test spaces must both be vector fields";
408  }
409 
410  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
411  const FiniteElement & test_fe,
412  ElementTransformation &Trans)
413  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
414 
415 
416  inline virtual void CalcTestShape(const FiniteElement & test_fe,
417  ElementTransformation &Trans,
418  DenseMatrix & shape)
419  { test_fe.CalcVShape(Trans, shape); }
420 
421  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
422  ElementTransformation &Trans,
423  DenseMatrix & shape)
424  { trial_fe.CalcVShape(Trans, shape); }
425 
430 
431 private:
432 
433 #ifndef MFEM_THREAD_SAFE
434  Vector V;
435  Vector D;
436  DenseMatrix M;
437  DenseMatrix test_shape;
438  DenseMatrix trial_shape;
439  DenseMatrix test_shape_tmp;
440 #endif
441 
442 };
443 
444 /** An abstract class for integrating the product of a scalar basis function and
445  the inner product of a vector basis function with a vector coefficient. In
446  2D the inner product can be replaced with a cross product. */
448 {
449 public:
450 
451  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
452  const FiniteElement &test_fe,
454  DenseMatrix &elmat);
455 
456 protected:
457 
458  MixedScalarVectorIntegrator(VectorCoefficient &vq, bool _transpose = false,
459  bool _cross_2d = false)
460  : VQ(&vq), transpose(_transpose), cross_2d(_cross_2d) {}
461 
462  inline virtual bool VerifyFiniteElementTypes(
463  const FiniteElement & trial_fe,
464  const FiniteElement & test_fe) const
465  {
466  return ((transpose &&
468  test_fe.GetRangeType() == mfem::FiniteElement::SCALAR ) ||
469  (!transpose &&
472  );
473  }
474 
475  inline virtual const char * FiniteElementTypeFailureMessage() const
476  {
477  if ( transpose )
478  {
479  return "MixedScalarVectorIntegrator: "
480  "Trial space must be a vector field "
481  "and the test space must be a scalar field";
482  }
483  else
484  {
485  return "MixedScalarVectorIntegrator: "
486  "Trial space must be a scalar field "
487  "and the test space must be a vector field";
488  }
489  }
490 
491  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
492  const FiniteElement & test_fe,
493  ElementTransformation &Trans)
494  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
495 
496 
497  inline virtual void CalcVShape(const FiniteElement & vector_fe,
498  ElementTransformation &Trans,
499  DenseMatrix & shape)
500  { vector_fe.CalcVShape(Trans, shape); }
501 
502  inline virtual void CalcShape(const FiniteElement & scalar_fe,
503  ElementTransformation &Trans,
504  Vector & shape)
505  { scalar_fe.CalcPhysShape(Trans, shape); }
506 
508  bool transpose;
509  bool cross_2d; // In 2D use a cross product rather than a dot product
510 
511 private:
512 
513 #ifndef MFEM_THREAD_SAFE
514  Vector V;
515  DenseMatrix vshape;
516  Vector shape;
517  Vector vshape_tmp;
518 #endif
519 
520 };
521 
522 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 1D, 2D,
523  or 3D and where Q is an optional scalar coefficient, u and v are each in H1
524  or L2. */
526 {
527 public:
530  : MixedScalarIntegrator(q) { same_calc_shape = true; }
531 };
532 
533 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 2D, or
534  3D and where Q is a vector coefficient, u is in H1 or L2 and v is in H(Curl)
535  or H(Div). */
537 {
538 public:
541 };
542 
543 /** Class for integrating the bilinear form a(u,v) := (Q D u, v) in 1D where Q
544  is an optional scalar coefficient, u is in H1, and v is in L2. */
546 {
547 public:
550  : MixedScalarIntegrator(q) {}
551 
552 protected:
553  inline virtual bool VerifyFiniteElementTypes(
554  const FiniteElement & trial_fe,
555  const FiniteElement & test_fe) const
556  {
557  return (trial_fe.GetDim() == 1 && test_fe.GetDim() == 1 &&
558  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
560  }
561 
562  inline virtual const char * FiniteElementTypeFailureMessage() const
563  {
564  return "MixedScalarDerivativeIntegrator: "
565  "Trial and test spaces must both be scalar fields in 1D "
566  "and the trial space must implement CaldDShape.";
567  }
568 
569  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
571  Vector & shape)
572  {
573  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
574  trial_fe.CalcPhysDShape(Trans, dshape);
575  }
576 };
577 
578 /** Class for integrating the bilinear form a(u,v) := -(Q u, D v) in 1D where Q
579  is an optional scalar coefficient, u is in L2, and v is in H1. */
581 {
582 public:
585  : MixedScalarIntegrator(q) {}
586 
587 protected:
588  inline virtual bool VerifyFiniteElementTypes(
589  const FiniteElement & trial_fe,
590  const FiniteElement & test_fe) const
591  {
592  return (trial_fe.GetDim() == 1 && test_fe.GetDim() == 1 &&
595  }
596 
597  inline virtual const char * FiniteElementTypeFailureMessage() const
598  {
599  return "MixedScalarWeakDerivativeIntegrator: "
600  "Trial and test spaces must both be scalar fields in 1D "
601  "and the test space must implement CalcDShape with "
602  "map type \"VALUE\".";
603  }
604 
605  inline virtual void CalcTestShape(const FiniteElement & test_fe,
607  Vector & shape)
608  {
609  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
610  test_fe.CalcPhysDShape(Trans, dshape);
611  shape *= -1.0;
612  }
613 };
614 
615 /** Class for integrating the bilinear form a(u,v) := (Q div u, v) in either 2D
616  or 3D where Q is an optional scalar coefficient, u is in H(Div), and v is a
617  scalar field. */
619 {
620 public:
623  : MixedScalarIntegrator(q) {}
624 
625 protected:
626  inline virtual bool VerifyFiniteElementTypes(
627  const FiniteElement & trial_fe,
628  const FiniteElement & test_fe) const
629  {
630  return (trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
632  }
633 
634  inline virtual const char * FiniteElementTypeFailureMessage() const
635  {
636  return "MixedScalarDivergenceIntegrator: "
637  "Trial must be H(Div) and the test space must be a "
638  "scalar field";
639  }
640 
641  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
642  const FiniteElement & test_fe,
644  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
645 
646  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
648  Vector & shape)
649  { trial_fe.CalcPhysDivShape(Trans, shape); }
650 };
651 
652 /** Class for integrating the bilinear form a(u,v) := (V div u, v) in either 2D
653  or 3D where V is a vector coefficient, u is in H(Div), and v is a vector
654  field. */
656 {
657 public:
660 
661 protected:
662  inline virtual bool VerifyFiniteElementTypes(
663  const FiniteElement & trial_fe,
664  const FiniteElement & test_fe) const
665  {
666  return (trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
668  }
669 
670  inline virtual const char * FiniteElementTypeFailureMessage() const
671  {
672  return "MixedVectorDivergenceIntegrator: "
673  "Trial must be H(Div) and the test space must be a "
674  "vector field";
675  }
676 
677  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
678  const FiniteElement & test_fe,
680  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
681 
682  inline virtual void CalcShape(const FiniteElement & scalar_fe,
684  Vector & shape)
685  { scalar_fe.CalcPhysDivShape(Trans, shape); }
686 };
687 
688 /** Class for integrating the bilinear form a(u,v) := -(Q u, div v) in either 2D
689  or 3D where Q is an optional scalar coefficient, u is in L2 or H1, and v is
690  in H(Div). */
692 {
693 public:
696  : MixedScalarIntegrator(q) {}
697 
698 protected:
699  inline virtual bool VerifyFiniteElementTypes(
700  const FiniteElement & trial_fe,
701  const FiniteElement & test_fe) const
702  {
703  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
704  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
705  }
706 
707  inline virtual const char * FiniteElementTypeFailureMessage() const
708  {
709  return "MixedScalarWeakGradientIntegrator: "
710  "Trial space must be a scalar field "
711  "and the test space must be H(Div)";
712  }
713 
714  virtual void CalcTestShape(const FiniteElement & test_fe,
716  Vector & shape)
717  {
718  test_fe.CalcPhysDivShape(Trans, shape);
719  shape *= -1.0;
720  }
721 };
722 
723 /** Class for integrating the bilinear form a(u,v) := (Q curl u, v) in 2D where
724  Q is an optional scalar coefficient, u is in H(Curl), and v is in L2 or
725  H1. */
727 {
728 public:
731  : MixedScalarIntegrator(q) {}
732 
733 protected:
734  inline virtual bool VerifyFiniteElementTypes(
735  const FiniteElement & trial_fe,
736  const FiniteElement & test_fe) const
737  {
738  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
739  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
741  }
742 
743  inline virtual const char * FiniteElementTypeFailureMessage() const
744  {
745  return "MixedScalarCurlIntegrator: "
746  "Trial must be H(Curl) and the test space must be a "
747  "scalar field";
748  }
749 
750  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
751  const FiniteElement & test_fe,
753  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
754 
755  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
757  Vector & shape)
758  {
759  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
760  trial_fe.CalcPhysCurlShape(Trans, dshape);
761  }
762 };
763 
764 /** Class for integrating the bilinear form a(u,v) := (Q u, curl v) in 2D where
765  Q is an optional scalar coefficient, u is in L2 or H1, and v is in
766  H(Curl). */
768 {
769 public:
772  : MixedScalarIntegrator(q) {}
773 
774 protected:
775  inline virtual bool VerifyFiniteElementTypes(
776  const FiniteElement & trial_fe,
777  const FiniteElement & test_fe) const
778  {
779  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
782  }
783 
784  inline virtual const char * FiniteElementTypeFailureMessage() const
785  {
786  return "MixedScalarWeakCurlIntegrator: "
787  "Trial space must be a scalar field "
788  "and the test space must be H(Curl)";
789  }
790 
791  inline virtual void CalcTestShape(const FiniteElement & test_fe,
793  Vector & shape)
794  {
795  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
796  test_fe.CalcPhysCurlShape(Trans, dshape);
797  }
798 };
799 
800 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 2D or
801  3D and where Q is an optional coefficient (of type scalar, matrix, or
802  diagonal matrix) u and v are each in H(Curl) or H(Div). */
804 {
805 public:
808  : MixedVectorIntegrator(q) { same_calc_shape = true; }
810  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
812  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
813 };
814 
815 /** Class for integrating the bilinear form a(u,v) := (V x u, v) in 3D and where
816  V is a vector coefficient u and v are each in H(Curl) or H(Div). */
818 {
819 public:
821  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
822 };
823 
824 /** Class for integrating the bilinear form a(u,v) := (V . u, v) in 2D or 3D and
825  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in H1 or
826  L2. */
828 {
829 public:
831  : MixedScalarVectorIntegrator(vq, true) {}
832 
833  inline virtual bool VerifyFiniteElementTypes(
834  const FiniteElement & trial_fe,
835  const FiniteElement & test_fe) const
836  {
837  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
839  }
840 
841  inline virtual const char * FiniteElementTypeFailureMessage() const
842  {
843  return "MixedDotProductIntegrator: "
844  "Trial space must be a vector field "
845  "and the test space must be a scalar field";
846  }
847 };
848 
849 /** Class for integrating the bilinear form a(u,v) := (-V . u, Div v) in 2D or
850  3D and where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
851  RT. */
853 {
854 public:
856  : MixedScalarVectorIntegrator(vq, true) {}
857 
858  inline virtual bool VerifyFiniteElementTypes(
859  const FiniteElement & trial_fe,
860  const FiniteElement & test_fe) const
861  {
862  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
864  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
865  }
866 
867  inline virtual const char * FiniteElementTypeFailureMessage() const
868  {
869  return "MixedWeakGradDotIntegrator: "
870  "Trial space must be a vector field "
871  "and the test space must be a vector field with a divergence";
872  }
873 
874  inline virtual void CalcShape(const FiniteElement & scalar_fe,
876  Vector & shape)
877  { scalar_fe.CalcPhysDivShape(Trans, shape); shape *= -1.0; }
878 };
879 
880 /** Class for integrating the bilinear form a(u,v) := (V x u, Grad v) in 3D and
881  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in H1. */
883 {
884 public:
886  : MixedVectorIntegrator(vq, false) {}
887 
888  inline virtual bool VerifyFiniteElementTypes(
889  const FiniteElement & trial_fe,
890  const FiniteElement & test_fe) const
891  {
892  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
896  }
897 
898  inline virtual const char * FiniteElementTypeFailureMessage() const
899  {
900  return "MixedWeakDivCrossIntegrator: "
901  "Trial space must be a vector field in 3D "
902  "and the test space must be a scalar field with a gradient";
903  }
904 
905  inline virtual void CalcTestShape(const FiniteElement & test_fe,
907  DenseMatrix & shape)
908  { test_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
909 };
910 
911 /** Class for integrating the bilinear form a(u,v) := (Q Grad u, Grad v) in 3D
912  or in 2D and where Q is a scalar or matrix coefficient u and v are both in
913  H1. */
915 {
916 public:
919  : MixedVectorIntegrator(q) { same_calc_shape = true; }
921  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
923  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
924 
925  inline virtual bool VerifyFiniteElementTypes(
926  const FiniteElement & trial_fe,
927  const FiniteElement & test_fe) const
928  {
929  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
930  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
933  }
934 
935  inline virtual const char * FiniteElementTypeFailureMessage() const
936  {
937  return "MixedGradGradIntegrator: "
938  "Trial and test spaces must both be scalar fields "
939  "with a gradient operator.";
940  }
941 
942  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
943  const FiniteElement & test_fe,
945  {
946  // Same as DiffusionIntegrator
947  return test_fe.Space() == FunctionSpace::Pk ?
948  trial_fe.GetOrder() + test_fe.GetOrder() - 2 :
949  trial_fe.GetOrder() + test_fe.GetOrder() + test_fe.GetDim() - 1;
950  }
951 
952  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
954  DenseMatrix & shape)
955  { trial_fe.CalcPhysDShape(Trans, shape); }
956 
957  inline virtual void CalcTestShape(const FiniteElement & test_fe,
959  DenseMatrix & shape)
960  { test_fe.CalcPhysDShape(Trans, shape); }
961 };
962 
963 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, Grad v) in 3D
964  or in 2D and where V is a vector coefficient u and v are both in H1. */
966 {
967 public:
969  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
970 
971  inline virtual bool VerifyFiniteElementTypes(
972  const FiniteElement & trial_fe,
973  const FiniteElement & test_fe) const
974  {
975  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
976  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
979  }
980 
981  inline virtual const char * FiniteElementTypeFailureMessage() const
982  {
983  return "MixedCrossGradGradIntegrator: "
984  "Trial and test spaces must both be scalar fields "
985  "with a gradient operator.";
986  }
987 
988  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
990  DenseMatrix & shape)
991  { trial_fe.CalcPhysDShape(Trans, shape); }
992 
993  inline virtual void CalcTestShape(const FiniteElement & test_fe,
995  DenseMatrix & shape)
996  { test_fe.CalcPhysDShape(Trans, shape); }
997 };
998 
999 /** Class for integrating the bilinear form a(u,v) := (Q Curl u, Curl v) in 3D
1000  and where Q is a scalar or matrix coefficient u and v are both in
1001  H(Curl). */
1003 {
1004 public:
1007  : MixedVectorIntegrator(q) { same_calc_shape = true; }
1009  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
1011  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
1012 
1013  inline virtual bool VerifyFiniteElementTypes(
1014  const FiniteElement & trial_fe,
1015  const FiniteElement & test_fe) const
1016  {
1017  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1018  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1019  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1021  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1022  }
1023 
1024  inline virtual const char * FiniteElementTypeFailureMessage() const
1025  {
1026  return "MixedCurlCurlIntegrator"
1027  "Trial and test spaces must both be vector fields in 3D "
1028  "with a curl.";
1029  }
1030 
1031  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1033  DenseMatrix & shape)
1034  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1035 
1036  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1038  DenseMatrix & shape)
1039  { test_fe.CalcPhysCurlShape(Trans, shape); }
1040 };
1041 
1042 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, Curl v) in 3D
1043  and where V is a vector coefficient u and v are both in H(Curl). */
1045 {
1046 public:
1048  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
1049 
1050  inline virtual bool VerifyFiniteElementTypes(
1051  const FiniteElement & trial_fe,
1052  const FiniteElement & test_fe) const
1053  {
1054  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1055  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1056  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1058  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1059  }
1060 
1061  inline virtual const char * FiniteElementTypeFailureMessage() const
1062  {
1063  return "MixedCrossCurlCurlIntegrator: "
1064  "Trial and test spaces must both be vector fields in 3D "
1065  "with a curl.";
1066  }
1067 
1068  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1070  DenseMatrix & shape)
1071  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1072 
1073  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1075  DenseMatrix & shape)
1076  { test_fe.CalcPhysCurlShape(Trans, shape); }
1077 };
1078 
1079 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, Grad v) in 3D
1080  and where V is a vector coefficient u is in H(Curl) and v is in H1. */
1082 {
1083 public:
1085  : MixedVectorIntegrator(vq, false) {}
1086 
1087  inline virtual bool VerifyFiniteElementTypes(
1088  const FiniteElement & trial_fe,
1089  const FiniteElement & test_fe) const
1090  {
1091  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1092  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1093  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1095  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1096  }
1097 
1098  inline virtual const char * FiniteElementTypeFailureMessage() const
1099  {
1100  return "MixedCrossCurlGradIntegrator"
1101  "Trial space must be a vector field in 3D with a curl"
1102  "and the test space must be a scalar field with a gradient";
1103  }
1104 
1105  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1107  DenseMatrix & shape)
1108  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1109 
1110  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1112  DenseMatrix & shape)
1113  { test_fe.CalcPhysDShape(Trans, shape); }
1114 };
1115 
1116 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, Curl v) in 3D
1117  and where V is a scalar coefficient u is in H1 and v is in H(Curl). */
1119 {
1120 public:
1122  : MixedVectorIntegrator(vq, false) {}
1123 
1124  inline virtual bool VerifyFiniteElementTypes(
1125  const FiniteElement & trial_fe,
1126  const FiniteElement & test_fe) const
1127  {
1128  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1129  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1130  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1132  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1133  }
1134 
1135  inline virtual const char * FiniteElementTypeFailureMessage() const
1136  {
1137  return "MixedCrossGradCurlIntegrator"
1138  "Trial space must be a scalar field in 3D with a gradient"
1139  "and the test space must be a vector field with a curl";
1140  }
1141 
1142  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1144  DenseMatrix & shape)
1145  { trial_fe.CalcPhysDShape(Trans, shape); }
1146 
1147  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1149  DenseMatrix & shape)
1150  { test_fe.CalcPhysCurlShape(Trans, shape); }
1151 };
1152 
1153 /** Class for integrating the bilinear form a(u,v) := (V x u, Curl v) in 3D and
1154  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
1155  H(Curl). */
1157 {
1158 public:
1160  : MixedVectorIntegrator(vq, false) {}
1161 
1162  inline virtual bool VerifyFiniteElementTypes(
1163  const FiniteElement & trial_fe,
1164  const FiniteElement & test_fe) const
1165  {
1166  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1167  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1169  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1170  }
1171 
1172  inline virtual const char * FiniteElementTypeFailureMessage() const
1173  {
1174  return "MixedWeakCurlCrossIntegrator: "
1175  "Trial space must be a vector field in 3D "
1176  "and the test space must be a vector field with a curl";
1177  }
1178 
1179  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1181  DenseMatrix & shape)
1182  { test_fe.CalcPhysCurlShape(Trans, shape); }
1183 };
1184 
1185 /** Class for integrating the bilinear form a(u,v) := (V x u, Curl v) in 2D and
1186  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
1187  H(Curl). */
1189 {
1190 public:
1192  : MixedScalarVectorIntegrator(vq, true, true) {}
1193 
1194  inline virtual bool VerifyFiniteElementTypes(
1195  const FiniteElement & trial_fe,
1196  const FiniteElement & test_fe) const
1197  {
1198  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1199  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1201  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1202  }
1203 
1204  inline virtual const char * FiniteElementTypeFailureMessage() const
1205  {
1206  return "MixedScalarWeakCurlCrossIntegrator: "
1207  "Trial space must be a vector field in 2D "
1208  "and the test space must be a vector field with a curl";
1209  }
1210 
1211  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1213  Vector & shape)
1214  {
1215  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
1216  scalar_fe.CalcPhysCurlShape(Trans, dshape);
1217  }
1218 };
1219 
1220 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, v) in 3D or
1221  in 2D and where V is a vector coefficient u is in H1 and v is in H(Curl) or
1222  H(Div). */
1224 {
1225 public:
1227  : MixedVectorIntegrator(vq, false) {}
1228 
1229  inline virtual bool VerifyFiniteElementTypes(
1230  const FiniteElement & trial_fe,
1231  const FiniteElement & test_fe) const
1232  {
1233  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1234  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1235  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1237  }
1238 
1239  inline virtual const char * FiniteElementTypeFailureMessage() const
1240  {
1241  return "MixedCrossGradIntegrator: "
1242  "Trial space must be a scalar field with a gradient operator"
1243  " and the test space must be a vector field both in 3D.";
1244  }
1245 
1246  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1248  DenseMatrix & shape)
1249  { trial_fe.CalcPhysDShape(Trans, shape); }
1250 
1251  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1253  DenseMatrix & shape)
1254  { test_fe.CalcVShape(Trans, shape); }
1255 };
1256 
1257 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, v) in 3D and
1258  where V is a vector coefficient u is in H(Curl) and v is in H(Curl) or
1259  H(Div). */
1261 {
1262 public:
1264  : MixedVectorIntegrator(vq, false) {}
1265 
1266  inline virtual bool VerifyFiniteElementTypes(
1267  const FiniteElement & trial_fe,
1268  const FiniteElement & test_fe) const
1269  {
1270  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1271  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1272  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1274  }
1275 
1276  inline virtual const char * FiniteElementTypeFailureMessage() const
1277  {
1278  return "MixedCrossCurlIntegrator: "
1279  "Trial space must be a vector field in 3D with a curl "
1280  "and the test space must be a vector field";
1281  }
1282 
1283  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1285  DenseMatrix & shape)
1286  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1287 };
1288 
1289 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, v) in 2D and
1290  where V is a vector coefficient u is in H(Curl) and v is in H(Curl) or
1291  H(Div). */
1293 {
1294 public:
1296  : MixedScalarVectorIntegrator(vq, false, true) {}
1297 
1298  inline virtual bool VerifyFiniteElementTypes(
1299  const FiniteElement & trial_fe,
1300  const FiniteElement & test_fe) const
1301  {
1302  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1303  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1304  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1306  }
1307 
1308  inline virtual const char * FiniteElementTypeFailureMessage() const
1309  {
1310  return "MixedCrossCurlIntegrator: "
1311  "Trial space must be a vector field in 2D with a curl "
1312  "and the test space must be a vector field";
1313  }
1314 
1315  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1317  Vector & shape)
1318  {
1319  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
1320  scalar_fe.CalcPhysCurlShape(Trans, dshape); shape *= -1.0;
1321  }
1322 };
1323 
1324 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, v) in 2D and
1325  where V is a vector coefficient u is in H1 and v is in H1 or L2. */
1327 {
1328 public:
1330  : MixedScalarVectorIntegrator(vq, true, true) {}
1331 
1332  inline virtual bool VerifyFiniteElementTypes(
1333  const FiniteElement & trial_fe,
1334  const FiniteElement & test_fe) const
1335  {
1336  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1337  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1338  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1340  }
1341 
1342  inline virtual const char * FiniteElementTypeFailureMessage() const
1343  {
1344  return "MixedScalarCrossGradIntegrator: "
1345  "Trial space must be a scalar field in 2D with a gradient "
1346  "and the test space must be a scalar field";
1347  }
1348 
1349  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1351  DenseMatrix & shape)
1352  { vector_fe.CalcPhysDShape(Trans, shape); }
1353 };
1354 
1355 /** Class for integrating the bilinear form a(u,v) := (V x u, v) in 2D and where
1356  V is a vector coefficient u is in ND or RT and v is in H1 or L2. */
1358 {
1359 public:
1361  : MixedScalarVectorIntegrator(vq, true, true) {}
1362 
1363  inline virtual bool VerifyFiniteElementTypes(
1364  const FiniteElement & trial_fe,
1365  const FiniteElement & test_fe) const
1366  {
1367  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1368  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1370  }
1371 
1372  inline virtual const char * FiniteElementTypeFailureMessage() const
1373  {
1374  return "MixedScalarCrossProductIntegrator: "
1375  "Trial space must be a vector field in 2D "
1376  "and the test space must be a scalar field";
1377  }
1378 };
1379 
1380 /** Class for integrating the bilinear form a(u,v) := (V x z u, v) in 2D and
1381  where V is a vector coefficient u is in H1 or L2 and v is in ND or RT. */
1383 {
1384 public:
1386  : MixedScalarVectorIntegrator(vq, false, true) {}
1387 
1388  inline virtual bool VerifyFiniteElementTypes(
1389  const FiniteElement & trial_fe,
1390  const FiniteElement & test_fe) const
1391  {
1392  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1393  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1395  }
1396 
1397  inline virtual const char * FiniteElementTypeFailureMessage() const
1398  {
1399  return "MixedScalarWeakCrossProductIntegrator: "
1400  "Trial space must be a scalar field in 2D "
1401  "and the test space must be a vector field";
1402  }
1403 
1404  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1406  Vector & shape)
1407  { scalar_fe.CalcPhysShape(Trans, shape); shape *= -1.0; }
1408 };
1409 
1410 /** Class for integrating the bilinear form a(u,v) := (V . Grad u, v) in 2D or
1411  3D and where V is a vector coefficient, u is in H1 and v is in H1 or L2. */
1413 {
1414 public:
1416  : MixedScalarVectorIntegrator(vq, true) {}
1417 
1418  inline virtual bool VerifyFiniteElementTypes(
1419  const FiniteElement & trial_fe,
1420  const FiniteElement & test_fe) const
1421  {
1422  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1423  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1425  }
1426 
1427  inline virtual const char * FiniteElementTypeFailureMessage() const
1428  {
1429  return "MixedDirectionalDerivativeIntegrator: "
1430  "Trial space must be a scalar field with a gradient "
1431  "and the test space must be a scalar field";
1432  }
1433 
1434  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1436  DenseMatrix & shape)
1437  { vector_fe.CalcPhysDShape(Trans, shape); }
1438 };
1439 
1440 /** Class for integrating the bilinear form a(u,v) := (-V . Grad u, Div v) in 2D
1441  or 3D and where V is a vector coefficient, u is in H1 and v is in RT. */
1443 {
1444 public:
1446  : MixedScalarVectorIntegrator(vq, true) {}
1447 
1448  inline virtual bool VerifyFiniteElementTypes(
1449  const FiniteElement & trial_fe,
1450  const FiniteElement & test_fe) const
1451  {
1452  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1453  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1455  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
1456  }
1457 
1458  inline virtual const char * FiniteElementTypeFailureMessage() const
1459  {
1460  return "MixedGradDivIntegrator: "
1461  "Trial space must be a scalar field with a gradient"
1462  "and the test space must be a vector field with a divergence";
1463  }
1464 
1465  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1467  DenseMatrix & shape)
1468  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1469 
1470  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1472  Vector & shape)
1473  { scalar_fe.CalcPhysDivShape(Trans, shape); }
1474 };
1475 
1476 /** Class for integrating the bilinear form a(u,v) := (-V Div u, Grad v) in 2D
1477  or 3D and where V is a vector coefficient, u is in RT and v is in H1. */
1479 {
1480 public:
1482  : MixedScalarVectorIntegrator(vq, false) {}
1483 
1484  inline virtual bool VerifyFiniteElementTypes(
1485  const FiniteElement & trial_fe,
1486  const FiniteElement & test_fe) const
1487  {
1488  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1489  trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
1492  );
1493  }
1494 
1495  inline virtual const char * FiniteElementTypeFailureMessage() const
1496  {
1497  return "MixedDivGradIntegrator: "
1498  "Trial space must be a vector field with a divergence"
1499  "and the test space must be a scalar field with a gradient";
1500  }
1501 
1502  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1504  DenseMatrix & shape)
1505  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1506 
1507  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1509  Vector & shape)
1510  { scalar_fe.CalcPhysDivShape(Trans, shape); }
1511 };
1512 
1513 /** Class for integrating the bilinear form a(u,v) := (-V u, Grad v) in 2D or 3D
1514  and where V is a vector coefficient, u is in H1 and v is in H1. */
1516 {
1517 public:
1519  : MixedScalarVectorIntegrator(vq, false) {}
1520 
1521  inline virtual bool VerifyFiniteElementTypes(
1522  const FiniteElement & trial_fe,
1523  const FiniteElement & test_fe) const
1524  {
1525  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1527  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1528  }
1529 
1530  inline virtual const char * FiniteElementTypeFailureMessage() const
1531  {
1532  return "MixedScalarWeakDivergenceIntegrator: "
1533  "Trial space must be a scalar field "
1534  "and the test space must be a scalar field with a gradient";
1535  }
1536 
1537  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1539  DenseMatrix & shape)
1540  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1541 };
1542 
1543 /** Class for integrating the bilinear form a(u,v) := (Q grad u, v) in either 2D
1544  or 3D and where Q is an optional coefficient (of type scalar, matrix, or
1545  diagonal matrix) u is in H1 and v is in H(Curl) or H(Div). */
1547 {
1548 public:
1551  : MixedVectorIntegrator(q) {}
1553  : MixedVectorIntegrator(dq, true) {}
1555  : MixedVectorIntegrator(mq) {}
1556 
1557 protected:
1558  inline virtual bool VerifyFiniteElementTypes(
1559  const FiniteElement & trial_fe,
1560  const FiniteElement & test_fe) const
1561  {
1562  return (trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1564  }
1565 
1566  inline virtual const char * FiniteElementTypeFailureMessage() const
1567  {
1568  return "MixedVectorGradientIntegrator: "
1569  "Trial spaces must be H1 and the test space must be a "
1570  "vector field in 2D or 3D";
1571  }
1572 
1573  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1575  DenseMatrix & shape)
1576  {
1577  trial_fe.CalcPhysDShape(Trans, shape);
1578  }
1579 
1581  virtual void AssemblePA(const FiniteElementSpace &trial_fes,
1582  const FiniteElementSpace &test_fes);
1583 
1584  virtual void AddMultPA(const Vector&, Vector&) const;
1585 
1586 private:
1587  DenseMatrix Jinv;
1588 
1589  // PA extension
1590  Vector pa_data;
1591  const DofToQuad *mapsO; ///< Not owned. DOF-to-quad map, open.
1592  const DofToQuad *mapsC; ///< Not owned. DOF-to-quad map, closed.
1593  const GeometricFactors *geom; ///< Not owned
1594  int dim, ne, dofs1D, quad1D;
1595 };
1596 
1597 /** Class for integrating the bilinear form a(u,v) := (Q curl u, v) in 3D and
1598  where Q is an optional coefficient (of type scalar, matrix, or diagonal
1599  matrix) u is in H(Curl) and v is in H(Div) or H(Curl). */
1601 {
1602 public:
1605  : MixedVectorIntegrator(q) {}
1607  : MixedVectorIntegrator(dq, true) {}
1609  : MixedVectorIntegrator(mq) {}
1610 
1611 protected:
1612  inline virtual bool VerifyFiniteElementTypes(
1613  const FiniteElement & trial_fe,
1614  const FiniteElement & test_fe) const
1615  {
1616  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1617  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1619  }
1620 
1621  inline virtual const char * FiniteElementTypeFailureMessage() const
1622  {
1623  return "MixedVectorCurlIntegrator: "
1624  "Trial space must be H(Curl) and the test space must be a "
1625  "vector field in 3D";
1626  }
1627 
1628  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1630  DenseMatrix & shape)
1631  {
1632  trial_fe.CalcPhysCurlShape(Trans, shape);
1633  }
1634 };
1635 
1636 /** Class for integrating the bilinear form a(u,v) := (Q u, curl v) in 3D and
1637  where Q is an optional coefficient (of type scalar, matrix, or diagonal
1638  matrix) u is in H(Div) or H(Curl) and v is in H(Curl). */
1640 {
1641 public:
1644  : MixedVectorIntegrator(q) {}
1646  : MixedVectorIntegrator(dq, true) {}
1648  : MixedVectorIntegrator(mq) {}
1649 
1650 protected:
1651  inline virtual bool VerifyFiniteElementTypes(
1652  const FiniteElement & trial_fe,
1653  const FiniteElement & test_fe) const
1654  {
1655  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1656  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1657  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1658  }
1659 
1660  inline virtual const char * FiniteElementTypeFailureMessage() const
1661  {
1662  return "MixedVectorWeakCurlIntegrator: "
1663  "Trial space must be vector field in 3D and the "
1664  "test space must be H(Curl)";
1665  }
1666 
1667  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1669  DenseMatrix & shape)
1670  {
1671  test_fe.CalcPhysCurlShape(Trans, shape);
1672  }
1673 };
1674 
1675 /** Class for integrating the bilinear form a(u,v) := - (Q u, grad v) in either
1676  2D or 3D and where Q is an optional coefficient (of type scalar, matrix, or
1677  diagonal matrix) u is in H(Div) or H(Curl) and v is in H1. */
1679 {
1680 public:
1683  : MixedVectorIntegrator(q) {}
1685  : MixedVectorIntegrator(dq, true) {}
1687  : MixedVectorIntegrator(mq) {}
1688 
1689 protected:
1690  inline virtual bool VerifyFiniteElementTypes(
1691  const FiniteElement & trial_fe,
1692  const FiniteElement & test_fe) const
1693  {
1694  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1695  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1696  }
1697 
1698  inline virtual const char * FiniteElementTypeFailureMessage() const
1699  {
1700  return "MixedVectorWeakDivergenceIntegrator: "
1701  "Trial space must be vector field and the "
1702  "test space must be H1";
1703  }
1704 
1705  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1707  DenseMatrix & shape)
1708  {
1709  test_fe.CalcPhysDShape(Trans, shape);
1710  shape *= -1.0;
1711  }
1712 };
1713 
1714 /** Class for integrating the bilinear form a(u,v) := (Q grad u, v) where Q is a
1715  scalar coefficient, and v is a vector with components v_i in the same space
1716  as u. */
1718 {
1719 protected:
1721 
1722 private:
1723  Vector shape;
1724  DenseMatrix dshape;
1725  DenseMatrix gshape;
1726  DenseMatrix Jadj;
1727  DenseMatrix elmat_comp;
1728  // PA extension
1729  Vector pa_data;
1730  const DofToQuad *trial_maps, *test_maps; ///< Not owned
1731  const GeometricFactors *geom; ///< Not owned
1732  int dim, ne, nq;
1733  int trial_dofs1D, test_dofs1D, quad1D;
1734 
1735 public:
1737  Q{NULL}, trial_maps{NULL}, test_maps{NULL}, geom{NULL}
1738  { }
1740  Q{_q}, trial_maps{NULL}, test_maps{NULL}, geom{NULL}
1741  { }
1743  Q{&q}, trial_maps{NULL}, test_maps{NULL}, geom{NULL}
1744  { }
1745 
1746  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1747  const FiniteElement &test_fe,
1748  ElementTransformation &Trans,
1749  DenseMatrix &elmat);
1750 
1752  virtual void AssemblePA(const FiniteElementSpace &trial_fes,
1753  const FiniteElementSpace &test_fes);
1754 
1755  virtual void AddMultPA(const Vector &x, Vector &y) const;
1756  virtual void AddMultTransposePA(const Vector &x, Vector &y) const;
1757 
1758  static const IntegrationRule &GetRule(const FiniteElement &trial_fe,
1759  const FiniteElement &test_fe,
1760  ElementTransformation &Trans);
1761 };
1762 
1763 /** Class for integrating the bilinear form a(u,v) := (Q grad u, grad v) where Q
1764  can be a scalar or a matrix coefficient. */
1766 {
1767 protected:
1770 
1771 private:
1772  Vector vec, pointflux, shape;
1773 #ifndef MFEM_THREAD_SAFE
1774  DenseMatrix dshape, dshapedxt, invdfdx, mq;
1775  DenseMatrix te_dshape, te_dshapedxt;
1776 #endif
1777 
1778  // PA extension
1779  const FiniteElementSpace *fespace;
1780  const DofToQuad *maps; ///< Not owned
1781  const GeometricFactors *geom; ///< Not owned
1782  int dim, ne, dofs1D, quad1D;
1783  Vector pa_data;
1784 
1785 #ifdef MFEM_USE_CEED
1786  // CEED extension
1787  CeedData* ceedDataPtr;
1788 #endif
1789 
1790 public:
1791  /// Construct a diffusion integrator with coefficient Q = 1
1793  {
1794  Q = NULL;
1795  MQ = NULL;
1796  maps = NULL;
1797  geom = NULL;
1798 #ifdef MFEM_USE_CEED
1799  ceedDataPtr = NULL;
1800 #endif
1801  }
1802 
1803  /// Construct a diffusion integrator with a scalar coefficient q
1805  : Q(&q)
1806  {
1807  MQ = NULL;
1808  maps = NULL;
1809  geom = NULL;
1810 #ifdef MFEM_USE_CEED
1811  ceedDataPtr = NULL;
1812 #endif
1813  }
1814 
1815  /// Construct a diffusion integrator with a matrix coefficient q
1817  : MQ(&q)
1818  {
1819  Q = NULL;
1820  maps = NULL;
1821  geom = NULL;
1822 #ifdef MFEM_USE_CEED
1823  ceedDataPtr = NULL;
1824 #endif
1825  }
1826 
1828  {
1829 #ifdef MFEM_USE_CEED
1830  delete ceedDataPtr;
1831 #endif
1832  }
1833 
1834  /** Given a particular Finite Element computes the element stiffness matrix
1835  elmat. */
1836  virtual void AssembleElementMatrix(const FiniteElement &el,
1838  DenseMatrix &elmat);
1839  /** Given a trial and test Finite Element computes the element stiffness
1840  matrix elmat. */
1841  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1842  const FiniteElement &test_fe,
1844  DenseMatrix &elmat);
1845 
1846  /// Perform the local action of the BilinearFormIntegrator
1847  virtual void AssembleElementVector(const FiniteElement &el,
1849  const Vector &elfun, Vector &elvect);
1850 
1851  virtual void ComputeElementFlux(const FiniteElement &el,
1853  Vector &u, const FiniteElement &fluxelem,
1854  Vector &flux, bool with_coef = true);
1855 
1856  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
1858  Vector &flux, Vector *d_energy = NULL);
1859 
1861 
1862  virtual void AssemblePA(const FiniteElementSpace &fes);
1863 
1864  virtual void AssembleDiagonalPA(Vector &diag);
1865 
1866  virtual void AddMultPA(const Vector&, Vector&) const;
1867 
1868  static const IntegrationRule &GetRule(const FiniteElement &trial_fe,
1869  const FiniteElement &test_fe);
1870 
1871  void SetupPA(const FiniteElementSpace &fes, const bool force = false);
1872 };
1873 
1874 /** Class for local mass matrix assembling a(u,v) := (Q u, v) */
1876 {
1877 protected:
1878 #ifndef MFEM_THREAD_SAFE
1880 #endif
1882  // PA extension
1885  const DofToQuad *maps; ///< Not owned
1886  const GeometricFactors *geom; ///< Not owned
1887  int dim, ne, nq, dofs1D, quad1D;
1888 
1889 #ifdef MFEM_USE_CEED
1890  // CEED extension
1891  CeedData* ceedDataPtr;
1892 #endif
1893 
1894 public:
1897  {
1898  Q = NULL;
1899  maps = NULL;
1900  geom = NULL;
1901 #ifdef MFEM_USE_CEED
1902  ceedDataPtr = NULL;
1903 #endif
1904  }
1905 
1906  /// Construct a mass integrator with coefficient q
1908  : BilinearFormIntegrator(ir), Q(&q)
1909  {
1910  maps = NULL;
1911  geom = NULL;
1912 #ifdef MFEM_USE_CEED
1913  ceedDataPtr = NULL;
1914 #endif
1915  }
1916 
1918  {
1919 #ifdef MFEM_USE_CEED
1920  delete ceedDataPtr;
1921 #endif
1922  }
1923  /** Given a particular Finite Element computes the element mass matrix
1924  elmat. */
1925  virtual void AssembleElementMatrix(const FiniteElement &el,
1927  DenseMatrix &elmat);
1928  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1929  const FiniteElement &test_fe,
1931  DenseMatrix &elmat);
1932 
1934 
1935  virtual void AssemblePA(const FiniteElementSpace &fes);
1936 
1937  virtual void AssembleDiagonalPA(Vector &diag);
1938 
1939  virtual void AddMultPA(const Vector&, Vector&) const;
1940 
1941  static const IntegrationRule &GetRule(const FiniteElement &trial_fe,
1942  const FiniteElement &test_fe,
1944 
1945  void SetupPA(const FiniteElementSpace &fes, const bool force = false);
1946 };
1947 
1949 {
1950 public:
1952 
1954 
1955  virtual void AssembleFaceMatrix(const FiniteElement &el1,
1956  const FiniteElement &el2,
1958  DenseMatrix &elmat);
1959 };
1960 
1961 /// alpha (q . grad u, v)
1963 {
1964 protected:
1966  double alpha;
1967  // PA extension
1969  const DofToQuad *maps; ///< Not owned
1970  const GeometricFactors *geom; ///< Not owned
1971  int dim, ne, nq, dofs1D, quad1D;
1972 
1973 private:
1974 #ifndef MFEM_THREAD_SAFE
1975  DenseMatrix dshape, adjJ, Q_ir;
1976  Vector shape, vec2, BdFidxT;
1977 #endif
1978 
1979 public:
1981  : Q(&q) { alpha = a; }
1982  virtual void AssembleElementMatrix(const FiniteElement &,
1984  DenseMatrix &);
1985 
1987 
1988  virtual void AssemblePA(const FiniteElementSpace&);
1989 
1990  virtual void AddMultPA(const Vector&, Vector&) const;
1991 
1992  static const IntegrationRule &GetRule(const FiniteElement &el,
1994 
1995  static const IntegrationRule &GetRule(const FiniteElement &trial_fe,
1996  const FiniteElement &test_fe,
1998 };
1999 
2000 /// alpha (q . grad u, v) using the "group" FE discretization
2002 {
2003 protected:
2005  double alpha;
2006 
2007 private:
2008  DenseMatrix dshape, adjJ, Q_nodal, grad;
2009  Vector shape;
2010 
2011 public:
2013  : Q(&q) { alpha = a; }
2014  virtual void AssembleElementMatrix(const FiniteElement &,
2016  DenseMatrix &);
2017 };
2018 
2019 /** Class for integrating the bilinear form a(u,v) := (Q u, v),
2020  where u=(u1,...,un) and v=(v1,...,vn); ui and vi are defined
2021  by scalar FE through standard transformation. */
2023 {
2024 private:
2025  int vdim;
2026  Vector shape, te_shape, vec;
2027  DenseMatrix partelmat;
2028  DenseMatrix mcoeff;
2029  int Q_order;
2030 
2031 protected:
2035  // PA extension
2037  const DofToQuad *maps; ///< Not owned
2038  const GeometricFactors *geom; ///< Not owned
2039  int dim, ne, nq, dofs1D, quad1D;
2040 
2041 public:
2042  /// Construct an integrator with coefficient 1.0
2044  : vdim(-1), Q_order(0), Q(NULL), VQ(NULL), MQ(NULL) { }
2045  /** Construct an integrator with scalar coefficient q. If possible, save
2046  memory by using a scalar integrator since the resulting matrix is block
2047  diagonal with the same diagonal block repeated. */
2049  : vdim(-1), Q(&q) { VQ = NULL; MQ = NULL; Q_order = qo; }
2051  : BilinearFormIntegrator(ir), vdim(-1), Q(&q)
2052  { VQ = NULL; MQ = NULL; Q_order = 0; }
2053  /// Construct an integrator with diagonal coefficient q
2055  : vdim(q.GetVDim()), VQ(&q) { Q = NULL; MQ = NULL; Q_order = qo; }
2056  /// Construct an integrator with matrix coefficient q
2058  : vdim(q.GetVDim()), MQ(&q) { Q = NULL; VQ = NULL; Q_order = qo; }
2059 
2060  int GetVDim() const { return vdim; }
2061  void SetVDim(int vdim) { this->vdim = vdim; }
2062 
2063  virtual void AssembleElementMatrix(const FiniteElement &el,
2065  DenseMatrix &elmat);
2066  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2067  const FiniteElement &test_fe,
2069  DenseMatrix &elmat);
2071  virtual void AssemblePA(const FiniteElementSpace &fes);
2072  virtual void AssembleDiagonalPA(Vector &diag);
2073  virtual void AddMultPA(const Vector &x, Vector &y) const;
2074 };
2075 
2076 
2077 /** Class for integrating (div u, p) where u is a vector field given by
2078  VectorFiniteElement through Piola transformation (for RT elements); p is
2079  scalar function given by FiniteElement through standard transformation.
2080  Here, u is the trial function and p is the test function.
2081 
2082  Note: the element matrix returned by AssembleElementMatrix2 does NOT depend
2083  on the ElementTransformation Trans. */
2085 {
2086 protected:
2088 
2089 private:
2090 #ifndef MFEM_THREAD_SAFE
2091  Vector divshape, shape;
2092 #endif
2093 
2094 public:
2097  virtual void AssembleElementMatrix(const FiniteElement &el,
2099  DenseMatrix &elmat) { }
2100  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2101  const FiniteElement &test_fe,
2103  DenseMatrix &elmat);
2104 };
2105 
2106 
2107 /** Integrator for `(-Q u, grad v)` for Nedelec (`u`) and H1 (`v`) elements.
2108  This is equivalent to a weak divergence of the Nedelec basis functions. */
2110 {
2111 protected:
2113 
2114 private:
2115 #ifndef MFEM_THREAD_SAFE
2116  DenseMatrix dshape;
2117  DenseMatrix dshapedxt;
2118  DenseMatrix vshape;
2119  DenseMatrix invdfdx;
2120 #endif
2121 
2122 public:
2125  virtual void AssembleElementMatrix(const FiniteElement &el,
2127  DenseMatrix &elmat) { }
2128  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2129  const FiniteElement &test_fe,
2131  DenseMatrix &elmat);
2132 };
2133 
2134 /** Integrator for (curl u, v) for Nedelec and RT elements. If the trial and
2135  test spaces are switched, assembles the form (u, curl v). */
2137 {
2138 protected:
2140 
2141 private:
2142 #ifndef MFEM_THREAD_SAFE
2143  DenseMatrix curlshapeTrial;
2144  DenseMatrix vshapeTest;
2145  DenseMatrix curlshapeTrial_dFT;
2146 #endif
2147 
2148 public:
2151  virtual void AssembleElementMatrix(const FiniteElement &el,
2153  DenseMatrix &elmat) { }
2154  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2155  const FiniteElement &test_fe,
2157  DenseMatrix &elmat);
2158 };
2159 
2160 /// Class for integrating (Q D_i(u), v); u and v are scalars
2162 {
2163 protected:
2165 
2166 private:
2167  int xi;
2168  DenseMatrix dshape, dshapedxt, invdfdx;
2169  Vector shape, dshapedxi;
2170 
2171 public:
2172  DerivativeIntegrator(Coefficient &q, int i) : Q(&q), xi(i) { }
2173  virtual void AssembleElementMatrix(const FiniteElement &el,
2175  DenseMatrix &elmat)
2176  { AssembleElementMatrix2(el,el,Trans,elmat); }
2177  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2178  const FiniteElement &test_fe,
2180  DenseMatrix &elmat);
2181 };
2182 
2183 /// Integrator for (curl u, curl v) for Nedelec elements
2185 {
2186 private:
2187  Vector vec, pointflux;
2188 #ifndef MFEM_THREAD_SAFE
2189  DenseMatrix curlshape, curlshape_dFt, M;
2190  DenseMatrix vshape, projcurl;
2191 #endif
2192 
2193 protected:
2196 
2197  // PA extension
2199  const DofToQuad *mapsO; ///< Not owned. DOF-to-quad map, open.
2200  const DofToQuad *mapsC; ///< Not owned. DOF-to-quad map, closed.
2201  const GeometricFactors *geom; ///< Not owned
2202  int dim, ne, nq, dofs1D, quad1D;
2203 
2204 public:
2205  CurlCurlIntegrator() { Q = NULL; MQ = NULL; }
2206  /// Construct a bilinear form integrator for Nedelec elements
2207  CurlCurlIntegrator(Coefficient &q) : Q(&q) { MQ = NULL; }
2209 
2210  /* Given a particular Finite Element, compute the
2211  element curl-curl matrix elmat */
2212  virtual void AssembleElementMatrix(const FiniteElement &el,
2214  DenseMatrix &elmat);
2215 
2216  virtual void ComputeElementFlux(const FiniteElement &el,
2218  Vector &u, const FiniteElement &fluxelem,
2219  Vector &flux, bool with_coef);
2220 
2221  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
2223  Vector &flux, Vector *d_energy = NULL);
2224 
2226  virtual void AssemblePA(const FiniteElementSpace &fes);
2227  virtual void AddMultPA(const Vector &x, Vector &y) const;
2228  virtual void AssembleDiagonalPA(Vector& diag);
2229 };
2230 
2231 /** Integrator for (curl u, curl v) for FE spaces defined by 'dim' copies of a
2232  scalar FE space. */
2234 {
2235 private:
2236 #ifndef MFEM_THREAD_SAFE
2237  DenseMatrix dshape_hat, dshape, curlshape, Jadj, grad_hat, grad;
2238 #endif
2239 
2240 protected:
2242 
2243 public:
2245 
2247 
2248  /// Assemble an element matrix
2249  virtual void AssembleElementMatrix(const FiniteElement &el,
2251  DenseMatrix &elmat);
2252  /// Compute element energy: (1/2) (curl u, curl u)_E
2253  virtual double GetElementEnergy(const FiniteElement &el,
2255  const Vector &elfun);
2256 };
2257 
2258 /** Integrator for (Q u, v), where Q is an optional coefficient (of type scalar,
2259  vector (diagonal matrix), or matrix), trial function u is in H(Curl) or
2260  H(Div), and test function v is in H(Curl), H(Div), or v=(v1,...,vn), where
2261  vi are in H1. */
2263 {
2264 private:
2266  { Q = q; VQ = vq; MQ = mq; }
2267 
2268 #ifndef MFEM_THREAD_SAFE
2269  Vector shape;
2270  Vector D;
2271  DenseMatrix K;
2272  DenseMatrix partelmat;
2273  DenseMatrix test_vshape;
2274  DenseMatrix trial_vshape;
2275 #endif
2276 
2277 protected:
2281 
2282  // PA extension
2284  const DofToQuad *mapsO; ///< Not owned. DOF-to-quad map, open.
2285  const DofToQuad *mapsC; ///< Not owned. DOF-to-quad map, closed.
2286  const GeometricFactors *geom; ///< Not owned
2287  int dim, ne, nq, dofs1D, quad1D;
2288 
2289 public:
2290  VectorFEMassIntegrator() { Init(NULL, NULL, NULL); }
2291  VectorFEMassIntegrator(Coefficient *_q) { Init(_q, NULL, NULL); }
2292  VectorFEMassIntegrator(Coefficient &q) { Init(&q, NULL, NULL); }
2293  VectorFEMassIntegrator(VectorCoefficient *_vq) { Init(NULL, _vq, NULL); }
2294  VectorFEMassIntegrator(VectorCoefficient &vq) { Init(NULL, &vq, NULL); }
2295  VectorFEMassIntegrator(MatrixCoefficient *_mq) { Init(NULL, NULL, _mq); }
2296  VectorFEMassIntegrator(MatrixCoefficient &mq) { Init(NULL, NULL, &mq); }
2297 
2298  virtual void AssembleElementMatrix(const FiniteElement &el,
2300  DenseMatrix &elmat);
2301  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2302  const FiniteElement &test_fe,
2304  DenseMatrix &elmat);
2305 
2307  virtual void AssemblePA(const FiniteElementSpace &fes);
2308  virtual void AddMultPA(const Vector &x, Vector &y) const;
2309  virtual void AssembleDiagonalPA(Vector& diag);
2310 };
2311 
2312 /** Integrator for (Q div u, p) where u=(v1,...,vn) and all vi are in the same
2313  scalar FE space; p is also in a (different) scalar FE space. */
2315 {
2316 protected:
2318 
2319 private:
2320  Vector shape;
2321  Vector divshape;
2322  DenseMatrix dshape;
2323  DenseMatrix gshape;
2324  DenseMatrix Jadj;
2325  // PA extension
2326  Vector pa_data;
2327  const DofToQuad *trial_maps, *test_maps; ///< Not owned
2328  const GeometricFactors *geom; ///< Not owned
2329  int dim, ne, nq;
2330  int trial_dofs1D, test_dofs1D, quad1D;
2331 
2332 public:
2334  Q(NULL), trial_maps(NULL), test_maps(NULL), geom(NULL)
2335  { }
2337  Q(_q), trial_maps(NULL), test_maps(NULL), geom(NULL)
2338  { }
2340  Q(&q), trial_maps(NULL), test_maps(NULL), geom(NULL)
2341  { }
2342 
2343  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2344  const FiniteElement &test_fe,
2346  DenseMatrix &elmat);
2347 
2349  virtual void AssemblePA(const FiniteElementSpace &trial_fes,
2350  const FiniteElementSpace &test_fes);
2351 
2352  virtual void AddMultPA(const Vector &x, Vector &y) const;
2353  virtual void AddMultTransposePA(const Vector &x, Vector &y) const;
2354 
2355  static const IntegrationRule &GetRule(const FiniteElement &trial_fe,
2356  const FiniteElement &test_fe,
2358 };
2359 
2360 /// (Q div u, div v) for RT elements
2362 {
2363 protected:
2365 
2366 private:
2367 #ifndef MFEM_THREAD_SAFE
2368  Vector divshape;
2369 #endif
2370 
2371 public:
2372  DivDivIntegrator() { Q = NULL; }
2374 
2375  virtual void AssembleElementMatrix(const FiniteElement &el,
2377  DenseMatrix &elmat);
2378 };
2379 
2380 /** Integrator for
2381  (Q grad u, grad v) = sum_i (Q grad u_i, grad v_i) e_i e_i^T
2382  for FE spaces defined by 'dim' copies of a scalar FE space. Where e_i
2383  is the unit vector in the i-th direction. The resulting local element
2384  matrix is a block-diagonal matrix consisting of 'dim' copies of a scalar
2385  diffusion matrix in each diagonal block. */
2387 {
2388 protected:
2390 
2391  // PA extension
2392  const DofToQuad *maps; ///< Not owned
2393  const GeometricFactors *geom; ///< Not owned
2396 
2397 private:
2398  DenseMatrix Jinv;
2399  DenseMatrix dshape;
2400  DenseMatrix gshape;
2401  DenseMatrix pelmat;
2402 
2403 public:
2406 
2407  virtual void AssembleElementMatrix(const FiniteElement &el,
2409  DenseMatrix &elmat);
2410  virtual void AssembleElementVector(const FiniteElement &el,
2412  const Vector &elfun, Vector &elvect);
2414  virtual void AssemblePA(const FiniteElementSpace &fes);
2415  virtual void AssembleDiagonalPA(Vector &diag);
2416  virtual void AddMultPA(const Vector &x, Vector &y) const;
2417 };
2418 
2419 /** Integrator for the linear elasticity form:
2420  a(u,v) = (lambda div(u), div(v)) + (2 mu e(u), e(v)),
2421  where e(v) = (1/2) (grad(v) + grad(v)^T).
2422  This is a 'Vector' integrator, i.e. defined for FE spaces
2423  using multiple copies of a scalar FE space. */
2425 {
2426 protected:
2427  double q_lambda, q_mu;
2429 
2430 private:
2431 #ifndef MFEM_THREAD_SAFE
2432  Vector shape;
2433  DenseMatrix dshape, gshape, pelmat;
2434  Vector divshape;
2435 #endif
2436 
2437 public:
2439  { lambda = &l; mu = &m; }
2440  /** With this constructor lambda = q_l * m and mu = q_m * m;
2441  if dim * q_l + 2 * q_m = 0 then trace(sigma) = 0. */
2442  ElasticityIntegrator(Coefficient &m, double q_l, double q_m)
2443  { lambda = NULL; mu = &m; q_lambda = q_l; q_mu = q_m; }
2444 
2445  virtual void AssembleElementMatrix(const FiniteElement &,
2447  DenseMatrix &);
2448 
2449  /** Compute the stress corresponding to the local displacement @a u and
2450  interpolate it at the nodes of the given @a fluxelem. Only the symmetric
2451  part of the stress is stored, so that the size of @a flux is equal to
2452  the number of DOFs in @a fluxelem times dim*(dim+1)/2. In 2D, the order
2453  of the stress components is: s_xx, s_yy, s_xy. In 3D, it is: s_xx, s_yy,
2454  s_zz, s_xy, s_xz, s_yz. In other words, @a flux is the local vector for
2455  a FE space with dim*(dim+1)/2 vector components, based on the finite
2456  element @a fluxelem. */
2457  virtual void ComputeElementFlux(const FiniteElement &el,
2459  Vector &u,
2460  const FiniteElement &fluxelem,
2461  Vector &flux, bool with_coef = true);
2462 
2463  /** Compute the element energy (integral of the strain energy density)
2464  corresponding to the stress represented by @a flux which is a vector of
2465  coefficients multiplying the basis functions defined by @a fluxelem. In
2466  other words, @a flux is the local vector for a FE space with
2467  dim*(dim+1)/2 vector components, based on the finite element @a fluxelem.
2468  The number of components, dim*(dim+1)/2 is such that it represents the
2469  symmetric part of the (symmetric) stress tensor. The order of the
2470  components is: s_xx, s_yy, s_xy in 2D, and s_xx, s_yy, s_zz, s_xy, s_xz,
2471  s_yz in 3D. */
2472  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
2474  Vector &flux, Vector *d_energy = NULL);
2475 };
2476 
2477 /** Integrator for the DG form:
2478  alpha < rho_u (u.n) {v},[w] > + beta < rho_u |u.n| [v],[w] >,
2479  where v and w are the trial and test variables, respectively, and rho/u are
2480  given scalar/vector coefficients. The vector coefficient, u, is assumed to
2481  be continuous across the faces and when given the scalar coefficient, rho,
2482  is assumed to be discontinuous. The integrator uses the upwind value of rho,
2483  rho_u, which is value from the side into which the vector coefficient, u,
2484  points. */
2486 {
2487 protected:
2490  double alpha, beta;
2491  // PA extension
2493  const DofToQuad *maps; ///< Not owned
2494  const FaceGeometricFactors *geom; ///< Not owned
2495  int dim, nf, nq, dofs1D, quad1D;
2496 
2497 private:
2498  Vector shape1, shape2;
2499 
2500 public:
2501  /// Construct integrator with rho = 1.
2503  { rho = NULL; u = &_u; alpha = a; beta = b; }
2504 
2506  double a, double b)
2507  { rho = &_rho; u = &_u; alpha = a; beta = b; }
2508 
2510  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2511  const FiniteElement &el2,
2513  DenseMatrix &elmat);
2514 
2516 
2517  virtual void AssemblePAInteriorFaces(const FiniteElementSpace &fes);
2518 
2519  virtual void AssemblePABoundaryFaces(const FiniteElementSpace &fes);
2520 
2521  virtual void AddMultTransposePA(const Vector &x, Vector &y) const;
2522 
2523  virtual void AddMultPA(const Vector&, Vector&) const;
2524 
2525  static const IntegrationRule &GetRule(Geometry::Type geom, int order,
2527 
2528 private:
2529  void SetupPA(const FiniteElementSpace &fes, FaceType type);
2530 };
2531 
2532 /** Integrator for the DG form:
2533 
2534  - < {(Q grad(u)).n}, [v] > + sigma < [u], {(Q grad(v)).n} >
2535  + kappa < {h^{-1} Q} [u], [v] >,
2536 
2537  where Q is a scalar or matrix diffusion coefficient and u, v are the trial
2538  and test spaces, respectively. The parameters sigma and kappa determine the
2539  DG method to be used (when this integrator is added to the "broken"
2540  DiffusionIntegrator):
2541  * sigma = -1, kappa >= kappa0: symm. interior penalty (IP or SIPG) method,
2542  * sigma = +1, kappa > 0: non-symmetric interior penalty (NIPG) method,
2543  * sigma = +1, kappa = 0: the method of Baumann and Oden. */
2545 {
2546 protected:
2549  double sigma, kappa;
2550 
2551  // these are not thread-safe!
2554 
2555 public:
2556  DGDiffusionIntegrator(const double s, const double k)
2557  : Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
2558  DGDiffusionIntegrator(Coefficient &q, const double s, const double k)
2559  : Q(&q), MQ(NULL), sigma(s), kappa(k) { }
2560  DGDiffusionIntegrator(MatrixCoefficient &q, const double s, const double k)
2561  : Q(NULL), MQ(&q), sigma(s), kappa(k) { }
2563  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2564  const FiniteElement &el2,
2566  DenseMatrix &elmat);
2567 };
2568 
2569 /** Integrator for the DG elasticity form, for the formulations see:
2570  - PhD Thesis of Jonas De Basabe, High-Order Finite %Element Methods for
2571  Seismic Wave Propagation, UT Austin, 2009, p. 23, and references therein
2572  - Peter Hansbo and Mats G. Larson, Discontinuous Galerkin and the
2573  Crouzeix-Raviart %Element: Application to Elasticity, PREPRINT 2000-09,
2574  p.3
2575 
2576  \f[
2577  - \left< \{ \tau(u) \}, [v] \right> + \alpha \left< \{ \tau(v) \}, [u]
2578  \right> + \kappa \left< h^{-1} \{ \lambda + 2 \mu \} [u], [v] \right>
2579  \f]
2580 
2581  where \f$ \left<u, v\right> = \int_{F} u \cdot v \f$, and \f$ F \f$ is a
2582  face which is either a boundary face \f$ F_b \f$ of an element \f$ K \f$ or
2583  an interior face \f$ F_i \f$ separating elements \f$ K_1 \f$ and \f$ K_2 \f$.
2584 
2585  In the bilinear form above \f$ \tau(u) \f$ is traction, and it's also
2586  \f$ \tau(u) = \sigma(u) \cdot \vec{n} \f$, where \f$ \sigma(u) \f$ is
2587  stress, and \f$ \vec{n} \f$ is the unit normal vector w.r.t. to \f$ F \f$.
2588 
2589  In other words, we have
2590  \f[
2591  - \left< \{ \sigma(u) \cdot \vec{n} \}, [v] \right> + \alpha \left< \{
2592  \sigma(v) \cdot \vec{n} \}, [u] \right> + \kappa \left< h^{-1} \{
2593  \lambda + 2 \mu \} [u], [v] \right>
2594  \f]
2595 
2596  For isotropic media
2597  \f[
2598  \begin{split}
2599  \sigma(u) &= \lambda \nabla \cdot u I + 2 \mu \varepsilon(u) \\
2600  &= \lambda \nabla \cdot u I + 2 \mu \frac{1}{2} (\nabla u + \nabla
2601  u^T) \\
2602  &= \lambda \nabla \cdot u I + \mu (\nabla u + \nabla u^T)
2603  \end{split}
2604  \f]
2605 
2606  where \f$ I \f$ is identity matrix, \f$ \lambda \f$ and \f$ \mu \f$ are Lame
2607  coefficients (see ElasticityIntegrator), \f$ u, v \f$ are the trial and test
2608  functions, respectively.
2609 
2610  The parameters \f$ \alpha \f$ and \f$ \kappa \f$ determine the DG method to
2611  use (when this integrator is added to the "broken" ElasticityIntegrator):
2612 
2613  - IIPG, \f$\alpha = 0\f$,
2614  C. Dawson, S. Sun, M. Wheeler, Compatible algorithms for coupled flow and
2615  transport, Comp. Meth. Appl. Mech. Eng., 193(23-26), 2565-2580, 2004.
2616 
2617  - SIPG, \f$\alpha = -1\f$,
2618  M. Grote, A. Schneebeli, D. Schotzau, Discontinuous Galerkin Finite
2619  %Element Method for the Wave Equation, SINUM, 44(6), 2408-2431, 2006.
2620 
2621  - NIPG, \f$\alpha = 1\f$,
2622  B. Riviere, M. Wheeler, V. Girault, A Priori Error Estimates for Finite
2623  %Element Methods Based on Discontinuous Approximation Spaces for Elliptic
2624  Problems, SINUM, 39(3), 902-931, 2001.
2625 
2626  This is a '%Vector' integrator, i.e. defined for FE spaces using multiple
2627  copies of a scalar FE space.
2628  */
2630 {
2631 public:
2632  DGElasticityIntegrator(double alpha_, double kappa_)
2633  : lambda(NULL), mu(NULL), alpha(alpha_), kappa(kappa_) { }
2634 
2636  double alpha_, double kappa_)
2637  : lambda(&lambda_), mu(&mu_), alpha(alpha_), kappa(kappa_) { }
2638 
2640  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2641  const FiniteElement &el2,
2643  DenseMatrix &elmat);
2644 
2645 protected:
2647  double alpha, kappa;
2648 
2649 #ifndef MFEM_THREAD_SAFE
2650  // values of all scalar basis functions for one component of u (which is a
2651  // vector) at the integration point in the reference space
2653  // values of derivatives of all scalar basis functions for one component
2654  // of u (which is a vector) at the integration point in the reference space
2656  // Adjugate of the Jacobian of the transformation: adjJ = det(J) J^{-1}
2658  // gradient of shape functions in the real (physical, not reference)
2659  // coordinates, scaled by det(J):
2660  // dshape_ps(jdof,jm) = sum_{t} adjJ(t,jm)*dshape(jdof,t)
2662  Vector nor; // nor = |weight(J_face)| n
2663  Vector nL1, nL2; // nL1 = (lambda1 * ip.weight / detJ1) nor
2664  Vector nM1, nM2; // nM1 = (mu1 * ip.weight / detJ1) nor
2665  Vector dshape1_dnM, dshape2_dnM; // dshape1_dnM = dshape1_ps . nM1
2666  // 'jmat' corresponds to the term: kappa <h^{-1} {lambda + 2 mu} [u], [v]>
2668 #endif
2669 
2670  static void AssembleBlock(
2671  const int dim, const int row_ndofs, const int col_ndofs,
2672  const int row_offset, const int col_offset,
2673  const double jmatcoef, const Vector &col_nL, const Vector &col_nM,
2674  const Vector &row_shape, const Vector &col_shape,
2675  const Vector &col_dshape_dnM, const DenseMatrix &col_dshape,
2676  DenseMatrix &elmat, DenseMatrix &jmat);
2677 };
2678 
2679 /** Integrator for the DPG form: < v, [w] > over all faces (the interface) where
2680  the trial variable v is defined on the interface and the test variable w is
2681  defined inside the elements, generally in a DG space. */
2683 {
2684 private:
2685  Vector face_shape, shape1, shape2;
2686 
2687 public:
2690  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
2691  const FiniteElement &test_fe1,
2692  const FiniteElement &test_fe2,
2694  DenseMatrix &elmat);
2695 };
2696 
2697 /** Integrator for the form: < v, [w.n] > over all faces (the interface) where
2698  the trial variable v is defined on the interface and the test variable w is
2699  in an H(div)-conforming space. */
2701 {
2702 private:
2703  Vector face_shape, normal, shape1_n, shape2_n;
2704  DenseMatrix shape1, shape2;
2705 
2706 public:
2709  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
2710  const FiniteElement &test_fe1,
2711  const FiniteElement &test_fe2,
2713  DenseMatrix &elmat);
2714 };
2715 
2716 /** Abstract class to serve as a base for local interpolators to be used in the
2717  DiscreteLinearOperator class. */
2719 
2720 
2721 /** Class for constructing the gradient as a DiscreteLinearOperator from an
2722  H1-conforming space to an H(curl)-conforming space. The range space can be
2723  vector L2 space as well. */
2725 {
2726 public:
2727  virtual void AssembleElementMatrix2(const FiniteElement &h1_fe,
2728  const FiniteElement &nd_fe,
2730  DenseMatrix &elmat)
2731  { nd_fe.ProjectGrad(h1_fe, Trans, elmat); }
2732 };
2733 
2734 
2735 /** Class for constructing the identity map as a DiscreteLinearOperator. This
2736  is the discrete embedding matrix when the domain space is a subspace of
2737  the range space. Otherwise, a dof projection matrix is constructed. */
2739 {
2740 public:
2741  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2742  const FiniteElement &ran_fe,
2744  DenseMatrix &elmat)
2745  { ran_fe.Project(dom_fe, Trans, elmat); }
2746 };
2747 
2748 
2749 /** Class for constructing the (local) discrete curl matrix which can be used
2750  as an integrator in a DiscreteLinearOperator object to assemble the global
2751  discrete curl matrix. */
2753 {
2754 public:
2755  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2756  const FiniteElement &ran_fe,
2758  DenseMatrix &elmat)
2759  { ran_fe.ProjectCurl(dom_fe, Trans, elmat); }
2760 };
2761 
2762 
2763 /** Class for constructing the (local) discrete divergence matrix which can
2764  be used as an integrator in a DiscreteLinearOperator object to assemble
2765  the global discrete divergence matrix.
2766 
2767  Note: Since the dofs in the L2_FECollection are nodal values, the local
2768  discrete divergence matrix (with an RT-type domain space) will depend on
2769  the transformation. On the other hand, the local matrix returned by
2770  VectorFEDivergenceIntegrator is independent of the transformation. */
2772 {
2773 public:
2774  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2775  const FiniteElement &ran_fe,
2777  DenseMatrix &elmat)
2778  { ran_fe.ProjectDiv(dom_fe, Trans, elmat); }
2779 };
2780 
2781 
2782 /** A trace face interpolator class for interpolating the normal component of
2783  the domain space, e.g. vector H1, into the range space, e.g. the trace of
2784  RT which uses FiniteElement::INTEGRAL map type. */
2786 {
2787 public:
2788  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2789  const FiniteElement &ran_fe,
2791  DenseMatrix &elmat);
2792 };
2793 
2794 /** Interpolator of a scalar coefficient multiplied by a scalar field onto
2795  another scalar field. Note that this can produce inaccurate fields unless
2796  the target is sufficiently high order. */
2798 {
2799 public:
2801 
2802  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2803  const FiniteElement &ran_fe,
2805  DenseMatrix &elmat);
2806 
2807 protected:
2809 };
2810 
2811 /** Interpolator of a scalar coefficient multiplied by a vector field onto
2812  another vector field. Note that this can produce inaccurate fields unless
2813  the target is sufficiently high order. */
2815 {
2816 public:
2818  : Q(&sc) { }
2819 
2820  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2821  const FiniteElement &ran_fe,
2823  DenseMatrix &elmat);
2824 protected:
2826 };
2827 
2828 /** Interpolator of a vector coefficient multiplied by a scalar field onto
2829  another vector field. Note that this can produce inaccurate fields unless
2830  the target is sufficiently high order. */
2832 {
2833 public:
2835  : VQ(&vc) { }
2836 
2837  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2838  const FiniteElement &ran_fe,
2840  DenseMatrix &elmat);
2841 protected:
2843 };
2844 
2845 /** Interpolator of the cross product between a vector coefficient and an
2846  H(curl)-conforming field onto an H(div)-conforming field. The range space
2847  can also be vector L2. */
2849 {
2850 public:
2852  : VQ(&vc) { }
2853 
2854  virtual void AssembleElementMatrix2(const FiniteElement &nd_fe,
2855  const FiniteElement &rt_fe,
2857  DenseMatrix &elmat);
2858 protected:
2860 };
2861 
2862 /** Interpolator of the inner product between a vector coefficient and an
2863  H(div)-conforming field onto an L2-conforming field. The range space can
2864  also be H1. */
2866 {
2867 public:
2869 
2870  virtual void AssembleElementMatrix2(const FiniteElement &rt_fe,
2871  const FiniteElement &l2_fe,
2873  DenseMatrix &elmat);
2874 protected:
2876 };
2877 
2878 }
2879 
2880 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:232
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:588
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
Definition: bilininteg.hpp:308
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:981
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
Definition: bilininteg.hpp:375
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
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:108
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:308
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:743
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const DofToQuad * mapsO
Not owned. DOF-to-quad map, open.
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:714
MixedCrossCurlGradIntegrator(VectorCoefficient &vq)
DGTraceIntegrator(VectorCoefficient &_u, double a, double b)
Construct integrator with rho = 1.
virtual void AssemblePAInteriorFaces(const FiniteElementSpace &fes)
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:276
SumIntegrator(int own_integs=1)
Definition: bilininteg.hpp:284
static const IntegrationRule & GetRule(Geometry::Type geom, int order, FaceElementTransformations &T)
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:40
virtual const char * FiniteElementTypeFailureMessage() const
MatrixCoefficient * MQ
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:396
MixedScalarCrossProductIntegrator(VectorCoefficient &vq)
CurlCurlIntegrator(MatrixCoefficient &m)
MatrixCoefficient * MQ
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:988
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Definition: fe.cpp:177
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)
VectorCoefficient * u
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:491
VectorCoefficient * VQ
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:130
MixedVectorWeakDivergenceIntegrator(VectorCoefficient &dq)
VectorCoefficient * VQ
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Virtual method required for Zienkiewicz-Zhu type error estimators.
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:61
virtual void AssemblePAInteriorFaces(const FiniteElementSpace &fes)
Definition: bilininteg.hpp:217
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:874
ElasticityIntegrator(Coefficient &l, Coefficient &m)
VectorFEMassIntegrator(VectorCoefficient *_vq)
VectorCoefficient * VQ
Definition: bilininteg.hpp:427
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:342
const GeometricFactors * geom
Not owned.
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:569
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:605
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
virtual void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
MixedScalarMassIntegrator(Coefficient &q)
Definition: bilininteg.hpp:529
VectorMassIntegrator(Coefficient &q, const IntegrationRule *ir)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:957
alpha (q . grad u, v) using the &quot;group&quot; FE discretization
const GeometricFactors * geom
Not owned.
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
void SetupPA(const FiniteElementSpace &fes, const bool force=false)
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:318
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:699
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
DivDivIntegrator(Coefficient &q)
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
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)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:125
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:134
const DofToQuad * maps
Not owned.
BilinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: bilininteg.hpp:27
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
Definition: bilininteg.hpp:212
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:328
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:408
void AddIntegrator(BilinearFormIntegrator *integ)
Definition: bilininteg.hpp:286
const DofToQuad * mapsC
Not owned. DOF-to-quad map, closed.
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:935
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1332
MixedVectorCurlIntegrator(VectorCoefficient &dq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:347
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef)
Virtual method required for Zienkiewicz-Zhu type error estimators.
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
Polynomials of order k.
Definition: fe.hpp:215
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Definition: bilininteg.hpp:175
MixedCrossProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:820
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
Definition: bilininteg.cpp:60
MixedScalarVectorIntegrator(VectorCoefficient &vq, bool _transpose=false, bool _cross_2d=false)
Definition: bilininteg.hpp:458
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:161
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:166
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Implements CalcDivShape methods.
Definition: fe.hpp:293
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:858
void SetupPA(const FiniteElementSpace &fes, const bool force=false)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:475
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
MixedScalarCrossGradIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:952
const DofToQuad * maps
Not owned.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:871
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:925
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:75
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:172
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
const GeometricFactors * geom
Not owned.
MixedCurlCurlIntegrator(MatrixCoefficient &mq)
MatrixCoefficient * MQ
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
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.
Structure for storing face geometric factors: coordinates, Jacobians, determinants of the Jacobians...
Definition: mesh.hpp:1376
VectorCrossProductInterpolator(VectorCoefficient &vc)
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:322
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:502
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
Definition: bilininteg.cpp:23
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:340
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
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:116
MixedGradGradIntegrator(VectorCoefficient &dq)
Definition: bilininteg.hpp:920
Integrator that inverts the matrix assembled by another integrator.
Definition: bilininteg.hpp:258
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:169
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
VectorCoefficient * DQ
Definition: bilininteg.hpp:428
MatrixCoefficient * MQ
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
Implements CalcCurlShape methods.
Definition: fe.hpp:294
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
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:787
VectorInnerProductInterpolator(VectorCoefficient &vc)
MixedWeakDivCrossIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:885
virtual const char * FiniteElementTypeFailureMessage() const
FaceType
Definition: mesh.hpp:42
virtual const char * FiniteElementTypeFailureMessage() const
VectorDivergenceIntegrator(Coefficient &q)
MixedCrossCurlCurlIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:141
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
MatrixCoefficient * MQ
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedVectorWeakCurlIntegrator(VectorCoefficient &dq)
MixedScalarWeakCurlIntegrator(Coefficient &q)
Definition: bilininteg.hpp:771
MixedCurlCurlIntegrator(Coefficient &q)
GradientIntegrator(Coefficient &q)
double b
Definition: lissajous.cpp:42
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:626
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:148
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
virtual const char * FiniteElementTypeFailureMessage() const
const DofToQuad * mapsC
Not owned. DOF-to-quad map, closed.
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: bilininteg.hpp:110
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.cpp:856
VectorFEDivergenceIntegrator(Coefficient &q)
virtual void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
Definition: nonlininteg.cpp:18
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:888
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
VectorFECurlIntegrator(Coefficient &q)
Implements CalcDShape methods.
Definition: fe.hpp:292
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 const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:404
DiffusionIntegrator()
Construct a diffusion integrator with coefficient Q = 1.
MixedVectorMassIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:811
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:622
MixedVectorIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:393
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:416
MixedVectorMassIntegrator(Coefficient &q)
Definition: bilininteg.hpp:807
virtual void AssemblePABoundaryFaces(const FiniteElementSpace &fes)
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:195
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:677
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
ScalarProductInterpolator(Coefficient &sc)
virtual const char * FiniteElementTypeFailureMessage() const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:775
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)
MixedVectorCurlIntegrator(MatrixCoefficient &mq)
MixedCurlCurlIntegrator(VectorCoefficient &dq)
virtual void AssemblePAInteriorFaces(const FiniteElementSpace &fes)
Definition: bilininteg.cpp:36
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:791
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:634
BoundaryMassIntegrator(Coefficient &q)
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:962
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
Definition: bilininteg.cpp:99
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
Definition: bilininteg.cpp:54
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
Definition: bilininteg.cpp:584
MixedWeakGradDotIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:855
ScalarVectorProductInterpolator(Coefficient &sc)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:841
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:682
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
Definition: bilininteg.hpp:232
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
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:421
const GeometricFactors * geom
Not owned.
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorDivergenceIntegrator(Coefficient *_q)
MixedVectorIntegrator(Coefficient &q)
Definition: bilininteg.hpp:388
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:784
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:922
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:914
VectorFEMassIntegrator(MatrixCoefficient &mq)
MixedVectorMassIntegrator(VectorCoefficient &dq)
Definition: bilininteg.hpp:809
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
MixedScalarWeakCurlCrossIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:734
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:66
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:662
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)
static const IntegrationRule & GetRule(const FiniteElement &el, ElementTransformation &Trans)
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssemblePABoundaryFaces(const FiniteElementSpace &fes)
Definition: bilininteg.hpp:222
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:185
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
MixedScalarWeakDivergenceIntegrator(VectorCoefficient &vq)
const FaceGeometricFactors * geom
Not owned.
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:525
DGDiffusionIntegrator(const double s, const double k)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:905
LumpedIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
Definition: bilininteg.hpp:247
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe)
Definition: bilininteg.cpp:764
MixedCrossGradCurlIntegrator(VectorCoefficient &vq)
MixedVectorDivergenceIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:658
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe.hpp:128
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
Definition: bilininteg.cpp:48
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:54
MassIntegrator(Coefficient &q, const IntegrationRule *ir=NULL)
Construct a mass integrator with coefficient q.
const FiniteElementSpace * fespace
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:330
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
MixedCrossGradGradIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:968
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:641
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedVectorIntegrator(VectorCoefficient &dq, bool diag=true)
Definition: bilininteg.hpp:390
MixedScalarCrossCurlIntegrator(VectorCoefficient &vq)
MatrixCoefficient * MQ
Definition: bilininteg.hpp:429
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Assemble an element matrix.
double a
Definition: lissajous.cpp:41
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:82
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:707
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:867
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:898
VectorMassIntegrator(Coefficient &q, int qo=0)
VectorFEMassIntegrator(Coefficient *_q)
DGElasticityIntegrator(double alpha_, double kappa_)
MixedScalarDerivativeIntegrator(Coefficient &q)
Definition: bilininteg.hpp:549
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorFEMassIntegrator(MatrixCoefficient *_mq)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Definition: bilininteg.cpp:640
MixedDotProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:830
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorCurlIntegrator(Coefficient &q)
const GeometricFactors * geom
Not owned.
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
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:597
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
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:819
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
int dim
Definition: ex24.cpp:43
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
const DofToQuad * maps
Not owned.
virtual void AssembleElementMatrix2(const FiniteElement &nd_fe, const FiniteElement &rt_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:231
VectorMassIntegrator(MatrixCoefficient &q, int qo=0)
Construct an integrator with matrix coefficient q.
GradientIntegrator(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:410
GroupConvectionIntegrator(VectorCoefficient &q, double a=1.0)
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:497
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:993
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:161
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 void AssembleElementMatrix2(const FiniteElement &rt_fe, const FiniteElement &l2_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:562
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
DGDiffusionIntegrator(MatrixCoefficient &q, const double s, const double k)
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
const DofToQuad * maps
Not owned.
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:755
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
Compute element energy: (1/2) (curl u, curl u)_E.
virtual void AssembleFaceGrad(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integr...
Definition: bilininteg.hpp:115
MixedCrossCurlIntegrator(VectorCoefficient &vq)
Vector data type.
Definition: vector.hpp:48
int GetDerivType() const
Definition: fe.hpp:336
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:670
MixedVectorProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:539
VectorScalarProductInterpolator(VectorCoefficient &vc)
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:730
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:74
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Definition: bilininteg.hpp:150
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleDiagonalPA(Vector &diag)
Assemble diagonal and add it to Vector diag.
DerivativeIntegrator(Coefficient &q, int i)
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.cpp:464
DGTraceIntegrator(Coefficient &_rho, VectorCoefficient &_u, double a, double b)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:833
CurlCurlIntegrator(Coefficient &q)
Construct a bilinear form integrator for Nedelec elements.
virtual const char * FiniteElementTypeFailureMessage() const
MixedGradGradIntegrator(Coefficient &q)
Definition: bilininteg.hpp:918
MixedVectorGradientIntegrator(MatrixCoefficient &mq)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:942
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:330
const DofToQuad * mapsO
Not owned. DOF-to-quad map, open.
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:646
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
Definition: bilininteg.hpp:227
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:475
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
MassIntegrator(const IntegrationRule *ir=NULL)
const GeometricFactors * geom
Not owned.
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:336
virtual double ComputeFluxEnergy(const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Definition: bilininteg.cpp:698
MixedVectorWeakCurlIntegrator(Coefficient &q)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:971
TransposeIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
Definition: bilininteg.hpp:192
MixedDivGradIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:553
virtual void AssemblePA(const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:750
MixedScalarIntegrator(Coefficient &q)
Definition: bilininteg.hpp:320
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:265
virtual void AssemblePABoundaryFaces(const FiniteElementSpace &fes)
Definition: bilininteg.cpp:42
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:462
const DofToQuad * maps
Not owned.
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
VectorCoefficient * Q
alpha (q . grad u, v)