MFEM  v3.4
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, 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_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  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 boundary integration L(v) := (g, v)
124 {
125  Vector shape;
126  Coefficient &Q;
127  int oa, ob;
128 public:
129  /// Constructs a boundary integrator with a given Coefficient QG
130  BoundaryLFIntegrator(Coefficient &QG, int a = 1, int b = 1)
131  : Q(QG), oa(a), ob(b) { }
132 
133  /** Given a particular boundary Finite Element and a transformation (Tr)
134  computes the element boundary vector, elvect. */
135  virtual void AssembleRHSElementVect(const FiniteElement &el,
137  Vector &elvect);
138 
140 };
141 
142 /// Class for boundary integration \f$ L(v) = (g \cdot n, v) \f$
144 {
145  Vector shape;
147  int oa, ob;
148 public:
149  /// Constructs a boundary integrator with a given Coefficient QG
151  : Q(QG), oa(a), ob(b) { }
152 
153  virtual void AssembleRHSElementVect(const FiniteElement &el,
155  Vector &elvect);
156 
158 };
159 
160 /// Class for boundary integration \f$ L(v) = (g \cdot \tau, v) \f$ in 2D
162 {
163  Vector shape;
165  int oa, ob;
166 public:
167  /// Constructs a boundary integrator with a given Coefficient QG
169  : Q(QG), oa(a), ob(b) { }
170 
171  virtual void AssembleRHSElementVect(const FiniteElement &el,
173  Vector &elvect);
174 
176 };
177 
178 /** Class for domain integration of L(v) := (f, v), where
179  f=(f1,...,fn) and v=(v1,...,vn). */
181 {
182 private:
183  Vector shape, Qvec;
185 
186 public:
187  /// Constructs a domain integrator with a given VectorCoefficient
189  : DeltaLFIntegrator(QF), Q(QF) { }
190 
191  /** Given a particular Finite Element and a transformation (Tr)
192  computes the element right hand side element vector, elvect. */
193  virtual void AssembleRHSElementVect(const FiniteElement &el,
195  Vector &elvect);
196 
197  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
199  Vector &elvect);
200 
202 };
203 
204 /** Class for boundary integration of L(v) := (g, v), where
205  f=(f1,...,fn) and v=(v1,...,vn). */
207 {
208 private:
209  Vector shape, vec;
211 
212 public:
213  /// Constructs a boundary integrator with a given VectorCoefficient QG
215 
216  /** Given a particular boundary Finite Element and a transformation (Tr)
217  computes the element boundary vector, elvect. */
218  virtual void AssembleRHSElementVect(const FiniteElement &el,
220  Vector &elvect);
221 
222  // For DG spaces
223  virtual void AssembleRHSElementVect(const FiniteElement &el,
225  Vector &elvect);
226 
228 };
229 
230 /// \f$ (f, v)_{\Omega} \f$ for VectorFiniteElements (Nedelec, Raviart-Thomas)
232 {
233 private:
234  VectorCoefficient &QF;
235  DenseMatrix vshape;
236  Vector vec;
237 
238 public:
240  : DeltaLFIntegrator(F), QF(F) { }
241 
242  virtual void AssembleRHSElementVect(const FiniteElement &el,
244  Vector &elvect);
245 
246  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
248  Vector &elvect);
249 
251 };
252 
253 
254 /** \f$ (f, v \cdot n)_{\partial\Omega} \f$ for vector test function
255  v=(v1,...,vn) where all vi are in the same scalar FE space and f is a
256  scalar function. */
258 {
259 private:
260  double Sign;
261  Coefficient *F;
262  Vector shape, nor;
263 
264 public:
266  const IntegrationRule *ir = NULL)
267  : LinearFormIntegrator(ir), Sign(s), F(&f) { }
268 
269  virtual void AssembleRHSElementVect(const FiniteElement &el,
271  Vector &elvect);
272 
274 };
275 
276 /** Class for boundary integration of (f, v.n) for scalar coefficient f and
277  RT vector test function v. This integrator works with RT spaces defined
278  using the RT_FECollection class. */
280 {
281 private:
282  Coefficient &F;
283  Vector shape;
284 
285 public:
287 
288  virtual void AssembleRHSElementVect(const FiniteElement &el,
290  Vector &elvect);
291 
293 };
294 
295 /// Class for boundary integration \f$ L(v) = (n \times f, v) \f$
297 {
298 private:
300 
301 public:
303 
304  virtual void AssembleRHSElementVect(const FiniteElement &el,
306  Vector &elvect);
307 
309 };
310 
311 
312 /** Class for boundary integration of the linear form:
313  (alpha/2) < (u.n) f, w > - beta < |u.n| f, w >,
314  where f and u are given scalar and vector coefficients, respectively,
315  and w is the scalar test function. */
317 {
318 private:
319  Coefficient *f;
321  double alpha, beta;
322 
323  Vector shape;
324 
325 public:
327  double a, double b)
328  { f = &_f; u = &_u; alpha = a; beta = b; }
329 
330  virtual void AssembleRHSElementVect(const FiniteElement &el,
332  Vector &elvect);
333  virtual void AssembleRHSElementVect(const FiniteElement &el,
335  Vector &elvect);
336 };
337 
338 /** Boundary linear integrator for imposing non-zero Dirichlet boundary
339  conditions, to be used in conjunction with DGDiffusionIntegrator.
340  Specifically, given the Dirichlet data u_D, the linear form assembles the
341  following integrals on the boundary:
342 
343  sigma < u_D, (Q grad(v)).n > + kappa < {h^{-1} Q} u_D, v >,
344 
345  where Q is a scalar or matrix diffusion coefficient and v is the test
346  function. The parameters sigma and kappa should be the same as the ones
347  used in the DGDiffusionIntegrator. */
349 {
350 protected:
353  double sigma, kappa;
354 
355  // these are not thread-safe!
358 
359 public:
360  DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
361  : uD(&u), Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
363  const double s, const double k)
364  : uD(&u), Q(&q), MQ(NULL), sigma(s), kappa(k) { }
366  const double s, const double k)
367  : uD(&u), Q(NULL), MQ(&q), sigma(s), kappa(k) { }
368 
369  virtual void AssembleRHSElementVect(const FiniteElement &el,
371  Vector &elvect);
372  virtual void AssembleRHSElementVect(const FiniteElement &el,
374  Vector &elvect);
375 };
376 
377 
378 /** Boundary linear form integrator for imposing non-zero Dirichlet boundary
379  conditions, in a DG elasticity formulation. Specifically, the linear form is
380  given by
381 
382  alpha < u_D, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
383  + kappa < h^{-1} (lambda + 2 mu) u_D, v >,
384 
385  where u_D is the given Dirichlet data. The parameters alpha, kappa, lambda
386  and mu, should match the parameters with the same names used in the bilinear
387  form integrator, DGElasticityIntegrator. */
389 {
390 protected:
393  double alpha, kappa;
394 
395 #ifndef MFEM_THREAD_SAFE
404 #endif
405 
406 public:
408  Coefficient &lambda_, Coefficient &mu_,
409  double alpha_, double kappa_)
410  : uD(uD_), lambda(&lambda_), mu(&mu_), alpha(alpha_), kappa(kappa_) { }
411 
412  virtual void AssembleRHSElementVect(const FiniteElement &el,
414  Vector &elvect);
415  virtual void AssembleRHSElementVect(const FiniteElement &el,
417  Vector &elvect);
418 };
419 
420 }
421 
422 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:140
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:387
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
void GetDeltaCenter(Vector &center)
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:328
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:452
VectorBoundaryLFIntegrator(VectorCoefficient &QG)
Constructs a boundary integrator with a given VectorCoefficient QG.
Definition: lininteg.hpp:214
BoundaryNormalLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:150
DGElasticityDirichletLFIntegrator(VectorCoefficient &uD_, Coefficient &lambda_, Coefficient &mu_, double alpha_, double kappa_)
Definition: lininteg.hpp:407
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:604
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
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.
const IntegrationRule * GetIntRule()
Definition: lininteg.hpp:40
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:168
Class for boundary integration L(v) := (g, v)
Definition: lininteg.hpp:123
void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:39
VectorFEBoundaryTangentLFIntegrator(VectorCoefficient &QG)
Definition: lininteg.hpp:302
const IntegrationRule * IntRule
Definition: lininteg.hpp:25
DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
Definition: lininteg.hpp:360
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:67
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:265
Class for boundary integration .
Definition: lininteg.hpp:143
DomainLFIntegrator(Coefficient &QF, int a=2, int b=0)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:100
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:228
VectorFEBoundaryFluxLFIntegrator(Coefficient &f)
Definition: lininteg.hpp:286
BoundaryLFIntegrator(Coefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:130
VectorDeltaCoefficient * vec_delta
Definition: lininteg.hpp:51
MatrixCoefficient * MQ
Definition: lininteg.hpp:352
bool IsDelta() const
Returns true if the derived class instance uses a delta coefficient.
Definition: lininteg.hpp:70
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:128
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:515
DomainLFIntegrator(Coefficient &QF, const IntegrationRule *ir)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:106
void GetDeltaCenter(Vector &center)
Definition: coefficient.cpp:77
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
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
DeltaCoefficient * delta
Definition: lininteg.hpp:50
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:188
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:354
VectorFEDomainLFIntegrator(VectorCoefficient &F)
Definition: lininteg.hpp:239
VectorDeltaCoefficient: DeltaCoefficient with a direction.
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:417
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:337
BoundaryTangentialLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:168
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:96
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:231
Vector data type.
Definition: vector.hpp:48
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:303
Class for boundary integration .
Definition: lininteg.hpp:296
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:210
DGDirichletLFIntegrator(Coefficient &u, MatrixCoefficient &q, const double s, const double k)
Definition: lininteg.hpp:365
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:26
Class for boundary integration in 2D.
Definition: lininteg.hpp:161
DGDirichletLFIntegrator(Coefficient &u, Coefficient &q, const double s, const double k)
Definition: lininteg.hpp:362
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:326