MFEM  v4.2.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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_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 
39  virtual void SetIntRule(const IntegrationRule *ir) { IntRule = ir; }
40  const IntegrationRule* GetIntRule() { return IntRule; }
41 
42  virtual ~LinearFormIntegrator() { }
43 };
44 
45 
46 /// Abstract class for integrators that support delta coefficients
48 {
49 protected:
52 
53  /** @brief This constructor should be used by derived classes that use a
54  scalar DeltaCoefficient. */
57  delta(dynamic_cast<DeltaCoefficient*>(&q)),
58  vec_delta(NULL) { }
59 
60  /** @brief This constructor should be used by derived classes that use a
61  VectorDeltaCoefficient. */
63  const IntegrationRule *ir = NULL)
65  delta(NULL),
66  vec_delta(dynamic_cast<VectorDeltaCoefficient*>(&vq)) { }
67 
68 public:
69  /// Returns true if the derived class instance uses a delta coefficient.
70  bool IsDelta() const { return (delta || vec_delta); }
71 
72  /// Returns the center of the delta coefficient.
73  void GetDeltaCenter(Vector &center)
74  {
75  if (delta) { delta->GetDeltaCenter(center); return; }
76  if (vec_delta) { vec_delta->GetDeltaCenter(center); return; }
77  center.SetSize(0);
78  }
79 
80  /** @brief Assemble the delta coefficient at the IntegrationPoint set in
81  @a Trans which is assumed to map to the delta coefficient center.
82 
83  @note This method should be called for one mesh element only, including
84  in parallel, even when the center of the delta coefficient is shared by
85  multiple elements. */
86  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
88  Vector &elvect) = 0;
89 };
90 
91 
92 /// Class for domain integration L(v) := (f, v)
94 {
95  Vector shape;
96  Coefficient &Q;
97  int oa, ob;
98 public:
99  /// Constructs a domain integrator with a given Coefficient
100  DomainLFIntegrator(Coefficient &QF, int a = 2, int b = 0)
101  // the old default was a = 1, b = 1
102  // for simple elliptic problems a = 2, b = -2 is OK
103  : DeltaLFIntegrator(QF), Q(QF), oa(a), ob(b) { }
104 
105  /// Constructs a domain integrator with a given Coefficient
107  : DeltaLFIntegrator(QF, ir), Q(QF), oa(1), ob(1) { }
108 
109  /** Given a particular Finite Element and a transformation (Tr)
110  computes the element right hand side element vector, elvect. */
111  virtual void AssembleRHSElementVect(const FiniteElement &el,
113  Vector &elvect);
114 
115  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
117  Vector &elvect);
118 
120 };
121 
122 /// Class for domain integrator L(v) := (f, grad v)
124 {
125 private:
126  Vector shape, Qvec;
128  DenseMatrix dshape;
129 
130 public:
131  /// Constructs the domain integrator (Q, grad v)
133  : DeltaLFIntegrator(QF), Q(QF) { }
134 
135  /** Given a particular Finite Element and a transformation (Tr)
136  computes the element right hand side element vector, elvect. */
137  virtual void AssembleRHSElementVect(const FiniteElement &el,
139  Vector &elvect);
140 
141  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
143  Vector &elvect);
144 
146 };
147 
148 
149 /// Class for boundary integration L(v) := (g, v)
151 {
152  Vector shape;
153  Coefficient &Q;
154  int oa, ob;
155 public:
156  /** @brief Constructs a boundary integrator with a given Coefficient @a QG.
157  Integration order will be @a a * basis_order + @a b. */
158  BoundaryLFIntegrator(Coefficient &QG, int a = 1, int b = 1)
159  : Q(QG), oa(a), ob(b) { }
160 
161  /** Given a particular boundary Finite Element and a transformation (Tr)
162  computes the element boundary vector, elvect. */
163  virtual void AssembleRHSElementVect(const FiniteElement &el,
165  Vector &elvect);
166  virtual void AssembleRHSElementVect(const FiniteElement &el,
168  Vector &elvect);
169 };
170 
171 /// Class for boundary integration \f$ L(v) = (g \cdot n, v) \f$
173 {
174  Vector shape;
176  int oa, ob;
177 public:
178  /// Constructs a boundary integrator with a given Coefficient QG
180  : Q(QG), oa(a), ob(b) { }
181 
182  virtual void AssembleRHSElementVect(const FiniteElement &el,
184  Vector &elvect);
185 
187 };
188 
189 /// Class for boundary integration \f$ L(v) = (g \cdot \tau, v) \f$ in 2D
191 {
192  Vector shape;
194  int oa, ob;
195 public:
196  /// Constructs a boundary integrator with a given Coefficient QG
198  : Q(QG), oa(a), ob(b) { }
199 
200  virtual void AssembleRHSElementVect(const FiniteElement &el,
202  Vector &elvect);
203 
205 };
206 
207 /** Class for domain integration of L(v) := (f, v), where
208  f=(f1,...,fn) and v=(v1,...,vn). */
210 {
211 private:
212  Vector shape, Qvec;
214 
215 public:
216  /// Constructs a domain integrator with a given VectorCoefficient
218  : DeltaLFIntegrator(QF), Q(QF) { }
219 
220  /** Given a particular Finite Element and a transformation (Tr)
221  computes the element right hand side element vector, elvect. */
222  virtual void AssembleRHSElementVect(const FiniteElement &el,
224  Vector &elvect);
225 
226  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
228  Vector &elvect);
229 
231 };
232 
233 /** Class for boundary integration of L(v) := (g, v), where
234  f=(f1,...,fn) and v=(v1,...,vn). */
236 {
237 private:
238  Vector shape, vec;
240 
241 public:
242  /// Constructs a boundary integrator with a given VectorCoefficient QG
244 
245  /** Given a particular boundary Finite Element and a transformation (Tr)
246  computes the element boundary vector, elvect. */
247  virtual void AssembleRHSElementVect(const FiniteElement &el,
249  Vector &elvect);
250 
251  // For DG spaces
252  virtual void AssembleRHSElementVect(const FiniteElement &el,
254  Vector &elvect);
255 
257 };
258 
259 /// \f$ (f, v)_{\Omega} \f$ for VectorFiniteElements (Nedelec, Raviart-Thomas)
261 {
262 private:
263  VectorCoefficient &QF;
264  DenseMatrix vshape;
265  Vector vec;
266 
267 public:
269  : DeltaLFIntegrator(F), QF(F) { }
270 
271  virtual void AssembleRHSElementVect(const FiniteElement &el,
273  Vector &elvect);
274 
275  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
277  Vector &elvect);
278 
280 };
281 
282 /// \f$ (Q, curl v)_{\Omega} \f$ for Nedelec Elements)
284 {
285 private:
286  VectorCoefficient *QF=nullptr;
287  DenseMatrix curlshape;
288  Vector vec;
289 
290 public:
291  /// Constructs the domain integrator (Q, curl v)
293  : DeltaLFIntegrator(F), QF(&F) { }
294 
295  virtual void AssembleRHSElementVect(const FiniteElement &el,
297  Vector &elvect);
298 
299  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
301  Vector &elvect);
302 
304 };
305 
306 /// \f$ (Q, div v)_{\Omega} \f$ for RT Elements)
308 {
309 private:
310  Vector divshape;
311  Coefficient &Q;
312 public:
313  /// Constructs the domain integrator (Q, div v)
315  : DeltaLFIntegrator(QF), Q(QF) { }
316 
317  /** Given a particular Finite Element and a transformation (Tr)
318  computes the element right hand side element vector, elvect. */
319  virtual void AssembleRHSElementVect(const FiniteElement &el,
321  Vector &elvect);
322 
323  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
325  Vector &elvect);
326 
328 };
329 
330 /** \f$ (f, v \cdot n)_{\partial\Omega} \f$ for vector test function
331  v=(v1,...,vn) where all vi are in the same scalar FE space and f is a
332  scalar function. */
334 {
335 private:
336  double Sign;
337  Coefficient *F;
338  Vector shape, nor;
339 
340 public:
342  const IntegrationRule *ir = NULL)
343  : LinearFormIntegrator(ir), Sign(s), F(&f) { }
344 
345  virtual void AssembleRHSElementVect(const FiniteElement &el,
347  Vector &elvect);
348 
350 };
351 
352 /** Class for boundary integration of (f, v.n) for scalar coefficient f and
353  RT vector test function v. This integrator works with RT spaces defined
354  using the RT_FECollection class. */
356 {
357 private:
358  Coefficient *F;
359  Vector shape;
360  int oa, ob; // these control the quadrature order, see DomainLFIntegrator
361 
362 public:
364  : F(NULL), oa(a), ob(b) { }
366  : F(&f), oa(a), ob(b) { }
367 
368  virtual void AssembleRHSElementVect(const FiniteElement &el,
370  Vector &elvect);
371 
373 };
374 
375 /// Class for boundary integration \f$ L(v) = (n \times f, v) \f$
377 {
378 private:
380  int oa, ob;
381 
382 public:
384  int a = 2, int b = 0)
385  : f(QG), oa(a), ob(b) { }
386 
387  virtual void AssembleRHSElementVect(const FiniteElement &el,
389  Vector &elvect);
390 
392 };
393 
394 
395 /** Class for boundary integration of the linear form:
396  (alpha/2) < (u.n) f, w > - beta < |u.n| f, w >,
397  where f and u are given scalar and vector coefficients, respectively,
398  and w is the scalar test function. */
400 {
401 private:
402  Coefficient *f;
404  double alpha, beta;
405 
406  Vector shape;
407 
408 public:
410  double a, double b)
411  { f = &_f; u = &_u; alpha = a; beta = b; }
412 
413  virtual void AssembleRHSElementVect(const FiniteElement &el,
415  Vector &elvect);
416  virtual void AssembleRHSElementVect(const FiniteElement &el,
418  Vector &elvect);
419 };
420 
421 /** Boundary linear integrator for imposing non-zero Dirichlet boundary
422  conditions, to be used in conjunction with DGDiffusionIntegrator.
423  Specifically, given the Dirichlet data u_D, the linear form assembles the
424  following integrals on the boundary:
425 
426  sigma < u_D, (Q grad(v)).n > + kappa < {h^{-1} Q} u_D, v >,
427 
428  where Q is a scalar or matrix diffusion coefficient and v is the test
429  function. The parameters sigma and kappa should be the same as the ones
430  used in the DGDiffusionIntegrator. */
432 {
433 protected:
436  double sigma, kappa;
437 
438  // these are not thread-safe!
441 
442 public:
443  DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
444  : uD(&u), Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
446  const double s, const double k)
447  : uD(&u), Q(&q), MQ(NULL), sigma(s), kappa(k) { }
449  const double s, const double k)
450  : uD(&u), Q(NULL), MQ(&q), sigma(s), kappa(k) { }
451 
452  virtual void AssembleRHSElementVect(const FiniteElement &el,
454  Vector &elvect);
455  virtual void AssembleRHSElementVect(const FiniteElement &el,
457  Vector &elvect);
458 };
459 
460 
461 /** Boundary linear form integrator for imposing non-zero Dirichlet boundary
462  conditions, in a DG elasticity formulation. Specifically, the linear form is
463  given by
464 
465  alpha < u_D, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
466  + kappa < h^{-1} (lambda + 2 mu) u_D, v >,
467 
468  where u_D is the given Dirichlet data. The parameters alpha, kappa, lambda
469  and mu, should match the parameters with the same names used in the bilinear
470  form integrator, DGElasticityIntegrator. */
472 {
473 protected:
476  double alpha, kappa;
477 
478 #ifndef MFEM_THREAD_SAFE
487 #endif
488 
489 public:
491  Coefficient &lambda_, Coefficient &mu_,
492  double alpha_, double kappa_)
493  : uD(uD_), lambda(&lambda_), mu(&mu_), alpha(alpha_), kappa(kappa_) { }
494 
495  virtual void AssembleRHSElementVect(const FiniteElement &el,
497  Vector &elvect);
498  virtual void AssembleRHSElementVect(const FiniteElement &el,
500  Vector &elvect);
501 };
502 
503 /** Class for domain integration of L(v) := (f, v), where
504  f=(f1,...,fn) and v=(v1,...,vn). that makes use of
505  VectorQuadratureFunctionCoefficient*/
507 {
508 private:
510 
511 public:
513  const IntegrationRule *ir)
514  : LinearFormIntegrator(ir), vqfc(vqfc)
515  {
516  if (ir)
517  {
518  MFEM_WARNING("Integration rule not used in this class. "
519  "The QuadratureFunction integration rules are used instead");
520  }
521  }
522 
524  virtual void AssembleRHSElementVect(const FiniteElement &fe,
526  Vector &elvect);
527 
528  virtual void SetIntRule(const IntegrationRule *ir)
529  {
530  MFEM_WARNING("Integration rule not used in this class. "
531  "The QuadratureFunction integration rules are used instead");
532  }
533 };
534 
535 /** Class for domain integration L(v) := (f, v) that makes use
536  of QuadratureFunctionCoefficient. */
538 {
539 private:
541 
542 public:
544  const IntegrationRule *ir)
545  : LinearFormIntegrator(ir), qfc(qfc)
546  {
547  if (ir)
548  {
549  MFEM_WARNING("Integration rule not used in this class. "
550  "The QuadratureFunction integration rules are used instead");
551  }
552  }
553 
555  virtual void AssembleRHSElementVect(const FiniteElement &fe,
557  Vector &elvect);
558 
559  virtual void SetIntRule(const IntegrationRule *ir)
560  {
561  MFEM_WARNING("Integration rule not used in this class. "
562  "The QuadratureFunction integration rules are used instead");
563  }
564 };
565 
566 }
567 
568 #endif
Abstract class for all finite elements.
Definition: fe.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:98
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:554
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
void GetDeltaCenter(Vector &center)
VectorFEDomainLFCurlIntegrator(VectorCoefficient &F)
Constructs the domain integrator (Q, curl v)
Definition: lininteg.hpp:292
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:512
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:484
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:620
VectorBoundaryLFIntegrator(VectorCoefficient &QG)
Constructs a boundary integrator with a given VectorCoefficient QG.
Definition: lininteg.hpp:243
Class for domain integrator L(v) := (f, grad v)
Definition: lininteg.hpp:123
BoundaryNormalLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:179
DGElasticityDirichletLFIntegrator(VectorCoefficient &uD_, Coefficient &lambda_, Coefficient &mu_, double alpha_, double kappa_)
Definition: lininteg.hpp:490
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:777
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:66
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:57
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
DeltaLFIntegrator(VectorCoefficient &vq, const IntegrationRule *ir=NULL)
This constructor should be used by derived classes that use a VectorDeltaCoefficient.
Definition: lininteg.hpp:62
virtual ~LinearFormIntegrator()
Definition: lininteg.hpp:42
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:40
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:249
Class for boundary integration L(v) := (g, v)
Definition: lininteg.hpp:150
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:39
const IntegrationRule * IntRule
Definition: lininteg.hpp:25
DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
Definition: lininteg.hpp:443
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:114
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:22
void GetDeltaCenter(Vector &center)
Returns the center of the delta coefficient.
Definition: lininteg.hpp:73
LinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:27
VectorBoundaryFluxLFIntegrator(Coefficient &f, double s=1.0, const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:341
DomainLFGradIntegrator(VectorCoefficient &QF)
Constructs the domain integrator (Q, grad v)
Definition: lininteg.hpp:132
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
Class for boundary integration .
Definition: lininteg.hpp:172
DomainLFIntegrator(Coefficient &QF, int a=2, int b=0)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:100
double b
Definition: lissajous.cpp:42
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:308
VectorFEBoundaryFluxLFIntegrator(int a=1, int b=-1)
Definition: lininteg.hpp:363
VectorFEBoundaryTangentLFIntegrator(VectorCoefficient &QG, int a=2, int b=0)
Definition: lininteg.hpp:383
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:158
VectorDeltaCoefficient * vec_delta
Definition: lininteg.hpp:51
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:467
MatrixCoefficient * MQ
Definition: lininteg.hpp:435
VectorFEBoundaryFluxLFIntegrator(Coefficient &f, int a=2, int b=0)
Definition: lininteg.hpp:365
bool IsDelta() const
Returns true if the derived class instance uses a delta coefficient.
Definition: lininteg.hpp:70
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:528
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:209
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:686
DomainLFIntegrator(Coefficient &QF, const IntegrationRule *ir)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:106
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
Definition: coefficient.cpp:77
VectorFEDomainLFDivIntegrator(Coefficient &QF)
Constructs the domain integrator (Q, div v)
Definition: lininteg.hpp:314
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:55
QuadratureLFIntegrator(QuadratureFunctionCoefficient &qfc, const IntegrationRule *ir)
Definition: lininteg.hpp:543
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:50
Base class for Matrix Coefficients that optionally depend on time and space.
Abstract class for integrators that support delta coefficients.
Definition: lininteg.hpp:47
VectorDomainLFIntegrator(VectorCoefficient &QF)
Constructs a domain integrator with a given VectorCoefficient.
Definition: lininteg.hpp:217
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:521
VectorFEDomainLFIntegrator(VectorCoefficient &F)
Definition: lininteg.hpp:268
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:586
double a
Definition: lissajous.cpp:41
for Nedelec Elements)
Definition: lininteg.hpp:283
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:419
BoundaryTangentialLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:197
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:559
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:913
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:177
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:435
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:260
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:943
Vector data type.
Definition: vector.hpp:51
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:386
Class for boundary integration .
Definition: lininteg.hpp:376
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:291
DGDirichletLFIntegrator(Coefficient &u, MatrixCoefficient &q, const double s, const double k)
Definition: lininteg.hpp:448
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:26
Class for boundary integration in 2D.
Definition: lininteg.hpp:190
VectorQuadratureLFIntegrator(VectorQuadratureFunctionCoefficient &vqfc, const IntegrationRule *ir)
Definition: lininteg.hpp:512
DGDirichletLFIntegrator(Coefficient &u, Coefficient &q, const double s, const double k)
Definition: lininteg.hpp:445
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...
BoundaryFlowIntegrator(Coefficient &_f, VectorCoefficient &_u, double a, double b)
Definition: lininteg.hpp:409
double f(const Vector &p)