MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lininteg.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_LININTEG
13 #define MFEM_LININTEG
14 
15 #include "../config/config.hpp"
16 #include "coefficient.hpp"
17 
18 namespace mfem
19 {
20 
21 /// Abstract base class LinearFormIntegrator
23 {
24 protected:
26 
27  LinearFormIntegrator(const IntegrationRule *ir = NULL) { IntRule = ir; }
28 
29 public:
30  /** Given a particular Finite Element and a transformation (Tr)
31  computes the element vector, elvect. */
32  virtual void AssembleRHSElementVect(const FiniteElement &el,
34  Vector &elvect) = 0;
35  virtual void AssembleRHSElementVect(const FiniteElement &el,
37  Vector &elvect);
38  virtual void AssembleRHSElementVect(const FiniteElement &el1,
39  const FiniteElement &el2,
41  Vector &elvect);
42 
43  virtual void SetIntRule(const IntegrationRule *ir) { IntRule = ir; }
44  const IntegrationRule* GetIntRule() { return IntRule; }
45 
46  virtual ~LinearFormIntegrator() { }
47 };
48 
49 
50 /// Abstract class for integrators that support delta coefficients
52 {
53 protected:
56 
57  /** @brief This constructor should be used by derived classes that use a
58  scalar DeltaCoefficient. */
61  delta(dynamic_cast<DeltaCoefficient*>(&q)),
62  vec_delta(NULL) { }
63 
64  /** @brief This constructor should be used by derived classes that use a
65  VectorDeltaCoefficient. */
67  const IntegrationRule *ir = NULL)
69  delta(NULL),
70  vec_delta(dynamic_cast<VectorDeltaCoefficient*>(&vq)) { }
71 
72 public:
73  /// Returns true if the derived class instance uses a delta coefficient.
74  bool IsDelta() const { return (delta || vec_delta); }
75 
76  /// Returns the center of the delta coefficient.
77  void GetDeltaCenter(Vector &center)
78  {
79  if (delta) { delta->GetDeltaCenter(center); return; }
80  if (vec_delta) { vec_delta->GetDeltaCenter(center); return; }
81  center.SetSize(0);
82  }
83 
84  /** @brief Assemble the delta coefficient at the IntegrationPoint set in
85  @a Trans which is assumed to map to the delta coefficient center.
86 
87  @note This method should be called for one mesh element only, including
88  in parallel, even when the center of the delta coefficient is shared by
89  multiple elements. */
90  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
92  Vector &elvect) = 0;
93 };
94 
95 
96 /// Class for domain integration L(v) := (f, v)
98 {
99  Vector shape;
100  Coefficient &Q;
101  int oa, ob;
102 public:
103  /// Constructs a domain integrator with a given Coefficient
104  DomainLFIntegrator(Coefficient &QF, int a = 2, int b = 0)
105  // the old default was a = 1, b = 1
106  // for simple elliptic problems a = 2, b = -2 is OK
107  : DeltaLFIntegrator(QF), Q(QF), oa(a), ob(b) { }
108 
109  /// Constructs a domain integrator with a given Coefficient
111  : DeltaLFIntegrator(QF, ir), Q(QF), oa(1), ob(1) { }
112 
113  /** Given a particular Finite Element and a transformation (Tr)
114  computes the element right hand side element vector, elvect. */
115  virtual void AssembleRHSElementVect(const FiniteElement &el,
117  Vector &elvect);
118 
119  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
121  Vector &elvect);
122 
124 };
125 
126 /// Class for domain integrator L(v) := (f, grad v)
128 {
129 private:
130  Vector shape, Qvec;
132  DenseMatrix dshape;
133 
134 public:
135  /// Constructs the domain integrator (Q, grad v)
137  : DeltaLFIntegrator(QF), Q(QF) { }
138 
139  /** Given a particular Finite Element and a transformation (Tr)
140  computes the element right hand side element vector, elvect. */
141  virtual void AssembleRHSElementVect(const FiniteElement &el,
143  Vector &elvect);
144 
145  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
147  Vector &elvect);
148 
150 };
151 
152 
153 /// Class for boundary integration L(v) := (g, v)
155 {
156  Vector shape;
157  Coefficient &Q;
158  int oa, ob;
159 public:
160  /** @brief Constructs a boundary integrator with a given Coefficient @a QG.
161  Integration order will be @a a * basis_order + @a b. */
162  BoundaryLFIntegrator(Coefficient &QG, int a = 1, int b = 1)
163  : Q(QG), oa(a), ob(b) { }
164 
165  /** Given a particular boundary Finite Element and a transformation (Tr)
166  computes the element boundary vector, elvect. */
167  virtual void AssembleRHSElementVect(const FiniteElement &el,
169  Vector &elvect);
170  virtual void AssembleRHSElementVect(const FiniteElement &el,
172  Vector &elvect);
173 
175 };
176 
177 /// Class for boundary integration \f$ L(v) = (g \cdot n, v) \f$
179 {
180  Vector shape;
182  int oa, ob;
183 public:
184  /// Constructs a boundary integrator with a given Coefficient QG
186  : Q(QG), oa(a), ob(b) { }
187 
188  virtual void AssembleRHSElementVect(const FiniteElement &el,
190  Vector &elvect);
191 
193 };
194 
195 /// Class for boundary integration \f$ L(v) = (g \cdot \tau, v) \f$ in 2D
197 {
198  Vector shape;
200  int oa, ob;
201 public:
202  /// Constructs a boundary integrator with a given Coefficient QG
204  : Q(QG), oa(a), ob(b) { }
205 
206  virtual void AssembleRHSElementVect(const FiniteElement &el,
208  Vector &elvect);
209 
211 };
212 
213 /** Class for domain integration of L(v) := (f, v), where
214  f=(f1,...,fn) and v=(v1,...,vn). */
216 {
217 private:
218  Vector shape, Qvec;
220 
221 public:
222  /// Constructs a domain integrator with a given VectorCoefficient
224  : DeltaLFIntegrator(QF), Q(QF) { }
225 
226  /** Given a particular Finite Element and a transformation (Tr)
227  computes the element right hand side element vector, elvect. */
228  virtual void AssembleRHSElementVect(const FiniteElement &el,
230  Vector &elvect);
231 
232  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
234  Vector &elvect);
235 
237 };
238 
239 /** Class for boundary integration of L(v) := (g, v), where
240  f=(f1,...,fn) and v=(v1,...,vn). */
242 {
243 private:
244  Vector shape, vec;
246 
247 public:
248  /// Constructs a boundary integrator with a given VectorCoefficient QG
250 
251  /** Given a particular boundary Finite Element and a transformation (Tr)
252  computes the element boundary vector, elvect. */
253  virtual void AssembleRHSElementVect(const FiniteElement &el,
255  Vector &elvect);
256 
257  // For DG spaces
258  virtual void AssembleRHSElementVect(const FiniteElement &el,
260  Vector &elvect);
261 
263 };
264 
265 /// \f$ (f, v)_{\Omega} \f$ for VectorFiniteElements (Nedelec, Raviart-Thomas)
267 {
268 private:
269  VectorCoefficient &QF;
270  DenseMatrix vshape;
271  Vector vec;
272 
273 public:
275  : DeltaLFIntegrator(F), QF(F) { }
276 
277  virtual void AssembleRHSElementVect(const FiniteElement &el,
279  Vector &elvect);
280 
281  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
283  Vector &elvect);
284 
286 };
287 
288 /// \f$ (Q, curl v)_{\Omega} \f$ for Nedelec Elements)
290 {
291 private:
292  VectorCoefficient *QF=nullptr;
293  DenseMatrix curlshape;
294  Vector vec;
295 
296 public:
297  /// Constructs the domain integrator (Q, curl v)
299  : DeltaLFIntegrator(F), QF(&F) { }
300 
301  virtual void AssembleRHSElementVect(const FiniteElement &el,
303  Vector &elvect);
304 
305  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
307  Vector &elvect);
308 
310 };
311 
312 /// \f$ (Q, div v)_{\Omega} \f$ for RT Elements)
314 {
315 private:
316  Vector divshape;
317  Coefficient &Q;
318 public:
319  /// Constructs the domain integrator (Q, div v)
321  : DeltaLFIntegrator(QF), Q(QF) { }
322 
323  /** Given a particular Finite Element and a transformation (Tr)
324  computes the element right hand side element vector, elvect. */
325  virtual void AssembleRHSElementVect(const FiniteElement &el,
327  Vector &elvect);
328 
329  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
331  Vector &elvect);
332 
334 };
335 
336 /** \f$ (f, v \cdot n)_{\partial\Omega} \f$ for vector test function
337  v=(v1,...,vn) where all vi are in the same scalar FE space and f is a
338  scalar function. */
340 {
341 private:
342  double Sign;
343  Coefficient *F;
344  Vector shape, nor;
345 
346 public:
348  const IntegrationRule *ir = NULL)
349  : LinearFormIntegrator(ir), Sign(s), F(&f) { }
350 
351  virtual void AssembleRHSElementVect(const FiniteElement &el,
353  Vector &elvect);
354 
356 };
357 
358 /** Class for boundary integration of (f, v.n) for scalar coefficient f and
359  RT vector test function v. This integrator works with RT spaces defined
360  using the RT_FECollection class. */
362 {
363 private:
364  Coefficient *F;
365  Vector shape;
366  int oa, ob; // these control the quadrature order, see DomainLFIntegrator
367 
368 public:
370  : F(NULL), oa(a), ob(b) { }
372  : F(&f), oa(a), ob(b) { }
373 
374  virtual void AssembleRHSElementVect(const FiniteElement &el,
376  Vector &elvect);
377 
379 };
380 
381 /// Class for boundary integration \f$ L(v) = (n \times f, v) \f$
383 {
384 private:
386  int oa, ob;
387 
388 public:
390  int a = 2, int b = 0)
391  : f(QG), oa(a), ob(b) { }
392 
393  virtual void AssembleRHSElementVect(const FiniteElement &el,
395  Vector &elvect);
396 
398 };
399 
400 
401 /** Class for boundary integration of the linear form:
402  (alpha/2) < (u.n) f, w > - beta < |u.n| f, w >,
403  where f and u are given scalar and vector coefficients, respectively,
404  and w is the scalar test function. */
406 {
407 private:
408  Coefficient *f;
410  double alpha, beta;
411 
412  Vector shape;
413 
414 public:
416  double a)
417  { f = &f_; u = &u_; alpha = a; beta = 0.5*a; }
418 
420  double a, double b)
421  { f = &f_; u = &u_; alpha = a; beta = b; }
422 
423  virtual void AssembleRHSElementVect(const FiniteElement &el,
425  Vector &elvect);
426  virtual void AssembleRHSElementVect(const FiniteElement &el,
428  Vector &elvect);
429 
431 };
432 
433 
434 /** Boundary linear integrator for imposing non-zero Dirichlet boundary
435  conditions, to be used in conjunction with DGDiffusionIntegrator.
436  Specifically, given the Dirichlet data u_D, the linear form assembles the
437  following integrals on the boundary:
438 
439  sigma < u_D, (Q grad(v)).n > + kappa < {h^{-1} Q} u_D, v >,
440 
441  where Q is a scalar or matrix diffusion coefficient and v is the test
442  function. The parameters sigma and kappa should be the same as the ones
443  used in the DGDiffusionIntegrator. */
445 {
446 protected:
449  double sigma, kappa;
450 
451  // these are not thread-safe!
454 
455 public:
456  DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
457  : uD(&u), Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
459  const double s, const double k)
460  : uD(&u), Q(&q), MQ(NULL), sigma(s), kappa(k) { }
462  const double s, const double k)
463  : uD(&u), Q(NULL), MQ(&q), sigma(s), kappa(k) { }
464 
465  virtual void AssembleRHSElementVect(const FiniteElement &el,
467  Vector &elvect);
468  virtual void AssembleRHSElementVect(const FiniteElement &el,
470  Vector &elvect);
471 
473 };
474 
475 
476 /** Boundary linear form integrator for imposing non-zero Dirichlet boundary
477  conditions, in a DG elasticity formulation. Specifically, the linear form is
478  given by
479 
480  alpha < u_D, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
481  + kappa < h^{-1} (lambda + 2 mu) u_D, v >,
482 
483  where u_D is the given Dirichlet data. The parameters alpha, kappa, lambda
484  and mu, should match the parameters with the same names used in the bilinear
485  form integrator, DGElasticityIntegrator. */
487 {
488 protected:
491  double alpha, kappa;
492 
493 #ifndef MFEM_THREAD_SAFE
502 #endif
503 
504 public:
506  Coefficient &lambda_, Coefficient &mu_,
507  double alpha_, double kappa_)
508  : uD(uD_), lambda(&lambda_), mu(&mu_), alpha(alpha_), kappa(kappa_) { }
509 
510  virtual void AssembleRHSElementVect(const FiniteElement &el,
512  Vector &elvect);
513  virtual void AssembleRHSElementVect(const FiniteElement &el,
515  Vector &elvect);
516 
518 };
519 
520 /** Class for domain integration of L(v) := (f, v), where
521  f=(f1,...,fn) and v=(v1,...,vn). that makes use of
522  VectorQuadratureFunctionCoefficient*/
524 {
525 private:
527 
528 public:
530  const IntegrationRule *ir)
531  : LinearFormIntegrator(ir), vqfc(vqfc)
532  {
533  if (ir)
534  {
535  MFEM_WARNING("Integration rule not used in this class. "
536  "The QuadratureFunction integration rules are used instead");
537  }
538  }
539 
541  virtual void AssembleRHSElementVect(const FiniteElement &fe,
543  Vector &elvect);
544 
545  virtual void SetIntRule(const IntegrationRule *ir)
546  {
547  MFEM_WARNING("Integration rule not used in this class. "
548  "The QuadratureFunction integration rules are used instead");
549  }
550 };
551 
552 /** Class for domain integration L(v) := (f, v) that makes use
553  of QuadratureFunctionCoefficient. */
555 {
556 private:
558 
559 public:
561  const IntegrationRule *ir)
562  : LinearFormIntegrator(ir), qfc(qfc)
563  {
564  if (ir)
565  {
566  MFEM_WARNING("Integration rule not used in this class. "
567  "The QuadratureFunction integration rules are used instead");
568  }
569  }
570 
572  virtual void AssembleRHSElementVect(const FiniteElement &fe,
574  Vector &elvect);
575 
576  virtual void SetIntRule(const IntegrationRule *ir)
577  {
578  MFEM_WARNING("Integration rule not used in this class. "
579  "The QuadratureFunction integration rules are used instead");
580  }
581 };
582 
583 }
584 
585 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:235
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:104
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:568
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
void GetDeltaCenter(Vector &center)
VectorFEDomainLFCurlIntegrator(VectorCoefficient &F)
Constructs the domain integrator (Q, curl v)
Definition: lininteg.hpp:298
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:526
Base class for vector Coefficients that optionally depend on time and space.
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:498
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:656
VectorBoundaryLFIntegrator(VectorCoefficient &QG)
Constructs a boundary integrator with a given VectorCoefficient QG.
Definition: lininteg.hpp:249
Class for domain integrator L(v) := (f, grad v)
Definition: lininteg.hpp:127
BoundaryNormalLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:185
DGElasticityDirichletLFIntegrator(VectorCoefficient &uD_, Coefficient &lambda_, Coefficient &mu_, double alpha_, double kappa_)
Definition: lininteg.hpp:505
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:813
Vector quadrature function coefficient which requires that the quadrature rules used for this vector ...
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:72
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:63
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:480
DeltaLFIntegrator(VectorCoefficient &vq, const IntegrationRule *ir=NULL)
This constructor should be used by derived classes that use a VectorDeltaCoefficient.
Definition: lininteg.hpp:66
virtual ~LinearFormIntegrator()
Definition: lininteg.hpp:46
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
const IntegrationRule * GetIntRule()
Definition: lininteg.hpp:44
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:262
Class for boundary integration L(v) := (g, v)
Definition: lininteg.hpp:154
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:43
const IntegrationRule * IntRule
Definition: lininteg.hpp:25
BoundaryFlowIntegrator(Coefficient &f_, VectorCoefficient &u_, double a, double b)
Definition: lininteg.hpp:419
DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
Definition: lininteg.hpp:456
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:120
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:22
void GetDeltaCenter(Vector &center)
Returns the center of the delta coefficient.
Definition: lininteg.hpp:77
LinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:27
VectorBoundaryFluxLFIntegrator(Coefficient &f, double s=1.0, const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:347
DomainLFGradIntegrator(VectorCoefficient &QF)
Constructs the domain integrator (Q, grad v)
Definition: lininteg.hpp:136
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
BoundaryFlowIntegrator(Coefficient &f_, VectorCoefficient &u_, double a)
Definition: lininteg.hpp:415
Class for boundary integration .
Definition: lininteg.hpp:178
DomainLFIntegrator(Coefficient &QF, int a=2, int b=0)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:104
double b
Definition: lissajous.cpp:42
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:321
VectorFEBoundaryFluxLFIntegrator(int a=1, int b=-1)
Definition: lininteg.hpp:369
VectorFEBoundaryTangentLFIntegrator(VectorCoefficient &QG, int a=2, int b=0)
Definition: lininteg.hpp:389
BoundaryLFIntegrator(Coefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG. Integration order will be a * basis_ord...
Definition: lininteg.hpp:162
VectorDeltaCoefficient * vec_delta
Definition: lininteg.hpp:55
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:481
MatrixCoefficient * MQ
Definition: lininteg.hpp:448
VectorFEBoundaryFluxLFIntegrator(Coefficient &f, int a=2, int b=0)
Definition: lininteg.hpp:371
bool IsDelta() const
Returns true if the derived class instance uses a delta coefficient.
Definition: lininteg.hpp:74
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:545
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:222
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:722
DomainLFIntegrator(Coefficient &QF, const IntegrationRule *ir)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:110
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
VectorFEDomainLFDivIntegrator(Coefficient &QF)
Constructs the domain integrator (Q, div v)
Definition: lininteg.hpp:320
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)=0
DeltaLFIntegrator(Coefficient &q, const IntegrationRule *ir=NULL)
This constructor should be used by derived classes that use a scalar DeltaCoefficient.
Definition: lininteg.hpp:59
QuadratureLFIntegrator(QuadratureFunctionCoefficient &qfc, const IntegrationRule *ir)
Definition: lininteg.hpp:560
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
DeltaCoefficient * delta
Definition: lininteg.hpp:54
Base class for Matrix Coefficients that optionally depend on time and space.
Abstract class for integrators that support delta coefficients.
Definition: lininteg.hpp:51
VectorDomainLFIntegrator(VectorCoefficient &QF)
Constructs a domain integrator with a given VectorCoefficient.
Definition: lininteg.hpp:223
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:535
VectorFEDomainLFIntegrator(VectorCoefficient &F)
Definition: lininteg.hpp:274
Vector coefficient defined by a scalar DeltaCoefficient and a constant vector direction.
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:600
double a
Definition: lissajous.cpp:41
for Nedelec Elements)
Definition: lininteg.hpp:289
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:433
BoundaryTangentialLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:203
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:576
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:949
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:183
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:449
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:979
Vector data type.
Definition: vector.hpp:60
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:399
Class for boundary integration .
Definition: lininteg.hpp:382
RefCoord s[3]
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:304
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
DGDirichletLFIntegrator(Coefficient &u, MatrixCoefficient &q, const double s, const double k)
Definition: lininteg.hpp:461
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:32
Class for boundary integration in 2D.
Definition: lininteg.hpp:196
VectorQuadratureLFIntegrator(VectorQuadratureFunctionCoefficient &vqfc, const IntegrationRule *ir)
Definition: lininteg.hpp:529
DGDirichletLFIntegrator(Coefficient &u, Coefficient &q, const double s, const double k)
Definition: lininteg.hpp:458
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect)=0
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...