MFEM  v4.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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_BILININTEG
13 #define MFEM_BILININTEG
14 
15 #include "../config/config.hpp"
16 #include "nonlininteg.hpp"
17 #include "fespace.hpp"
18 
19 namespace mfem
20 {
21 
22 /// Abstract base class BilinearFormIntegrator
24 {
25 protected:
27  : NonlinearFormIntegrator(ir) { }
28 
29 public:
30  // TODO: add support for other assembly levels (in addition to PA) and their
31  // actions.
32 
33  // TODO: for mixed meshes the quadrature rules to be used by methods like
34  // AssemblePA() can be given as a QuadratureSpace, e.g. using a new method:
35  // SetQuadratureSpace().
36 
37  // TODO: the methods for the various assembly levels make sense even in the
38  // base class NonlinearFormIntegrator, except that not all assembly levels
39  // make sense for the action of the nonlinear operator (but they all make
40  // sense for its Jacobian).
41 
42  /// Method defining partial assembly.
43  /** The result of the partial assembly is stored internally so that it can be
44  used later in the methods AddMultPA() and AddMultTransposePA(). */
45  virtual void AssemblePA(const FiniteElementSpace &fes);
46 
47  /// Method for partially assembled action.
48  /** Perform the action of integrator on the input @a x and add the result to
49  the output @a y. Both @a x and @a y are E-vectors, i.e. they represent
50  the element-wise discontinuous version of the FE space.
51 
52  This method can be called only after the method AssemblePA() has been
53  called. */
54  virtual void AddMultPA(const Vector &x, Vector &y) const;
55 
56  /// Method for partially assembled transposed action.
57  /** Perform the transpose action of integrator on the input @a x and add the
58  result to the output @a y. Both @a x and @a y are E-vectors, i.e. they
59  represent the element-wise discontinuous version of the FE space.
60 
61  This method can be called only after the method AssemblePA() has been
62  called. */
63  virtual void AddMultTransposePA(const Vector &x, Vector &y) const;
64 
65  /// Given a particular Finite Element computes the element matrix elmat.
66  virtual void AssembleElementMatrix(const FiniteElement &el,
68  DenseMatrix &elmat);
69 
70  /** Compute the local matrix representation of a bilinear form
71  a(u,v) defined on different trial (given by u) and test
72  (given by v) spaces. The rows in the local matrix correspond
73  to the test dofs and the columns -- to the trial dofs. */
74  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
75  const FiniteElement &test_fe,
77  DenseMatrix &elmat);
78 
79  virtual void AssembleFaceMatrix(const FiniteElement &el1,
80  const FiniteElement &el2,
82  DenseMatrix &elmat);
83 
84  /** Abstract method used for assembling TraceFaceIntegrators in a
85  MixedBilinearForm. */
86  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
87  const FiniteElement &test_fe1,
88  const FiniteElement &test_fe2,
90  DenseMatrix &elmat);
91 
92  /// Perform the local action of the BilinearFormIntegrator
93  virtual void AssembleElementVector(const FiniteElement &el,
95  const Vector &elfun, Vector &elvect);
96 
97  virtual void AssembleElementGrad(const FiniteElement &el,
99  const Vector &elfun, DenseMatrix &elmat)
100  { AssembleElementMatrix(el, Tr, elmat); }
101 
102  virtual void AssembleFaceGrad(const FiniteElement &el1,
103  const FiniteElement &el2,
105  const Vector &elfun, DenseMatrix &elmat)
106  { AssembleFaceMatrix(el1, el2, Tr, elmat); }
107 
108  /** @brief Virtual method required for Zienkiewicz-Zhu type error estimators.
109 
110  The purpose of the method is to compute a local "flux" finite element
111  function given a local finite element solution. The "flux" function has
112  to be computed in terms of its coefficients (represented by the Vector
113  @a flux) which multiply the basis functions defined by the FiniteElement
114  @a fluxelem. Typically, the "flux" function will have more than one
115  component and consequently @a flux should be store the coefficients of
116  all components: first all coefficient for component 0, then all
117  coefficients for component 1, etc. What the "flux" function represents
118  depends on the specific integrator. For example, in the case of
119  DiffusionIntegrator, the flux is the gradient of the solution multiplied
120  by the diffusion coefficient.
121 
122  @param[in] el FiniteElement of the solution.
123  @param[in] Trans The ElementTransformation describing the physical
124  position of the mesh element.
125  @param[in] u Solution coefficients representing the expansion of the
126  solution function in the basis of @a el.
127  @param[in] fluxelem FiniteElement of the "flux".
128  @param[out] flux "Flux" coefficients representing the expansion of the
129  "flux" function in the basis of @a fluxelem. The size
130  of @a flux as a Vector has to be set by this method,
131  e.g. using Vector::SetSize().
132  @param[in] with_coef If zero (the default value is 1) the implementation
133  of the method may choose not to scale the "flux"
134  function by any coefficients describing the
135  integrator.
136  */
137  virtual void ComputeElementFlux(const FiniteElement &el,
139  Vector &u,
140  const FiniteElement &fluxelem,
141  Vector &flux, int with_coef = 1) { }
142 
143  /** @brief Virtual method required for Zienkiewicz-Zhu type error estimators.
144 
145  The purpose of this method is to compute a local number that measures the
146  energy of a given "flux" function (see ComputeElementFlux() for a
147  description of the "flux" function). Typically, the energy of a "flux"
148  function should be equal to a_local(u,u), if the "flux" is defined from
149  a solution u; here a_local(.,.) denotes the element-local bilinear
150  form represented by the integrator.
151 
152  @param[in] fluxelem FiniteElement of the "flux".
153  @param[in] Trans The ElementTransformation describing the physical
154  position of the mesh element.
155  @param[in] flux "Flux" coefficients representing the expansion of the
156  "flux" function in the basis of @a fluxelem.
157  @param[out] d_energy If not NULL, the given Vector should be set to
158  represent directional energy split that can be used
159  for anisotropic error estimation.
160  @returns The computed energy.
161  */
162  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
164  Vector &flux, Vector *d_energy = NULL)
165  { return 0.0; }
166 
168 };
169 
171 {
172 private:
173  int own_bfi;
175 
176  DenseMatrix bfi_elmat;
177 
178 public:
179  TransposeIntegrator (BilinearFormIntegrator *_bfi, int _own_bfi = 1)
180  { bfi = _bfi; own_bfi = _own_bfi; }
181 
182  virtual void AssembleElementMatrix(const FiniteElement &el,
184  DenseMatrix &elmat);
185 
186  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
187  const FiniteElement &test_fe,
189  DenseMatrix &elmat);
190 
192  virtual void AssembleFaceMatrix(const FiniteElement &el1,
193  const FiniteElement &el2,
195  DenseMatrix &elmat);
196 
197  virtual ~TransposeIntegrator() { if (own_bfi) { delete bfi; } }
198 };
199 
201 {
202 private:
203  int own_bfi;
205 
206 public:
207  LumpedIntegrator (BilinearFormIntegrator *_bfi, int _own_bfi = 1)
208  { bfi = _bfi; own_bfi = _own_bfi; }
209 
210  virtual void AssembleElementMatrix(const FiniteElement &el,
212  DenseMatrix &elmat);
213 
214  virtual ~LumpedIntegrator() { if (own_bfi) { delete bfi; } }
215 };
216 
217 /// Integrator that inverts the matrix assembled by another integrator.
219 {
220 private:
221  int own_integrator;
222  BilinearFormIntegrator *integrator;
223 
224 public:
225  InverseIntegrator(BilinearFormIntegrator *integ, int own_integ = 1)
226  { integrator = integ; own_integrator = own_integ; }
227 
228  virtual void AssembleElementMatrix(const FiniteElement &el,
230  DenseMatrix &elmat);
231 
232  virtual ~InverseIntegrator() { if (own_integrator) { delete integrator; } }
233 };
234 
235 /// Integrator defining a sum of multiple Integrators.
237 {
238 private:
239  int own_integrators;
240  DenseMatrix elem_mat;
241  Array<BilinearFormIntegrator*> integrators;
242 
243 public:
244  SumIntegrator(int own_integs = 1) { own_integrators = own_integs; }
245 
247  { integrators.Append(integ); }
248 
249  virtual void AssembleElementMatrix(const FiniteElement &el,
251  DenseMatrix &elmat);
252 
253  virtual ~SumIntegrator();
254 };
255 
256 /** An abstract class for integrating the product of two scalar basis functions
257  with an optional scalar coefficient. */
259 {
260 public:
261 
262  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
263  const FiniteElement &test_fe,
265  DenseMatrix &elmat);
266 
267  /// Support for use in BilinearForm. Can be used only when appropriate.
268  virtual void AssembleElementMatrix(const FiniteElement &fe,
269  ElementTransformation &Trans,
270  DenseMatrix &elmat)
271  { AssembleElementMatrix2(fe, fe, Trans, elmat); }
272 
273 protected:
274  /// This parameter can be set by derived methods to enable single shape
275  /// evaluation in case CalcTestShape() and CalcTrialShape() return the same
276  /// result if given the same FiniteElement. The default is false.
278 
281 
282  inline virtual bool VerifyFiniteElementTypes(
283  const FiniteElement & trial_fe,
284  const FiniteElement & test_fe) const
285  {
286  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
288  }
289 
290  inline virtual const char * FiniteElementTypeFailureMessage() const
291  {
292  return "MixedScalarIntegrator: "
293  "Trial and test spaces must both be scalar fields.";
294  }
295 
296  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
297  const FiniteElement & test_fe,
298  ElementTransformation &Trans)
299  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
300 
301 
302  inline virtual void CalcTestShape(const FiniteElement & test_fe,
303  ElementTransformation &Trans,
304  Vector & shape)
305  { test_fe.CalcPhysShape(Trans, shape); }
306 
307  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
308  ElementTransformation &Trans,
309  Vector & shape)
310  { trial_fe.CalcPhysShape(Trans, shape); }
311 
313 
314 private:
315 
316 #ifndef MFEM_THREAD_SAFE
317  Vector test_shape;
318  Vector trial_shape;
319 #endif
320 
321 };
322 
323 /** An abstract class for integrating the inner product of two vector basis
324  functions with an optional scalar, vector, or matrix coefficient. */
326 {
327 public:
328 
329  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
330  const FiniteElement &test_fe,
332  DenseMatrix &elmat);
333 
334  /// Support for use in BilinearForm. Can be used only when appropriate.
335  virtual void AssembleElementMatrix(const FiniteElement &fe,
336  ElementTransformation &Trans,
337  DenseMatrix &elmat)
338  { AssembleElementMatrix2(fe, fe, Trans, elmat); }
339 
340 protected:
341  /// This parameter can be set by derived methods to enable single shape
342  /// evaluation in case CalcTestShape() and CalcTrialShape() return the same
343  /// result if given the same FiniteElement. The default is false.
345 
347  : same_calc_shape(false), Q(NULL), VQ(NULL), DQ(NULL), MQ(NULL) {}
349  : same_calc_shape(false), Q(&q), VQ(NULL), DQ(NULL), MQ(NULL) {}
351  : same_calc_shape(false), Q(NULL), VQ(diag?NULL:&dq), DQ(diag?&dq:NULL),
352  MQ(NULL) {}
354  : same_calc_shape(false), Q(NULL), VQ(NULL), DQ(NULL), MQ(&mq) {}
355 
356  inline virtual bool VerifyFiniteElementTypes(
357  const FiniteElement & trial_fe,
358  const FiniteElement & test_fe) const
359  {
360  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
362  }
363 
364  inline virtual const char * FiniteElementTypeFailureMessage() const
365  {
366  return "MixedVectorIntegrator: "
367  "Trial and test spaces must both be vector fields";
368  }
369 
370  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
371  const FiniteElement & test_fe,
372  ElementTransformation &Trans)
373  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
374 
375 
376  inline virtual void CalcTestShape(const FiniteElement & test_fe,
377  ElementTransformation &Trans,
378  DenseMatrix & shape)
379  { test_fe.CalcVShape(Trans, shape); }
380 
381  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
382  ElementTransformation &Trans,
383  DenseMatrix & shape)
384  { trial_fe.CalcVShape(Trans, shape); }
385 
390 
391 private:
392 
393 #ifndef MFEM_THREAD_SAFE
394  Vector V;
395  Vector D;
396  DenseMatrix M;
397  DenseMatrix test_shape;
398  DenseMatrix trial_shape;
399  DenseMatrix test_shape_tmp;
400 #endif
401 
402 };
403 
404 /** An abstract class for integrating the product of a scalar basis function and
405  the inner product of a vector basis function with a vector coefficient. In
406  2D the inner product can be replaced with a cross product. */
408 {
409 public:
410 
411  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
412  const FiniteElement &test_fe,
414  DenseMatrix &elmat);
415 
416 protected:
417 
418  MixedScalarVectorIntegrator(VectorCoefficient &vq, bool _transpose = false,
419  bool _cross_2d = false)
420  : VQ(&vq), transpose(_transpose), cross_2d(_cross_2d) {}
421 
422  inline virtual bool VerifyFiniteElementTypes(
423  const FiniteElement & trial_fe,
424  const FiniteElement & test_fe) const
425  {
426  return ((transpose &&
428  test_fe.GetRangeType() == mfem::FiniteElement::SCALAR ) ||
429  (!transpose &&
432  );
433  }
434 
435  inline virtual const char * FiniteElementTypeFailureMessage() const
436  {
437  if ( transpose )
438  {
439  return "MixedScalarVectorIntegrator: "
440  "Trial space must be a vector field "
441  "and the test space must be a scalar field";
442  }
443  else
444  {
445  return "MixedScalarVectorIntegrator: "
446  "Trial space must be a scalar field "
447  "and the test space must be a vector field";
448  }
449  }
450 
451  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
452  const FiniteElement & test_fe,
453  ElementTransformation &Trans)
454  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW(); }
455 
456 
457  inline virtual void CalcVShape(const FiniteElement & vector_fe,
458  ElementTransformation &Trans,
459  DenseMatrix & shape)
460  { vector_fe.CalcVShape(Trans, shape); }
461 
462  inline virtual void CalcShape(const FiniteElement & scalar_fe,
463  ElementTransformation &Trans,
464  Vector & shape)
465  { scalar_fe.CalcPhysShape(Trans, shape); }
466 
468  bool transpose;
469  bool cross_2d; // In 2D use a cross product rather than a dot product
470 
471 private:
472 
473 #ifndef MFEM_THREAD_SAFE
474  Vector V;
475  DenseMatrix vshape;
476  Vector shape;
477  Vector vshape_tmp;
478 #endif
479 
480 };
481 
482 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 1D, 2D,
483  or 3D and where Q is an optional scalar coefficient, u and v are each in H1
484  or L2. */
486 {
487 public:
490  : MixedScalarIntegrator(q) { same_calc_shape = true; }
491 };
492 
493 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 2D, or
494  3D and where Q is a vector coefficient, u is in H1 or L2 and v is in H(Curl)
495  or H(Div). */
497 {
498 public:
501 };
502 
503 /** Class for integrating the bilinear form a(u,v) := (Q D u, v) in 1D where Q
504  is an optional scalar coefficient, u is in H1, and v is in L2. */
506 {
507 public:
510  : MixedScalarIntegrator(q) {}
511 
512 protected:
513  inline virtual bool VerifyFiniteElementTypes(
514  const FiniteElement & trial_fe,
515  const FiniteElement & test_fe) const
516  {
517  return (trial_fe.GetDim() == 1 && test_fe.GetDim() == 1 &&
518  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
520  }
521 
522  inline virtual const char * FiniteElementTypeFailureMessage() const
523  {
524  return "MixedScalarDerivativeIntegrator: "
525  "Trial and test spaces must both be scalar fields in 1D "
526  "and the trial space must implement CaldDShape.";
527  }
528 
529  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
531  Vector & shape)
532  {
533  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
534  trial_fe.CalcPhysDShape(Trans, dshape);
535  }
536 };
537 
538 /** Class for integrating the bilinear form a(u,v) := -(Q u, D v) in 1D where Q
539  is an optional scalar coefficient, u is in L2, and v is in H1. */
541 {
542 public:
545  : MixedScalarIntegrator(q) {}
546 
547 protected:
548  inline virtual bool VerifyFiniteElementTypes(
549  const FiniteElement & trial_fe,
550  const FiniteElement & test_fe) const
551  {
552  return (trial_fe.GetDim() == 1 && test_fe.GetDim() == 1 &&
555  }
556 
557  inline virtual const char * FiniteElementTypeFailureMessage() const
558  {
559  return "MixedScalarWeakDerivativeIntegrator: "
560  "Trial and test spaces must both be scalar fields in 1D "
561  "and the test space must implement CalcDShape with "
562  "map type \"VALUE\".";
563  }
564 
565  inline virtual void CalcTestShape(const FiniteElement & test_fe,
567  Vector & shape)
568  {
569  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
570  test_fe.CalcPhysDShape(Trans, dshape);
571  shape *= -1.0;
572  }
573 };
574 
575 /** Class for integrating the bilinear form a(u,v) := (Q div u, v) in either 2D
576  or 3D where Q is an optional scalar coefficient, u is in H(Div), and v is a
577  scalar field. */
579 {
580 public:
583  : MixedScalarIntegrator(q) {}
584 
585 protected:
586  inline virtual bool VerifyFiniteElementTypes(
587  const FiniteElement & trial_fe,
588  const FiniteElement & test_fe) const
589  {
590  return (trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
592  }
593 
594  inline virtual const char * FiniteElementTypeFailureMessage() const
595  {
596  return "MixedScalarDivergenceIntegrator: "
597  "Trial must be H(Div) and the test space must be a "
598  "scalar field";
599  }
600 
601  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
602  const FiniteElement & test_fe,
604  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
605 
606  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
608  Vector & shape)
609  { trial_fe.CalcPhysDivShape(Trans, shape); }
610 };
611 
612 /** Class for integrating the bilinear form a(u,v) := (V div u, v) in either 2D
613  or 3D where V is a vector coefficient, u is in H(Div), and v is a vector
614  field. */
616 {
617 public:
620 
621 protected:
622  inline virtual bool VerifyFiniteElementTypes(
623  const FiniteElement & trial_fe,
624  const FiniteElement & test_fe) const
625  {
626  return (trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
628  }
629 
630  inline virtual const char * FiniteElementTypeFailureMessage() const
631  {
632  return "MixedVectorDivergenceIntegrator: "
633  "Trial must be H(Div) and the test space must be a "
634  "vector field";
635  }
636 
637  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
638  const FiniteElement & test_fe,
640  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
641 
642  inline virtual void CalcShape(const FiniteElement & scalar_fe,
644  Vector & shape)
645  { scalar_fe.CalcPhysDivShape(Trans, shape); }
646 };
647 
648 /** Class for integrating the bilinear form a(u,v) := -(Q u, div v) in either 2D
649  or 3D where Q is an optional scalar coefficient, u is in L2 or H1, and v is
650  in H(Div). */
652 {
653 public:
656  : MixedScalarIntegrator(q) {}
657 
658 protected:
659  inline virtual bool VerifyFiniteElementTypes(
660  const FiniteElement & trial_fe,
661  const FiniteElement & test_fe) const
662  {
663  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
664  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
665  }
666 
667  inline virtual const char * FiniteElementTypeFailureMessage() const
668  {
669  return "MixedScalarWeakGradientIntegrator: "
670  "Trial space must be a scalar field "
671  "and the test space must be H(Div)";
672  }
673 
674  virtual void CalcTestShape(const FiniteElement & test_fe,
676  Vector & shape)
677  {
678  test_fe.CalcPhysDivShape(Trans, shape);
679  shape *= -1.0;
680  }
681 };
682 
683 /** Class for integrating the bilinear form a(u,v) := (Q curl u, v) in 2D where
684  Q is an optional scalar coefficient, u is in H(Curl), and v is in L2 or
685  H1. */
687 {
688 public:
691  : MixedScalarIntegrator(q) {}
692 
693 protected:
694  inline virtual bool VerifyFiniteElementTypes(
695  const FiniteElement & trial_fe,
696  const FiniteElement & test_fe) const
697  {
698  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
699  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
701  }
702 
703  inline virtual const char * FiniteElementTypeFailureMessage() const
704  {
705  return "MixedScalarCurlIntegrator: "
706  "Trial must be H(Curl) and the test space must be a "
707  "scalar field";
708  }
709 
710  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
711  const FiniteElement & test_fe,
713  { return trial_fe.GetOrder() + test_fe.GetOrder() + Trans.OrderW() - 1; }
714 
715  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
717  Vector & shape)
718  {
719  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
720  trial_fe.CalcPhysCurlShape(Trans, dshape);
721  }
722 };
723 
724 /** Class for integrating the bilinear form a(u,v) := (Q u, curl v) in 2D where
725  Q is an optional scalar coefficient, u is in L2 or H1, and v is in
726  H(Curl). */
728 {
729 public:
732  : MixedScalarIntegrator(q) {}
733 
734 protected:
735  inline virtual bool VerifyFiniteElementTypes(
736  const FiniteElement & trial_fe,
737  const FiniteElement & test_fe) const
738  {
739  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
742  }
743 
744  inline virtual const char * FiniteElementTypeFailureMessage() const
745  {
746  return "MixedScalarWeakCurlIntegrator: "
747  "Trial space must be a scalar field "
748  "and the test space must be H(Curl)";
749  }
750 
751  inline virtual void CalcTestShape(const FiniteElement & test_fe,
753  Vector & shape)
754  {
755  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
756  test_fe.CalcPhysCurlShape(Trans, dshape);
757  }
758 };
759 
760 /** Class for integrating the bilinear form a(u,v) := (Q u, v) in either 2D or
761  3D and where Q is an optional coefficient (of type scalar, matrix, or
762  diagonal matrix) u and v are each in H(Curl) or H(Div). */
764 {
765 public:
768  : MixedVectorIntegrator(q) { same_calc_shape = true; }
770  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
772  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
773 };
774 
775 /** Class for integrating the bilinear form a(u,v) := (V x u, v) in 3D and where
776  V is a vector coefficient u and v are each in H(Curl) or H(Div). */
778 {
779 public:
781  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
782 };
783 
784 /** Class for integrating the bilinear form a(u,v) := (V . u, v) in 2D or 3D and
785  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in H1 or
786  L2. */
788 {
789 public:
791  : MixedScalarVectorIntegrator(vq, true) {}
792 
793  inline virtual bool VerifyFiniteElementTypes(
794  const FiniteElement & trial_fe,
795  const FiniteElement & test_fe) const
796  {
797  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
799  }
800 
801  inline virtual const char * FiniteElementTypeFailureMessage() const
802  {
803  return "MixedDotProductIntegrator: "
804  "Trial space must be a vector field "
805  "and the test space must be a scalar field";
806  }
807 };
808 
809 /** Class for integrating the bilinear form a(u,v) := (-V . u, Div v) in 2D or
810  3D and where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
811  RT. */
813 {
814 public:
816  : MixedScalarVectorIntegrator(vq, true) {}
817 
818  inline virtual bool VerifyFiniteElementTypes(
819  const FiniteElement & trial_fe,
820  const FiniteElement & test_fe) const
821  {
822  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
824  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
825  }
826 
827  inline virtual const char * FiniteElementTypeFailureMessage() const
828  {
829  return "MixedWeakGradDotIntegrator: "
830  "Trial space must be a vector field "
831  "and the test space must be a vector field with a divergence";
832  }
833 
834  inline virtual void CalcShape(const FiniteElement & scalar_fe,
836  Vector & shape)
837  { scalar_fe.CalcPhysDivShape(Trans, shape); shape *= -1.0; }
838 };
839 
840 /** Class for integrating the bilinear form a(u,v) := (V x u, Grad v) in 3D and
841  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in H1. */
843 {
844 public:
846  : MixedVectorIntegrator(vq, false) {}
847 
848  inline virtual bool VerifyFiniteElementTypes(
849  const FiniteElement & trial_fe,
850  const FiniteElement & test_fe) const
851  {
852  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
856  }
857 
858  inline virtual const char * FiniteElementTypeFailureMessage() const
859  {
860  return "MixedWeakDivCrossIntegrator: "
861  "Trial space must be a vector field in 3D "
862  "and the test space must be a scalar field with a gradient";
863  }
864 
865  inline virtual void CalcTestShape(const FiniteElement & test_fe,
867  DenseMatrix & shape)
868  { test_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
869 };
870 
871 /** Class for integrating the bilinear form a(u,v) := (Q Grad u, Grad v) in 3D
872  or in 2D and where Q is a scalar or matrix coefficient u and v are both in
873  H1. */
875 {
876 public:
879  : MixedVectorIntegrator(q) { same_calc_shape = true; }
881  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
883  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
884 
885  inline virtual bool VerifyFiniteElementTypes(
886  const FiniteElement & trial_fe,
887  const FiniteElement & test_fe) const
888  {
889  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
890  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
893  }
894 
895  inline virtual const char * FiniteElementTypeFailureMessage() const
896  {
897  return "MixedGradGradIntegrator: "
898  "Trial and test spaces must both be scalar fields "
899  "with a gradient operator.";
900  }
901 
902  inline virtual int GetIntegrationOrder(const FiniteElement & trial_fe,
903  const FiniteElement & test_fe,
905  {
906  // Same as DiffusionIntegrator
907  return test_fe.Space() == FunctionSpace::Pk ?
908  trial_fe.GetOrder() + test_fe.GetOrder() - 2 :
909  trial_fe.GetOrder() + test_fe.GetOrder() + test_fe.GetDim() - 1;
910  }
911 
912  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
914  DenseMatrix & shape)
915  { trial_fe.CalcPhysDShape(Trans, shape); }
916 
917  inline virtual void CalcTestShape(const FiniteElement & test_fe,
919  DenseMatrix & shape)
920  { test_fe.CalcPhysDShape(Trans, shape); }
921 };
922 
923 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, Grad v) in 3D
924  or in 2D and where V is a vector coefficient u and v are both in H1. */
926 {
927 public:
929  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
930 
931  inline virtual bool VerifyFiniteElementTypes(
932  const FiniteElement & trial_fe,
933  const FiniteElement & test_fe) const
934  {
935  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
936  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
939  }
940 
941  inline virtual const char * FiniteElementTypeFailureMessage() const
942  {
943  return "MixedCrossGradGradIntegrator: "
944  "Trial and test spaces must both be scalar fields "
945  "with a gradient operator.";
946  }
947 
948  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
950  DenseMatrix & shape)
951  { trial_fe.CalcPhysDShape(Trans, shape); }
952 
953  inline virtual void CalcTestShape(const FiniteElement & test_fe,
955  DenseMatrix & shape)
956  { test_fe.CalcPhysDShape(Trans, shape); }
957 };
958 
959 /** Class for integrating the bilinear form a(u,v) := (Q Curl u, Curl v) in 3D
960  and where Q is a scalar or matrix coefficient u and v are both in
961  H(Curl). */
963 {
964 public:
967  : MixedVectorIntegrator(q) { same_calc_shape = true; }
969  : MixedVectorIntegrator(dq, true) { same_calc_shape = true; }
971  : MixedVectorIntegrator(mq) { same_calc_shape = true; }
972 
973  inline virtual bool VerifyFiniteElementTypes(
974  const FiniteElement & trial_fe,
975  const FiniteElement & test_fe) const
976  {
977  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
979  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
982  }
983 
984  inline virtual const char * FiniteElementTypeFailureMessage() const
985  {
986  return "MixedCurlCurlIntegrator"
987  "Trial and test spaces must both be vector fields in 3D "
988  "with a curl.";
989  }
990 
991  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
993  DenseMatrix & shape)
994  { trial_fe.CalcPhysCurlShape(Trans, shape); }
995 
996  inline virtual void CalcTestShape(const FiniteElement & test_fe,
998  DenseMatrix & shape)
999  { test_fe.CalcPhysCurlShape(Trans, shape); }
1000 };
1001 
1002 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, Curl v) in 3D
1003  and where V is a vector coefficient u and v are both in H(Curl). */
1005 {
1006 public:
1008  : MixedVectorIntegrator(vq, false) { same_calc_shape = true; }
1009 
1010  inline virtual bool VerifyFiniteElementTypes(
1011  const FiniteElement & trial_fe,
1012  const FiniteElement & test_fe) const
1013  {
1014  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1015  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1016  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1018  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1019  }
1020 
1021  inline virtual const char * FiniteElementTypeFailureMessage() const
1022  {
1023  return "MixedCrossCurlCurlIntegrator: "
1024  "Trial and test spaces must both be vector fields in 3D "
1025  "with a curl.";
1026  }
1027 
1028  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1030  DenseMatrix & shape)
1031  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1032 
1033  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1035  DenseMatrix & shape)
1036  { test_fe.CalcPhysCurlShape(Trans, shape); }
1037 };
1038 
1039 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, Grad v) in 3D
1040  and where V is a vector coefficient u is in H(Curl) and v is in H1. */
1042 {
1043 public:
1045  : MixedVectorIntegrator(vq, false) {}
1046 
1047  inline virtual bool VerifyFiniteElementTypes(
1048  const FiniteElement & trial_fe,
1049  const FiniteElement & test_fe) const
1050  {
1051  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1052  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1053  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1055  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1056  }
1057 
1058  inline virtual const char * FiniteElementTypeFailureMessage() const
1059  {
1060  return "MixedCrossCurlGradIntegrator"
1061  "Trial space must be a vector field in 3D with a curl"
1062  "and the test space must be a scalar field with a gradient";
1063  }
1064 
1065  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1067  DenseMatrix & shape)
1068  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1069 
1070  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1072  DenseMatrix & shape)
1073  { test_fe.CalcPhysDShape(Trans, shape); }
1074 };
1075 
1076 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, Curl v) in 3D
1077  and where V is a scalar coefficient u is in H1 and v is in H(Curl). */
1079 {
1080 public:
1082  : MixedVectorIntegrator(vq, false) {}
1083 
1084  inline virtual bool VerifyFiniteElementTypes(
1085  const FiniteElement & trial_fe,
1086  const FiniteElement & test_fe) const
1087  {
1088  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1089  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1090  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1092  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1093  }
1094 
1095  inline virtual const char * FiniteElementTypeFailureMessage() const
1096  {
1097  return "MixedCrossGradCurlIntegrator"
1098  "Trial space must be a scalar field in 3D with a gradient"
1099  "and the test space must be a vector field with a curl";
1100  }
1101 
1102  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1104  DenseMatrix & shape)
1105  { trial_fe.CalcPhysDShape(Trans, shape); }
1106 
1107  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1109  DenseMatrix & shape)
1110  { test_fe.CalcPhysCurlShape(Trans, shape); }
1111 };
1112 
1113 /** Class for integrating the bilinear form a(u,v) := (V x u, Curl v) in 3D and
1114  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
1115  H(Curl). */
1117 {
1118 public:
1120  : MixedVectorIntegrator(vq, false) {}
1121 
1122  inline virtual bool VerifyFiniteElementTypes(
1123  const FiniteElement & trial_fe,
1124  const FiniteElement & test_fe) const
1125  {
1126  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1127  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1129  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1130  }
1131 
1132  inline virtual const char * FiniteElementTypeFailureMessage() const
1133  {
1134  return "MixedWeakCurlCrossIntegrator: "
1135  "Trial space must be a vector field in 3D "
1136  "and the test space must be a vector field with a curl";
1137  }
1138 
1139  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1141  DenseMatrix & shape)
1142  { test_fe.CalcPhysCurlShape(Trans, shape); }
1143 };
1144 
1145 /** Class for integrating the bilinear form a(u,v) := (V x u, Curl v) in 2D and
1146  where V is a vector coefficient u is in H(Curl) or H(Div) and v is in
1147  H(Curl). */
1149 {
1150 public:
1152  : MixedScalarVectorIntegrator(vq, true, true) {}
1153 
1154  inline virtual bool VerifyFiniteElementTypes(
1155  const FiniteElement & trial_fe,
1156  const FiniteElement & test_fe) const
1157  {
1158  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1159  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1161  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1162  }
1163 
1164  inline virtual const char * FiniteElementTypeFailureMessage() const
1165  {
1166  return "MixedScalarWeakCurlCrossIntegrator: "
1167  "Trial space must be a vector field in 2D "
1168  "and the test space must be a vector field with a curl";
1169  }
1170 
1171  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1173  Vector & shape)
1174  {
1175  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
1176  scalar_fe.CalcPhysCurlShape(Trans, dshape);
1177  }
1178 };
1179 
1180 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, v) in 3D or
1181  in 2D and where V is a vector coefficient u is in H1 and v is in H(Curl) or
1182  H(Div). */
1184 {
1185 public:
1187  : MixedVectorIntegrator(vq, false) {}
1188 
1189  inline virtual bool VerifyFiniteElementTypes(
1190  const FiniteElement & trial_fe,
1191  const FiniteElement & test_fe) const
1192  {
1193  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1194  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1195  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1197  }
1198 
1199  inline virtual const char * FiniteElementTypeFailureMessage() const
1200  {
1201  return "MixedCrossGradIntegrator: "
1202  "Trial space must be a scalar field with a gradient operator"
1203  " and the test space must be a vector field both in 3D.";
1204  }
1205 
1206  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1208  DenseMatrix & shape)
1209  { trial_fe.CalcPhysDShape(Trans, shape); }
1210 
1211  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1213  DenseMatrix & shape)
1214  { test_fe.CalcVShape(Trans, shape); }
1215 };
1216 
1217 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, v) in 3D and
1218  where V is a vector coefficient u is in H(Curl) and v is in H(Curl) or
1219  H(Div). */
1221 {
1222 public:
1224  : MixedVectorIntegrator(vq, false) {}
1225 
1226  inline virtual bool VerifyFiniteElementTypes(
1227  const FiniteElement & trial_fe,
1228  const FiniteElement & test_fe) const
1229  {
1230  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1231  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1232  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1234  }
1235 
1236  inline virtual const char * FiniteElementTypeFailureMessage() const
1237  {
1238  return "MixedCrossCurlIntegrator: "
1239  "Trial space must be a vector field in 3D with a curl "
1240  "and the test space must be a vector field";
1241  }
1242 
1243  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1245  DenseMatrix & shape)
1246  { trial_fe.CalcPhysCurlShape(Trans, shape); }
1247 };
1248 
1249 /** Class for integrating the bilinear form a(u,v) := (V x Curl u, v) in 2D and
1250  where V is a vector coefficient u is in H(Curl) and v is in H(Curl) or
1251  H(Div). */
1253 {
1254 public:
1256  : MixedScalarVectorIntegrator(vq, false, true) {}
1257 
1258  inline virtual bool VerifyFiniteElementTypes(
1259  const FiniteElement & trial_fe,
1260  const FiniteElement & test_fe) const
1261  {
1262  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1263  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1264  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1266  }
1267 
1268  inline virtual const char * FiniteElementTypeFailureMessage() const
1269  {
1270  return "MixedCrossCurlIntegrator: "
1271  "Trial space must be a vector field in 2D with a curl "
1272  "and the test space must be a vector field";
1273  }
1274 
1275  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1277  Vector & shape)
1278  {
1279  DenseMatrix dshape(shape.GetData(), shape.Size(), 1);
1280  scalar_fe.CalcPhysCurlShape(Trans, dshape); shape *= -1.0;
1281  }
1282 };
1283 
1284 /** Class for integrating the bilinear form a(u,v) := (V x Grad u, v) in 2D and
1285  where V is a vector coefficient u is in H1 and v is in H1 or L2. */
1287 {
1288 public:
1290  : MixedScalarVectorIntegrator(vq, true, true) {}
1291 
1292  inline virtual bool VerifyFiniteElementTypes(
1293  const FiniteElement & trial_fe,
1294  const FiniteElement & test_fe) const
1295  {
1296  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1297  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1298  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1300  }
1301 
1302  inline virtual const char * FiniteElementTypeFailureMessage() const
1303  {
1304  return "MixedScalarCrossGradIntegrator: "
1305  "Trial space must be a scalar field in 2D with a gradient "
1306  "and the test space must be a scalar field";
1307  }
1308 
1309  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1311  DenseMatrix & shape)
1312  { vector_fe.CalcPhysDShape(Trans, shape); }
1313 };
1314 
1315 /** Class for integrating the bilinear form a(u,v) := (V x u, v) in 2D and where
1316  V is a vector coefficient u is in ND or RT and v is in H1 or L2. */
1318 {
1319 public:
1321  : MixedScalarVectorIntegrator(vq, true, true) {}
1322 
1323  inline virtual bool VerifyFiniteElementTypes(
1324  const FiniteElement & trial_fe,
1325  const FiniteElement & test_fe) const
1326  {
1327  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1328  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1330  }
1331 
1332  inline virtual const char * FiniteElementTypeFailureMessage() const
1333  {
1334  return "MixedScalarCrossProductIntegrator: "
1335  "Trial space must be a vector field in 2D "
1336  "and the test space must be a scalar field";
1337  }
1338 };
1339 
1340 /** Class for integrating the bilinear form a(u,v) := (V x z u, v) in 2D and
1341  where V is a vector coefficient u is in H1 or L2 and v is in ND or RT. */
1343 {
1344 public:
1346  : MixedScalarVectorIntegrator(vq, false, true) {}
1347 
1348  inline virtual bool VerifyFiniteElementTypes(
1349  const FiniteElement & trial_fe,
1350  const FiniteElement & test_fe) const
1351  {
1352  return (trial_fe.GetDim() == 2 && test_fe.GetDim() == 2 &&
1353  trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1355  }
1356 
1357  inline virtual const char * FiniteElementTypeFailureMessage() const
1358  {
1359  return "MixedScalarWeakCrossProductIntegrator: "
1360  "Trial space must be a scalar field in 2D "
1361  "and the test space must be a vector field";
1362  }
1363 
1364  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1366  Vector & shape)
1367  { scalar_fe.CalcPhysShape(Trans, shape); shape *= -1.0; }
1368 };
1369 
1370 /** Class for integrating the bilinear form a(u,v) := (V . Grad u, v) in 2D or
1371  3D and where V is a vector coefficient, u is in H1 and v is in H1 or L2. */
1373 {
1374 public:
1376  : MixedScalarVectorIntegrator(vq, true) {}
1377 
1378  inline virtual bool VerifyFiniteElementTypes(
1379  const FiniteElement & trial_fe,
1380  const FiniteElement & test_fe) const
1381  {
1382  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1383  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1385  }
1386 
1387  inline virtual const char * FiniteElementTypeFailureMessage() const
1388  {
1389  return "MixedDirectionalDerivativeIntegrator: "
1390  "Trial space must be a scalar field with a gradient "
1391  "and the test space must be a scalar field";
1392  }
1393 
1394  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1396  DenseMatrix & shape)
1397  { vector_fe.CalcPhysDShape(Trans, shape); }
1398 };
1399 
1400 /** Class for integrating the bilinear form a(u,v) := (-V . Grad u, Div v) in 2D
1401  or 3D and where V is a vector coefficient, u is in H1 and v is in RT. */
1403 {
1404 public:
1406  : MixedScalarVectorIntegrator(vq, true) {}
1407 
1408  inline virtual bool VerifyFiniteElementTypes(
1409  const FiniteElement & trial_fe,
1410  const FiniteElement & test_fe) const
1411  {
1412  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1413  trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1415  test_fe.GetDerivType() == mfem::FiniteElement::DIV );
1416  }
1417 
1418  inline virtual const char * FiniteElementTypeFailureMessage() const
1419  {
1420  return "MixedGradDivIntegrator: "
1421  "Trial space must be a scalar field with a gradient"
1422  "and the test space must be a vector field with a divergence";
1423  }
1424 
1425  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1427  DenseMatrix & shape)
1428  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1429 
1430  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1432  Vector & shape)
1433  { scalar_fe.CalcPhysDivShape(Trans, shape); }
1434 };
1435 
1436 /** Class for integrating the bilinear form a(u,v) := (-V Div u, Grad v) in 2D
1437  or 3D and where V is a vector coefficient, u is in RT and v is in H1. */
1439 {
1440 public:
1442  : MixedScalarVectorIntegrator(vq, false) {}
1443 
1444  inline virtual bool VerifyFiniteElementTypes(
1445  const FiniteElement & trial_fe,
1446  const FiniteElement & test_fe) const
1447  {
1448  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1449  trial_fe.GetDerivType() == mfem::FiniteElement::DIV &&
1452  );
1453  }
1454 
1455  inline virtual const char * FiniteElementTypeFailureMessage() const
1456  {
1457  return "MixedDivGradIntegrator: "
1458  "Trial space must be a vector field with a divergence"
1459  "and the test space must be a scalar field with a gradient";
1460  }
1461 
1462  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1464  DenseMatrix & shape)
1465  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1466 
1467  inline virtual void CalcShape(const FiniteElement & scalar_fe,
1469  Vector & shape)
1470  { scalar_fe.CalcPhysDivShape(Trans, shape); }
1471 };
1472 
1473 /** Class for integrating the bilinear form a(u,v) := (-V u, Grad v) in 2D or 3D
1474  and where V is a vector coefficient, u is in H1 and v is in H1. */
1476 {
1477 public:
1479  : MixedScalarVectorIntegrator(vq, false) {}
1480 
1481  inline virtual bool VerifyFiniteElementTypes(
1482  const FiniteElement & trial_fe,
1483  const FiniteElement & test_fe) const
1484  {
1485  return (trial_fe.GetRangeType() == mfem::FiniteElement::SCALAR &&
1487  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1488  }
1489 
1490  inline virtual const char * FiniteElementTypeFailureMessage() const
1491  {
1492  return "MixedScalarWeakDivergenceIntegrator: "
1493  "Trial space must be a scalar field "
1494  "and the test space must be a scalar field with a gradient";
1495  }
1496 
1497  inline virtual void CalcVShape(const FiniteElement & vector_fe,
1499  DenseMatrix & shape)
1500  { vector_fe.CalcPhysDShape(Trans, shape); shape *= -1.0; }
1501 };
1502 
1503 /** Class for integrating the bilinear form a(u,v) := (Q grad u, v) in either 2D
1504  or 3D and where Q is an optional coefficient (of type scalar, matrix, or
1505  diagonal matrix) u is in H1 and v is in H(Curl) or H(Div). */
1507 {
1508 public:
1511  : MixedVectorIntegrator(q) {}
1513  : MixedVectorIntegrator(dq, true) {}
1515  : MixedVectorIntegrator(mq) {}
1516 
1517 protected:
1518  inline virtual bool VerifyFiniteElementTypes(
1519  const FiniteElement & trial_fe,
1520  const FiniteElement & test_fe) const
1521  {
1522  return (trial_fe.GetDerivType() == mfem::FiniteElement::GRAD &&
1524  }
1525 
1526  inline virtual const char * FiniteElementTypeFailureMessage() const
1527  {
1528  return "MixedVectorGradientIntegrator: "
1529  "Trial spaces must be H1 and the test space must be a "
1530  "vector field in 2D or 3D";
1531  }
1532 
1533  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1535  DenseMatrix & shape)
1536  {
1537  trial_fe.CalcPhysDShape(Trans, shape);
1538  }
1539 
1540 private:
1541  DenseMatrix Jinv;
1542 };
1543 
1544 /** Class for integrating the bilinear form a(u,v) := (Q curl u, v) in 3D and
1545  where Q is an optional coefficient (of type scalar, matrix, or diagonal
1546  matrix) u is in H(Curl) and v is in H(Div) or H(Curl). */
1548 {
1549 public:
1552  : MixedVectorIntegrator(q) {}
1554  : MixedVectorIntegrator(dq, true) {}
1556  : MixedVectorIntegrator(mq) {}
1557 
1558 protected:
1559  inline virtual bool VerifyFiniteElementTypes(
1560  const FiniteElement & trial_fe,
1561  const FiniteElement & test_fe) const
1562  {
1563  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1564  trial_fe.GetDerivType() == mfem::FiniteElement::CURL &&
1566  }
1567 
1568  inline virtual const char * FiniteElementTypeFailureMessage() const
1569  {
1570  return "MixedVectorCurlIntegrator: "
1571  "Trial space must be H(Curl) and the test space must be a "
1572  "vector field in 3D";
1573  }
1574 
1575  inline virtual void CalcTrialShape(const FiniteElement & trial_fe,
1577  DenseMatrix & shape)
1578  {
1579  trial_fe.CalcPhysCurlShape(Trans, shape);
1580  }
1581 };
1582 
1583 /** Class for integrating the bilinear form a(u,v) := (Q u, curl v) in 3D and
1584  where Q is an optional coefficient (of type scalar, matrix, or diagonal
1585  matrix) u is in H(Div) or H(Curl) and v is in H(Curl). */
1587 {
1588 public:
1591  : MixedVectorIntegrator(q) {}
1593  : MixedVectorIntegrator(dq, true) {}
1595  : MixedVectorIntegrator(mq) {}
1596 
1597 protected:
1598  inline virtual bool VerifyFiniteElementTypes(
1599  const FiniteElement & trial_fe,
1600  const FiniteElement & test_fe) const
1601  {
1602  return (trial_fe.GetDim() == 3 && test_fe.GetDim() == 3 &&
1603  trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1604  test_fe.GetDerivType() == mfem::FiniteElement::CURL );
1605  }
1606 
1607  inline virtual const char * FiniteElementTypeFailureMessage() const
1608  {
1609  return "MixedVectorWeakCurlIntegrator: "
1610  "Trial space must be vector field in 3D and the "
1611  "test space must be H(Curl)";
1612  }
1613 
1614  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1616  DenseMatrix & shape)
1617  {
1618  test_fe.CalcPhysCurlShape(Trans, shape);
1619  }
1620 };
1621 
1622 /** Class for integrating the bilinear form a(u,v) := - (Q u, grad v) in either
1623  2D or 3D and where Q is an optional coefficient (of type scalar, matrix, or
1624  diagonal matrix) u is in H(Div) or H(Curl) and v is in H1. */
1626 {
1627 public:
1630  : MixedVectorIntegrator(q) {}
1632  : MixedVectorIntegrator(dq, true) {}
1634  : MixedVectorIntegrator(mq) {}
1635 
1636 protected:
1637  inline virtual bool VerifyFiniteElementTypes(
1638  const FiniteElement & trial_fe,
1639  const FiniteElement & test_fe) const
1640  {
1641  return (trial_fe.GetRangeType() == mfem::FiniteElement::VECTOR &&
1642  test_fe.GetDerivType() == mfem::FiniteElement::GRAD );
1643  }
1644 
1645  inline virtual const char * FiniteElementTypeFailureMessage() const
1646  {
1647  return "MixedVectorWeakDivergenceIntegrator: "
1648  "Trial space must be vector field and the "
1649  "test space must be H1";
1650  }
1651 
1652  inline virtual void CalcTestShape(const FiniteElement & test_fe,
1654  DenseMatrix & shape)
1655  {
1656  test_fe.CalcPhysDShape(Trans, shape);
1657  shape *= -1.0;
1658  }
1659 };
1660 
1661 /** Class for integrating the bilinear form a(u,v) := (Q grad u, grad v) where Q
1662  can be a scalar or a matrix coefficient. */
1664 {
1665 protected:
1668 
1669 private:
1670  Vector vec, pointflux, shape;
1671 #ifndef MFEM_THREAD_SAFE
1672  DenseMatrix dshape, dshapedxt, invdfdx, mq;
1673  DenseMatrix te_dshape, te_dshapedxt;
1674 #endif
1675 
1676  // PA extension
1677  const DofToQuad *maps; ///< Not owned
1678  const GeometricFactors *geom; ///< Not owned
1679  int dim, ne, dofs1D, quad1D;
1680  Vector pa_data;
1681 
1682 public:
1683  /// Construct a diffusion integrator with coefficient Q = 1
1684  DiffusionIntegrator() { Q = NULL; MQ = NULL; maps = NULL; geom = NULL; }
1685 
1686  /// Construct a diffusion integrator with a scalar coefficient q
1688  : Q(&q) { MQ = NULL; maps = NULL; geom = NULL; }
1689 
1690  /// Construct a diffusion integrator with a matrix coefficient q
1692  : MQ(&q) { Q = NULL; maps = NULL; geom = NULL; }
1693 
1694  /** Given a particular Finite Element
1695  computes the element stiffness matrix elmat. */
1696  virtual void AssembleElementMatrix(const FiniteElement &el,
1698  DenseMatrix &elmat);
1699  /** Given a trial and test Finite Element computes the element stiffness
1700  matrix elmat. */
1701  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1702  const FiniteElement &test_fe,
1704  DenseMatrix &elmat);
1705 
1706  /// Perform the local action of the BilinearFormIntegrator
1707  virtual void AssembleElementVector(const FiniteElement &el,
1709  const Vector &elfun, Vector &elvect);
1710 
1711  virtual void ComputeElementFlux(const FiniteElement &el,
1713  Vector &u, const FiniteElement &fluxelem,
1714  Vector &flux, int with_coef = 1);
1715 
1716  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
1718  Vector &flux, Vector *d_energy = NULL);
1719 
1720  virtual void AssemblePA(const FiniteElementSpace&);
1721 
1722  virtual void AddMultPA(const Vector&, Vector&) const;
1723 
1724  static const IntegrationRule &GetRule(const FiniteElement &trial_fe,
1725  const FiniteElement &test_fe);
1726 };
1727 
1728 /** Class for local mass matrix assembling a(u,v) := (Q u, v) */
1730 {
1731 protected:
1732 #ifndef MFEM_THREAD_SAFE
1734 #endif
1736  // PA extension
1738  const DofToQuad *maps; ///< Not owned
1739  const GeometricFactors *geom; ///< Not owned
1740  int dim, ne, nq, dofs1D, quad1D;
1741 
1742 public:
1744  : BilinearFormIntegrator(ir) { Q = NULL; maps = NULL; geom = NULL; }
1745 
1746  /// Construct a mass integrator with coefficient q
1748  : BilinearFormIntegrator(ir), Q(&q) { maps = NULL; geom = NULL; }
1749 
1750  /** Given a particular Finite Element
1751  computes the element mass matrix elmat. */
1752  virtual void AssembleElementMatrix(const FiniteElement &el,
1754  DenseMatrix &elmat);
1755  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1756  const FiniteElement &test_fe,
1758  DenseMatrix &elmat);
1759 
1760  virtual void AssemblePA(const FiniteElementSpace&);
1761 
1762  virtual void AddMultPA(const Vector&, Vector&) const;
1763 
1764  static const IntegrationRule &GetRule(const FiniteElement &trial_fe,
1765  const FiniteElement &test_fe,
1767 };
1768 
1770 {
1771 public:
1773 
1775 
1776  virtual void AssembleFaceMatrix(const FiniteElement &el1,
1777  const FiniteElement &el2,
1779  DenseMatrix &elmat);
1780 };
1781 
1782 /// alpha (q . grad u, v)
1784 {
1785 protected:
1787  double alpha;
1788 
1789 private:
1790 #ifndef MFEM_THREAD_SAFE
1791  DenseMatrix dshape, adjJ, Q_ir;
1792  Vector shape, vec2, BdFidxT;
1793 #endif
1794 
1795 public:
1797  : Q(&q) { alpha = a; }
1798  virtual void AssembleElementMatrix(const FiniteElement &,
1800  DenseMatrix &);
1801 };
1802 
1803 /// alpha (q . grad u, v) using the "group" FE discretization
1805 {
1806 protected:
1808  double alpha;
1809 
1810 private:
1811  DenseMatrix dshape, adjJ, Q_nodal, grad;
1812  Vector shape;
1813 
1814 public:
1816  : Q(&q) { alpha = a; }
1817  virtual void AssembleElementMatrix(const FiniteElement &,
1819  DenseMatrix &);
1820 };
1821 
1822 /** Class for integrating the bilinear form a(u,v) := (Q u, v),
1823  where u=(u1,...,un) and v=(v1,...,vn); ui and vi are defined
1824  by scalar FE through standard transformation. */
1826 {
1827 private:
1828  int vdim;
1829  Vector shape, te_shape, vec;
1830  DenseMatrix partelmat;
1831  DenseMatrix mcoeff;
1832  int Q_order;
1833 
1834 protected:
1838 
1839 public:
1840  /// Construct an integrator with coefficient 1.0
1842  : vdim(-1), Q_order(0), Q(NULL), VQ(NULL), MQ(NULL) { }
1843  /** Construct an integrator with scalar coefficient q.
1844  If possible, save memory by using a scalar integrator since
1845  the resulting matrix is block diagonal with the same diagonal
1846  block repeated. */
1848  : vdim(-1), Q(&q) { VQ = NULL; MQ = NULL; Q_order = qo; }
1850  : BilinearFormIntegrator(ir), vdim(-1), Q(&q)
1851  { VQ = NULL; MQ = NULL; Q_order = 0; }
1852  /// Construct an integrator with diagonal coefficient q
1854  : vdim(q.GetVDim()), VQ(&q) { Q = NULL; MQ = NULL; Q_order = qo; }
1855  /// Construct an integrator with matrix coefficient q
1857  : vdim(q.GetVDim()), MQ(&q) { Q = NULL; VQ = NULL; Q_order = qo; }
1858 
1859  int GetVDim() const { return vdim; }
1860  void SetVDim(int vdim) { this->vdim = vdim; }
1861 
1862  virtual void AssembleElementMatrix(const FiniteElement &el,
1864  DenseMatrix &elmat);
1865  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1866  const FiniteElement &test_fe,
1868  DenseMatrix &elmat);
1869 };
1870 
1871 
1872 /** Class for integrating (div u, p) where u is a vector field given
1873  by VectorFiniteElement through Piola transformation (for RT
1874  elements); p is scalar function given by FiniteElement through
1875  standard transformation. Here, u is the trial function and p is
1876  the test function.
1877  Note: the element matrix returned by AssembleElementMatrix2
1878  does NOT depend on the ElementTransformation Trans. */
1880 {
1881 protected:
1883 
1884 private:
1885 #ifndef MFEM_THREAD_SAFE
1886  Vector divshape, shape;
1887 #endif
1888 
1889 public:
1892  virtual void AssembleElementMatrix(const FiniteElement &el,
1894  DenseMatrix &elmat) { }
1895  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1896  const FiniteElement &test_fe,
1898  DenseMatrix &elmat);
1899 };
1900 
1901 
1902 /** Integrator for `(-Q u, grad v)` for Nedelec (`u`) and H1 (`v`) elements.
1903  This is equivalent to a weak divergence of the Nedelec basis functions. */
1905 {
1906 protected:
1908 
1909 private:
1910 #ifndef MFEM_THREAD_SAFE
1911  DenseMatrix dshape;
1912  DenseMatrix dshapedxt;
1913  DenseMatrix vshape;
1914  DenseMatrix invdfdx;
1915 #endif
1916 
1917 public:
1920  virtual void AssembleElementMatrix(const FiniteElement &el,
1922  DenseMatrix &elmat) { }
1923  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1924  const FiniteElement &test_fe,
1926  DenseMatrix &elmat);
1927 };
1928 
1929 /** Integrator for (curl u, v) for Nedelec and RT elements. If the trial and
1930  test spaces are switched, assembles the form (u, curl v). */
1932 {
1933 protected:
1935 
1936 private:
1937 #ifndef MFEM_THREAD_SAFE
1938  DenseMatrix curlshapeTrial;
1939  DenseMatrix vshapeTest;
1940  DenseMatrix curlshapeTrial_dFT;
1941 #endif
1942 
1943 public:
1946  virtual void AssembleElementMatrix(const FiniteElement &el,
1948  DenseMatrix &elmat) { }
1949  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1950  const FiniteElement &test_fe,
1952  DenseMatrix &elmat);
1953 };
1954 
1955 /// Class for integrating (Q D_i(u), v); u and v are scalars
1957 {
1958 protected:
1960 
1961 private:
1962  int xi;
1963  DenseMatrix dshape, dshapedxt, invdfdx;
1964  Vector shape, dshapedxi;
1965 
1966 public:
1967  DerivativeIntegrator(Coefficient &q, int i) : Q(&q), xi(i) { }
1968  virtual void AssembleElementMatrix(const FiniteElement &el,
1970  DenseMatrix &elmat)
1971  { AssembleElementMatrix2(el,el,Trans,elmat); }
1972  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
1973  const FiniteElement &test_fe,
1975  DenseMatrix &elmat);
1976 };
1977 
1978 /// Integrator for (curl u, curl v) for Nedelec elements
1980 {
1981 private:
1982  Vector vec, pointflux;
1983 #ifndef MFEM_THREAD_SAFE
1984  DenseMatrix curlshape, curlshape_dFt, M;
1985  DenseMatrix vshape, projcurl;
1986 #endif
1987 
1988 protected:
1991 
1992 public:
1993  CurlCurlIntegrator() { Q = NULL; MQ = NULL; }
1994  /// Construct a bilinear form integrator for Nedelec elements
1995  CurlCurlIntegrator(Coefficient &q) : Q(&q) { MQ = NULL; }
1997 
1998  /* Given a particular Finite Element, compute the
1999  element curl-curl matrix elmat */
2000  virtual void AssembleElementMatrix(const FiniteElement &el,
2002  DenseMatrix &elmat);
2003 
2004  virtual void ComputeElementFlux(const FiniteElement &el,
2006  Vector &u, const FiniteElement &fluxelem,
2007  Vector &flux, int with_coef);
2008 
2009  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
2011  Vector &flux, Vector *d_energy = NULL);
2012 };
2013 
2014 /** Integrator for (curl u, curl v) for FE spaces defined by 'dim' copies of a
2015  scalar FE space. */
2017 {
2018 private:
2019 #ifndef MFEM_THREAD_SAFE
2020  DenseMatrix dshape_hat, dshape, curlshape, Jadj, grad_hat, grad;
2021 #endif
2022 
2023 protected:
2025 
2026 public:
2028 
2030 
2031  /// Assemble an element matrix
2032  virtual void AssembleElementMatrix(const FiniteElement &el,
2034  DenseMatrix &elmat);
2035  /// Compute element energy: (1/2) (curl u, curl u)_E
2036  virtual double GetElementEnergy(const FiniteElement &el,
2038  const Vector &elfun);
2039 };
2040 
2041 /// Integrator for (Q u, v) for VectorFiniteElements
2043 {
2044 private:
2046  { Q = q; VQ = vq; MQ = mq; }
2047 
2048 #ifndef MFEM_THREAD_SAFE
2049  Vector shape;
2050  Vector D;
2051  DenseMatrix K;
2052  DenseMatrix test_vshape;
2053  DenseMatrix trial_vshape;
2054 #endif
2055 
2056 protected:
2060 
2061 public:
2062  VectorFEMassIntegrator() { Init(NULL, NULL, NULL); }
2063  VectorFEMassIntegrator(Coefficient *_q) { Init(_q, NULL, NULL); }
2064  VectorFEMassIntegrator(Coefficient &q) { Init(&q, NULL, NULL); }
2065  VectorFEMassIntegrator(VectorCoefficient *_vq) { Init(NULL, _vq, NULL); }
2066  VectorFEMassIntegrator(VectorCoefficient &vq) { Init(NULL, &vq, NULL); }
2067  VectorFEMassIntegrator(MatrixCoefficient *_mq) { Init(NULL, NULL, _mq); }
2068  VectorFEMassIntegrator(MatrixCoefficient &mq) { Init(NULL, NULL, &mq); }
2069 
2070  virtual void AssembleElementMatrix(const FiniteElement &el,
2072  DenseMatrix &elmat);
2073  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2074  const FiniteElement &test_fe,
2076  DenseMatrix &elmat);
2077 };
2078 
2079 /** Integrator for (Q div u, p) where u=(v1,...,vn) and all vi are in the same
2080  scalar FE space; p is also in a (different) scalar FE space. */
2082 {
2083 protected:
2085 
2086 private:
2087  Vector shape;
2088  Vector divshape;
2089  DenseMatrix dshape;
2090  DenseMatrix gshape;
2091  DenseMatrix Jadj;
2092 
2093 public:
2097 
2098  virtual void AssembleElementMatrix2(const FiniteElement &trial_fe,
2099  const FiniteElement &test_fe,
2101  DenseMatrix &elmat);
2102 };
2103 
2104 /// (Q div u, div v) for RT elements
2106 {
2107 protected:
2109 
2110 private:
2111 #ifndef MFEM_THREAD_SAFE
2112  Vector divshape;
2113 #endif
2114 
2115 public:
2116  DivDivIntegrator() { Q = NULL; }
2118 
2119  virtual void AssembleElementMatrix(const FiniteElement &el,
2121  DenseMatrix &elmat);
2122 };
2123 
2124 /** Integrator for
2125  (Q grad u, grad v) = sum_i (Q grad u_i, grad v_i) e_i e_i^T
2126  for FE spaces defined by 'dim' copies of a scalar FE space. Where e_i
2127  is the unit vector in the i-th direction. The resulting local element
2128  matrix is a block-diagonal matrix consisting of 'dim' copies of a scalar
2129  diffusion matrix in each diagonal block. */
2131 {
2132 protected:
2134 
2135 private:
2136  DenseMatrix Jinv;
2137  DenseMatrix dshape;
2138  DenseMatrix gshape;
2139  DenseMatrix pelmat;
2140 
2141 public:
2144 
2145  virtual void AssembleElementMatrix(const FiniteElement &el,
2147  DenseMatrix &elmat);
2148  virtual void AssembleElementVector(const FiniteElement &el,
2150  const Vector &elfun, Vector &elvect);
2151 };
2152 
2153 /** Integrator for the linear elasticity form:
2154  a(u,v) = (lambda div(u), div(v)) + (2 mu e(u), e(v)),
2155  where e(v) = (1/2) (grad(v) + grad(v)^T).
2156  This is a 'Vector' integrator, i.e. defined for FE spaces
2157  using multiple copies of a scalar FE space. */
2159 {
2160 protected:
2161  double q_lambda, q_mu;
2163 
2164 private:
2165 #ifndef MFEM_THREAD_SAFE
2166  Vector shape;
2167  DenseMatrix dshape, gshape, pelmat;
2168  Vector divshape;
2169 #endif
2170 
2171 public:
2173  { lambda = &l; mu = &m; }
2174  /** With this constructor lambda = q_l * m and mu = q_m * m;
2175  if dim * q_l + 2 * q_m = 0 then trace(sigma) = 0. */
2176  ElasticityIntegrator(Coefficient &m, double q_l, double q_m)
2177  { lambda = NULL; mu = &m; q_lambda = q_l; q_mu = q_m; }
2178 
2179  virtual void AssembleElementMatrix(const FiniteElement &,
2181  DenseMatrix &);
2182 
2183  /** Compute the stress corresponding to the local displacement @a u and
2184  interpolate it at the nodes of the given @a fluxelem. Only the symmetric
2185  part of the stress is stored, so that the size of @a flux is equal to
2186  the number of DOFs in @a fluxelem times dim*(dim+1)/2. In 2D, the order
2187  of the stress components is: s_xx, s_yy, s_xy. In 3D, it is: s_xx, s_yy,
2188  s_zz, s_xy, s_xz, s_yz. In other words, @a flux is the local vector for
2189  a FE space with dim*(dim+1)/2 vector components, based on the finite
2190  element @a fluxelem. */
2191  virtual void ComputeElementFlux(const FiniteElement &el,
2193  Vector &u,
2194  const FiniteElement &fluxelem,
2195  Vector &flux, int with_coef = 1);
2196 
2197  /** Compute the element energy (integral of the strain energy density)
2198  corresponding to the stress represented by @a flux which is a vector of
2199  coefficients multiplying the basis functions defined by @a fluxelem. In
2200  other words, @a flux is the local vector for a FE space with
2201  dim*(dim+1)/2 vector components, based on the finite element @a fluxelem.
2202  The number of components, dim*(dim+1)/2 is such that it represents the
2203  symmetric part of the (symmetric) stress tensor. The order of the
2204  components is: s_xx, s_yy, s_xy in 2D, and s_xx, s_yy, s_zz, s_xy, s_xz,
2205  s_yz in 3D. */
2206  virtual double ComputeFluxEnergy(const FiniteElement &fluxelem,
2208  Vector &flux, Vector *d_energy = NULL);
2209 };
2210 
2211 /** Integrator for the DG form:
2212  alpha < rho_u (u.n) {v},[w] > + beta < rho_u |u.n| [v],[w] >,
2213  where v and w are the trial and test variables, respectively, and rho/u are
2214  given scalar/vector coefficients. The vector coefficient, u, is assumed to
2215  be continuous across the faces and when given the scalar coefficient, rho,
2216  is assumed to be discontinuous. The integrator uses the upwind value of rho,
2217  rho_u, which is value from the side into which the vector coefficient, u,
2218  points. */
2220 {
2221 protected:
2224  double alpha, beta;
2225 
2226 private:
2227  Vector shape1, shape2;
2228 
2229 public:
2230  /// Construct integrator with rho = 1.
2231  DGTraceIntegrator(VectorCoefficient &_u, double a, double b)
2232  { rho = NULL; u = &_u; alpha = a; beta = b; }
2233 
2235  double a, double b)
2236  { rho = &_rho; u = &_u; alpha = a; beta = b; }
2237 
2239  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2240  const FiniteElement &el2,
2242  DenseMatrix &elmat);
2243 };
2244 
2245 /** Integrator for the DG form:
2246 
2247  - < {(Q grad(u)).n}, [v] > + sigma < [u], {(Q grad(v)).n} >
2248  + kappa < {h^{-1} Q} [u], [v] >,
2249 
2250  where Q is a scalar or matrix diffusion coefficient and u, v are the trial
2251  and test spaces, respectively. The parameters sigma and kappa determine the
2252  DG method to be used (when this integrator is added to the "broken"
2253  DiffusionIntegrator):
2254  * sigma = -1, kappa >= kappa0: symm. interior penalty (IP or SIPG) method,
2255  * sigma = +1, kappa > 0: non-symmetric interior penalty (NIPG) method,
2256  * sigma = +1, kappa = 0: the method of Baumann and Oden. */
2258 {
2259 protected:
2262  double sigma, kappa;
2263 
2264  // these are not thread-safe!
2267 
2268 public:
2269  DGDiffusionIntegrator(const double s, const double k)
2270  : Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
2271  DGDiffusionIntegrator(Coefficient &q, const double s, const double k)
2272  : Q(&q), MQ(NULL), sigma(s), kappa(k) { }
2273  DGDiffusionIntegrator(MatrixCoefficient &q, const double s, const double k)
2274  : Q(NULL), MQ(&q), sigma(s), kappa(k) { }
2276  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2277  const FiniteElement &el2,
2279  DenseMatrix &elmat);
2280 };
2281 
2282 /** Integrator for the DG elasticity form, for the formulations see:
2283  - PhD Thesis of Jonas De Basabe, High-Order Finite %Element Methods for
2284  Seismic Wave Propagation, UT Austin, 2009, p. 23, and references therein
2285  - Peter Hansbo and Mats G. Larson, Discontinuous Galerkin and the
2286  Crouzeix-Raviart %Element: Application to Elasticity, PREPRINT 2000-09,
2287  p.3
2288 
2289  \f[
2290  - \left< \{ \tau(u) \}, [v] \right> + \alpha \left< \{ \tau(v) \}, [u]
2291  \right> + \kappa \left< h^{-1} \{ \lambda + 2 \mu \} [u], [v] \right>
2292  \f]
2293 
2294  where \f$ \left<u, v\right> = \int_{F} u \cdot v \f$, and \f$ F \f$ is a
2295  face which is either a boundary face \f$ F_b \f$ of an element \f$ K \f$ or
2296  an interior face \f$ F_i \f$ separating elements \f$ K_1 \f$ and \f$ K_2 \f$.
2297 
2298  In the bilinear form above \f$ \tau(u) \f$ is traction, and it's also
2299  \f$ \tau(u) = \sigma(u) \cdot \vec{n} \f$, where \f$ \sigma(u) \f$ is
2300  stress, and \f$ \vec{n} \f$ is the unit normal vector w.r.t. to \f$ F \f$.
2301 
2302  In other words, we have
2303  \f[
2304  - \left< \{ \sigma(u) \cdot \vec{n} \}, [v] \right> + \alpha \left< \{
2305  \sigma(v) \cdot \vec{n} \}, [u] \right> + \kappa \left< h^{-1} \{
2306  \lambda + 2 \mu \} [u], [v] \right>
2307  \f]
2308 
2309  For isotropic media
2310  \f[
2311  \begin{split}
2312  \sigma(u) &= \lambda \nabla \cdot u I + 2 \mu \varepsilon(u) \\
2313  &= \lambda \nabla \cdot u I + 2 \mu \frac{1}{2} (\nabla u + \nabla
2314  u^T) \\
2315  &= \lambda \nabla \cdot u I + \mu (\nabla u + \nabla u^T)
2316  \end{split}
2317  \f]
2318 
2319  where \f$ I \f$ is identity matrix, \f$ \lambda \f$ and \f$ \mu \f$ are Lame
2320  coefficients (see ElasticityIntegrator), \f$ u, v \f$ are the trial and test
2321  functions, respectively.
2322 
2323  The parameters \f$ \alpha \f$ and \f$ \kappa \f$ determine the DG method to
2324  use (when this integrator is added to the "broken" ElasticityIntegrator):
2325 
2326  - IIPG, \f$\alpha = 0\f$,
2327  C. Dawson, S. Sun, M. Wheeler, Compatible algorithms for coupled flow and
2328  transport, Comp. Meth. Appl. Mech. Eng., 193(23-26), 2565-2580, 2004.
2329 
2330  - SIPG, \f$\alpha = -1\f$,
2331  M. Grote, A. Schneebeli, D. Schotzau, Discontinuous Galerkin Finite
2332  %Element Method for the Wave Equation, SINUM, 44(6), 2408-2431, 2006.
2333 
2334  - NIPG, \f$\alpha = 1\f$,
2335  B. Riviere, M. Wheeler, V. Girault, A Priori Error Estimates for Finite
2336  %Element Methods Based on Discontinuous Approximation Spaces for Elliptic
2337  Problems, SINUM, 39(3), 902-931, 2001.
2338 
2339  This is a '%Vector' integrator, i.e. defined for FE spaces using multiple
2340  copies of a scalar FE space.
2341  */
2343 {
2344 public:
2345  DGElasticityIntegrator(double alpha_, double kappa_)
2346  : lambda(NULL), mu(NULL), alpha(alpha_), kappa(kappa_) { }
2347 
2349  double alpha_, double kappa_)
2350  : lambda(&lambda_), mu(&mu_), alpha(alpha_), kappa(kappa_) { }
2351 
2353  virtual void AssembleFaceMatrix(const FiniteElement &el1,
2354  const FiniteElement &el2,
2356  DenseMatrix &elmat);
2357 
2358 protected:
2360  double alpha, kappa;
2361 
2362 #ifndef MFEM_THREAD_SAFE
2363  // values of all scalar basis functions for one component of u (which is a
2364  // vector) at the integration point in the reference space
2366  // values of derivatives of all scalar basis functions for one component
2367  // of u (which is a vector) at the integration point in the reference space
2369  // Adjugate of the Jacobian of the transformation: adjJ = det(J) J^{-1}
2371  // gradient of shape functions in the real (physical, not reference)
2372  // coordinates, scaled by det(J):
2373  // dshape_ps(jdof,jm) = sum_{t} adjJ(t,jm)*dshape(jdof,t)
2375  Vector nor; // nor = |weight(J_face)| n
2376  Vector nL1, nL2; // nL1 = (lambda1 * ip.weight / detJ1) nor
2377  Vector nM1, nM2; // nM1 = (mu1 * ip.weight / detJ1) nor
2378  Vector dshape1_dnM, dshape2_dnM; // dshape1_dnM = dshape1_ps . nM1
2379  // 'jmat' corresponds to the term: kappa <h^{-1} {lambda + 2 mu} [u], [v]>
2381 #endif
2382 
2383  static void AssembleBlock(
2384  const int dim, const int row_ndofs, const int col_ndofs,
2385  const int row_offset, const int col_offset,
2386  const double jmatcoef, const Vector &col_nL, const Vector &col_nM,
2387  const Vector &row_shape, const Vector &col_shape,
2388  const Vector &col_dshape_dnM, const DenseMatrix &col_dshape,
2389  DenseMatrix &elmat, DenseMatrix &jmat);
2390 };
2391 
2392 /** Integrator for the DPG form: < v, [w] > over all faces (the interface) where
2393  the trial variable v is defined on the interface and the test variable w is
2394  defined inside the elements, generally in a DG space. */
2396 {
2397 private:
2398  Vector face_shape, shape1, shape2;
2399 
2400 public:
2403  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
2404  const FiniteElement &test_fe1,
2405  const FiniteElement &test_fe2,
2407  DenseMatrix &elmat);
2408 };
2409 
2410 /** Integrator for the form: < v, [w.n] > over all faces (the interface) where
2411  the trial variable v is defined on the interface and the test variable w is
2412  in an H(div)-conforming space. */
2414 {
2415 private:
2416  Vector face_shape, normal, shape1_n, shape2_n;
2417  DenseMatrix shape1, shape2;
2418 
2419 public:
2422  virtual void AssembleFaceMatrix(const FiniteElement &trial_face_fe,
2423  const FiniteElement &test_fe1,
2424  const FiniteElement &test_fe2,
2426  DenseMatrix &elmat);
2427 };
2428 
2429 /** Abstract class to serve as a base for local interpolators to be used in the
2430  DiscreteLinearOperator class. */
2432 
2433 
2434 /** Class for constructing the gradient as a DiscreteLinearOperator from an
2435  H1-conforming space to an H(curl)-conforming space. The range space can be
2436  vector L2 space as well. */
2438 {
2439 public:
2440  virtual void AssembleElementMatrix2(const FiniteElement &h1_fe,
2441  const FiniteElement &nd_fe,
2443  DenseMatrix &elmat)
2444  { nd_fe.ProjectGrad(h1_fe, Trans, elmat); }
2445 };
2446 
2447 
2448 /** Class for constructing the identity map as a DiscreteLinearOperator. This
2449  is the discrete embedding matrix when the domain space is a subspace of
2450  the range space. Otherwise, a dof projection matrix is constructed. */
2452 {
2453 public:
2454  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2455  const FiniteElement &ran_fe,
2457  DenseMatrix &elmat)
2458  { ran_fe.Project(dom_fe, Trans, elmat); }
2459 };
2460 
2461 
2462 /** Class for constructing the (local) discrete curl matrix which can be used
2463  as an integrator in a DiscreteLinearOperator object to assemble the global
2464  discrete curl matrix. */
2466 {
2467 public:
2468  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2469  const FiniteElement &ran_fe,
2471  DenseMatrix &elmat)
2472  { ran_fe.ProjectCurl(dom_fe, Trans, elmat); }
2473 };
2474 
2475 
2476 /** Class for constructing the (local) discrete divergence matrix which can
2477  be used as an integrator in a DiscreteLinearOperator object to assemble
2478  the global discrete divergence matrix.
2479 
2480  Note: Since the dofs in the L2_FECollection are nodal values, the local
2481  discrete divergence matrix (with an RT-type domain space) will depend on
2482  the transformation. On the other hand, the local matrix returned by
2483  VectorFEDivergenceIntegrator is independent of the transformation. */
2485 {
2486 public:
2487  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2488  const FiniteElement &ran_fe,
2490  DenseMatrix &elmat)
2491  { ran_fe.ProjectDiv(dom_fe, Trans, elmat); }
2492 };
2493 
2494 
2495 /** A trace face interpolator class for interpolating the normal component of
2496  the domain space, e.g. vector H1, into the range space, e.g. the trace of
2497  RT which uses FiniteElement::INTEGRAL map type. */
2499 {
2500 public:
2501  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2502  const FiniteElement &ran_fe,
2504  DenseMatrix &elmat);
2505 };
2506 
2507 /** Interpolator of a scalar coefficient multiplied by a scalar field onto
2508  another scalar field. Note that this can produce inaccurate fields unless
2509  the target is sufficiently high order. */
2511 {
2512 public:
2514 
2515  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2516  const FiniteElement &ran_fe,
2518  DenseMatrix &elmat);
2519 
2520 protected:
2522 };
2523 
2524 /** Interpolator of a scalar coefficient multiplied by a vector field onto
2525  another vector field. Note that this can produce inaccurate fields unless
2526  the target is sufficiently high order. */
2528 {
2529 public:
2531  : Q(&sc) { }
2532 
2533  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2534  const FiniteElement &ran_fe,
2536  DenseMatrix &elmat);
2537 protected:
2539 };
2540 
2541 /** Interpolator of a vector coefficient multiplied by a scalar field onto
2542  another vector field. Note that this can produce inaccurate fields unless
2543  the target is sufficiently high order. */
2545 {
2546 public:
2548  : VQ(&vc) { }
2549 
2550  virtual void AssembleElementMatrix2(const FiniteElement &dom_fe,
2551  const FiniteElement &ran_fe,
2553  DenseMatrix &elmat);
2554 protected:
2556 };
2557 
2558 /** Interpolator of the cross product between a vector coefficient and an
2559  H(curl)-conforming field onto an H(div)-conforming field. The range space
2560  can also be vector L2. */
2562 {
2563 public:
2565  : VQ(&vc) { }
2566 
2567  virtual void AssembleElementMatrix2(const FiniteElement &nd_fe,
2568  const FiniteElement &rt_fe,
2570  DenseMatrix &elmat);
2571 protected:
2573 };
2574 
2575 /** Interpolator of the inner product between a vector coefficient and an
2576  H(div)-conforming field onto an L2-conforming field. The range space can
2577  also be H1. */
2579 {
2580 public:
2582 
2583  virtual void AssembleElementMatrix2(const FiniteElement &rt_fe,
2584  const FiniteElement &l2_fe,
2586  DenseMatrix &elmat);
2587 protected:
2589 };
2590 
2591 }
2592 
2593 #endif
Implements CalcDivShape methods.
Definition: fe.hpp:290
Abstract class for Finite Elements.
Definition: fe.hpp:229
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:548
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
Definition: bilininteg.hpp:268
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:941
virtual void AssembleElementMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
Support for use in BilinearForm. Can be used only when appropriate.
Definition: bilininteg.hpp:335
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:83
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:305
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:703
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:85
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:674
MixedCrossCurlGradIntegrator(VectorCoefficient &vq)
DGTraceIntegrator(VectorCoefficient &_u, double a, double b)
Construct integrator with rho = 1.
MixedGradDivIntegrator(VectorCoefficient &vq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
DGDiffusionIntegrator(Coefficient &q, const double s, const double k)
Integrator defining a sum of multiple Integrators.
Definition: bilininteg.hpp:236
SumIntegrator(int own_integs=1)
Definition: bilininteg.hpp:244
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:356
MixedScalarCrossProductIntegrator(VectorCoefficient &vq)
CurlCurlIntegrator(MatrixCoefficient &m)
MatrixCoefficient * MQ
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:948
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:451
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 CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:834
ElasticityIntegrator(Coefficient &l, Coefficient &m)
VectorFEMassIntegrator(VectorCoefficient *_vq)
VectorCoefficient * VQ
Definition: bilininteg.hpp:387
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:302
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:529
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:565
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
MixedScalarMassIntegrator(Coefficient &q)
Definition: bilininteg.hpp:489
VectorMassIntegrator(Coefficient &q, const IntegrationRule *ir)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:917
alpha (q . grad u, v) using the &quot;group&quot; FE discretization
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:927
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:315
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:659
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:150
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 ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:996
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:100
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:109
BilinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: bilininteg.hpp:26
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:325
void AddIntegrator(BilinearFormIntegrator *integ)
Definition: bilininteg.hpp:246
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:895
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1284
MixedVectorCurlIntegrator(VectorCoefficient &dq)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:307
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:162
MixedCrossProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:780
virtual void AddMultTransposePA(const Vector &x, Vector &y) const
Method for partially assembled transposed action.
Definition: bilininteg.cpp:35
MixedScalarVectorIntegrator(VectorCoefficient &vq, bool _transpose=false, bool _cross_2d=false)
Definition: bilininteg.hpp:418
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:159
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorGradientIntegrator(Coefficient &q)
(Q div u, div v) for RT elements
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:818
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:435
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:912
const DofToQuad * maps
Not owned.
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Definition: bilininteg.cpp:549
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:781
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:885
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:147
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
const GeometricFactors * geom
Not owned.
MixedCurlCurlIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:970
MatrixCoefficient * MQ
VectorFEMassIntegrator(VectorCoefficient &vq)
VectorMassIntegrator()
Construct an integrator with coefficient 1.0.
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
VectorMassIntegrator(VectorCoefficient &q, int qo=0)
Construct an integrator with diagonal coefficient q.
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:282
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:462
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:315
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
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:91
MixedGradGradIntegrator(VectorCoefficient &dq)
Definition: bilininteg.hpp:880
Integrator that inverts the matrix assembled by another integrator.
Definition: bilininteg.hpp:218
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:388
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
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:696
VectorInnerProductInterpolator(VectorCoefficient &vc)
MixedWeakDivCrossIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:845
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
int dim
Definition: ex3.cpp:48
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:116
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
MatrixCoefficient * MQ
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedVectorWeakCurlIntegrator(VectorCoefficient &dq)
MixedScalarWeakCurlIntegrator(Coefficient &q)
Definition: bilininteg.hpp:731
MixedCurlCurlIntegrator(Coefficient &q)
Definition: bilininteg.hpp:966
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:586
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:123
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: bilininteg.hpp:97
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.cpp:766
VectorFEDivergenceIntegrator(Coefficient &q)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:991
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:848
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
VectorFECurlIntegrator(Coefficient &q)
virtual void AssembleElementMatrix2(const FiniteElement &h1_fe, const FiniteElement &nd_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Definition: bilininteg.hpp:137
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:364
DiffusionIntegrator()
Construct a diffusion integrator with coefficient Q = 1.
MixedVectorMassIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:771
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:582
MixedVectorIntegrator(MatrixCoefficient &mq)
Definition: bilininteg.hpp:353
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:376
MixedVectorMassIntegrator(Coefficient &q)
Definition: bilininteg.hpp:767
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:637
ScalarProductInterpolator(Coefficient &sc)
virtual const char * FiniteElementTypeFailureMessage() const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:735
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)
Definition: bilininteg.hpp:968
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:751
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:594
BoundaryMassIntegrator(Coefficient &q)
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:873
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:973
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
Definition: bilininteg.cpp:74
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
Definition: bilininteg.cpp:29
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the BilinearFormIntegrator.
Definition: bilininteg.cpp:493
MixedWeakGradDotIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:815
ScalarVectorProductInterpolator(Coefficient &sc)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:801
virtual void CalcShape(const FiniteElement &scalar_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:642
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:23
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:381
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorDivergenceIntegrator(Coefficient *_q)
MixedVectorIntegrator(Coefficient &q)
Definition: bilininteg.hpp:348
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:744
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:882
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
Definition: bilininteg.cpp:825
VectorFEMassIntegrator(MatrixCoefficient &mq)
MixedVectorMassIntegrator(VectorCoefficient &dq)
Definition: bilininteg.hpp:769
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
MixedScalarWeakCurlCrossIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:694
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:41
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:622
Implements CalcDShape methods.
Definition: fe.hpp:289
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Definition: fe.cpp:185
virtual const char * FiniteElementTypeFailureMessage() const
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
MixedScalarWeakDivergenceIntegrator(VectorCoefficient &vq)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:434
DGDiffusionIntegrator(const double s, const double k)
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:865
LumpedIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
Definition: bilininteg.hpp:207
static const IntegrationRule & GetRule(const FiniteElement &trial_fe, const FiniteElement &test_fe)
Definition: bilininteg.cpp:673
MixedCrossGradCurlIntegrator(VectorCoefficient &vq)
MixedVectorDivergenceIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:618
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe.hpp:125
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:55
MassIntegrator(Coefficient &q, const IntegrationRule *ir=NULL)
Construct a mass integrator with coefficient q.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:290
virtual void AddMultPA(const Vector &, Vector &) const
Method for partially assembled action.
MixedCrossGradGradIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:928
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:601
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
MixedVectorIntegrator(VectorCoefficient &dq, bool diag=true)
Definition: bilininteg.hpp:350
MixedScalarCrossCurlIntegrator(VectorCoefficient &vq)
MatrixCoefficient * MQ
Definition: bilininteg.hpp:389
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Assemble an element matrix.
virtual void AssembleFaceMatrix(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:57
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:667
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:827
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:858
VectorMassIntegrator(Coefficient &q, int qo=0)
VectorFEMassIntegrator(Coefficient *_q)
DGElasticityIntegrator(double alpha_, double kappa_)
MixedScalarDerivativeIntegrator(Coefficient &q)
Definition: bilininteg.hpp:509
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
VectorFEMassIntegrator(MatrixCoefficient *_mq)
MixedDotProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:790
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
MixedVectorCurlIntegrator(Coefficient &q)
MixedWeakCurlCrossIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
virtual void AssemblePA(const FiniteElementSpace &)
Method defining partial assembly.
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:557
Implements CalcCurlShape methods.
Definition: fe.hpp:291
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:728
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix2(const FiniteElement &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:206
VectorMassIntegrator(MatrixCoefficient &q, int qo=0)
Construct an integrator with matrix coefficient q.
MixedScalarWeakCrossProductIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
virtual void AssembleElementMatrix(const FiniteElement &, ElementTransformation &, DenseMatrix &)
Given a particular Finite Element computes the element matrix elmat.
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:370
GroupConvectionIntegrator(VectorCoefficient &q, double a=1.0)
virtual void CalcVShape(const FiniteElement &vector_fe, ElementTransformation &Trans, DenseMatrix &shape)
Definition: bilininteg.hpp:457
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:953
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:136
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:522
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:984
virtual const char * FiniteElementTypeFailureMessage() const
DGDiffusionIntegrator(MatrixCoefficient &q, const double s, const double k)
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:715
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:102
MixedCrossCurlIntegrator(VectorCoefficient &vq)
Vector data type.
Definition: vector.hpp:48
int GetDerivType() const
Definition: fe.hpp:333
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual const char * FiniteElementTypeFailureMessage() const
Definition: bilininteg.hpp:630
MixedVectorProductIntegrator(VectorCoefficient &vq)
Definition: bilininteg.hpp:499
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:690
virtual void AssembleElementMatrix2(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:49
virtual void CalcTestShape(const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
DerivativeIntegrator(Coefficient &q, int i)
DGTraceIntegrator(Coefficient &_rho, VectorCoefficient &_u, double a, double b)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:793
CurlCurlIntegrator(Coefficient &q)
Construct a bilinear form integrator for Nedelec elements.
virtual const char * FiniteElementTypeFailureMessage() const
MixedGradGradIntegrator(Coefficient &q)
Definition: bilininteg.hpp:878
MixedVectorGradientIntegrator(MatrixCoefficient &mq)
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:902
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:327
virtual void CalcTrialShape(const FiniteElement &trial_fe, ElementTransformation &Trans, Vector &shape)
Definition: bilininteg.hpp:606
virtual void AssembleElementMatrix2(const FiniteElement &dom_fe, const FiniteElement &ran_fe, ElementTransformation &Trans, DenseMatrix &elmat)
virtual void AssembleElementMatrix(const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
Definition: bilininteg.cpp:384
virtual const char * FiniteElementTypeFailureMessage() const
virtual const char * FiniteElementTypeFailureMessage() const
MassIntegrator(const IntegrationRule *ir=NULL)
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef)
Virtual method required for Zienkiewicz-Zhu type error estimators.
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:296
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:607
Integrator for (Q u, v) for VectorFiniteElements.
MixedVectorWeakCurlIntegrator(Coefficient &q)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:931
TransposeIntegrator(BilinearFormIntegrator *_bfi, int _own_bfi=1)
Definition: bilininteg.hpp:179
MixedDivGradIntegrator(VectorCoefficient &vq)
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:513
virtual int GetIntegrationOrder(const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
Definition: bilininteg.hpp:710
MixedScalarIntegrator(Coefficient &q)
Definition: bilininteg.hpp:280
Polynomials of order k.
Definition: fe.hpp:212
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:225
virtual bool VerifyFiniteElementTypes(const FiniteElement &trial_fe, const FiniteElement &test_fe) const
Definition: bilininteg.hpp:422
VectorCoefficient * Q
alpha (q . grad u, v)