MFEM  v3.2
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 
23 {
24 protected:
26 
28  { IntRule = ir; }
29 
30 public:
33  virtual void AssembleRHSElementVect(const FiniteElement &el,
35  Vector &elvect) = 0;
36  virtual void AssembleRHSElementVect(const FiniteElement &el,
38  Vector &elvect);
39 
40  void SetIntRule(const IntegrationRule *ir) { IntRule = ir; }
41 
42  virtual ~LinearFormIntegrator() { }
43 };
44 
45 
48 {
49  Vector shape;
50  Coefficient &Q;
51  int oa, ob;
52 public:
54  DomainLFIntegrator(Coefficient &QF, int a = 2, int b = 0)
55  // the old default was a = 1, b = 1
56  // for simple elliptic problems a = 2, b = -2 is ok
57  : Q(QF), oa(a), ob(b) { }
58 
61  : LinearFormIntegrator(ir), Q(QF), oa(1), ob(1) { }
62 
65  virtual void AssembleRHSElementVect(const FiniteElement &el,
67  Vector &elvect);
68 
70 };
71 
74 {
75  Vector shape;
76  Coefficient &Q;
77  int oa, ob;
78 public:
80  BoundaryLFIntegrator(Coefficient &QG, int a = 1, int b = 1)
81  : Q(QG), oa(a), ob(b) {};
82 
85  virtual void AssembleRHSElementVect(const FiniteElement &el,
87  Vector &elvect);
88 
90 };
91 
94 {
95  Vector shape;
97  int oa, ob;
98 public:
101  : Q(QG), oa(a), ob(b) {};
102 
103  virtual void AssembleRHSElementVect(const FiniteElement &el,
105  Vector &elvect);
106 
108 };
109 
112 {
113  Vector shape;
115  int oa, ob;
116 public:
119  : Q(QG), oa(a), ob(b) {};
120 
121  virtual void AssembleRHSElementVect(const FiniteElement &el,
123  Vector &elvect);
124 
126 };
127 
131 {
132 private:
133  Vector shape, Qvec;
135 
136 public:
139 
142  virtual void AssembleRHSElementVect(const FiniteElement &el,
144  Vector &elvect);
145 
147 };
148 
152 {
153 private:
154  Vector shape, vec;
156 
157 public:
160 
163  virtual void AssembleRHSElementVect(const FiniteElement &el,
165  Vector &elvect);
166 
168 };
169 
172 {
173 private:
174  VectorCoefficient &QF;
175  DenseMatrix vshape;
176  Vector vec;
177 
178 public:
180 
181  virtual void AssembleRHSElementVect(const FiniteElement &el,
183  Vector &elvect);
184 
186 };
187 
188 
193 {
194 private:
195  double Sign;
196  Coefficient *F;
197  Vector shape, nor;
198 
199 public:
201  const IntegrationRule *ir = NULL)
202  : LinearFormIntegrator(ir), Sign(s), F(&f) { }
203 
204  virtual void AssembleRHSElementVect(const FiniteElement &el,
206  Vector &elvect);
207 
209 };
210 
215 {
216 private:
217  Coefficient &F;
218  Vector shape;
219 
220 public:
222 
223  virtual void AssembleRHSElementVect(const FiniteElement &el,
225  Vector &elvect);
226 
228 };
229 
232 {
233 private:
235 
236 public:
238 
239  virtual void AssembleRHSElementVect(const FiniteElement &el,
241  Vector &elvect);
242 
244 };
245 
246 
252 {
253 private:
254  Coefficient *f;
256  double alpha, beta;
257 
258  Vector shape;
259 
260 public:
262  double a, double b)
263  { f = &_f; u = &_u; alpha = a; beta = b; }
264 
265  virtual void AssembleRHSElementVect(const FiniteElement &el,
267  Vector &elvect);
268  virtual void AssembleRHSElementVect(const FiniteElement &el,
270  Vector &elvect);
271 };
272 
284 {
285 protected:
288  double sigma, kappa;
289 
290  // these are not thread-safe!
293 
294 public:
295  DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
296  : uD(&u), Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
298  const double s, const double k)
299  : uD(&u), Q(&q), MQ(NULL), sigma(s), kappa(k) { }
301  const double s, const double k)
302  : uD(&u), Q(NULL), MQ(&q), sigma(s), kappa(k) { }
303 
304  virtual void AssembleRHSElementVect(const FiniteElement &el,
306  Vector &elvect);
307  virtual void AssembleRHSElementVect(const FiniteElement &el,
309  Vector &elvect);
310 };
311 
312 }
313 
314 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:44
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:303
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:47
Class for integration rule.
Definition: intrules.hpp:83
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:368
VectorBoundaryLFIntegrator(VectorCoefficient &QG)
Constructs a boundary integrator with a given VectorCoefficient QG.
Definition: lininteg.hpp:159
BoundaryNormalLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:100
virtual ~LinearFormIntegrator()
Definition: lininteg.hpp:42
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:158
Class for boundary integration L(v) := (g, v)
Definition: lininteg.hpp:73
void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:40
VectorFEBoundaryTangentLFIntegrator(VectorCoefficient &QG)
Definition: lininteg.hpp:237
const IntegrationRule * IntRule
Definition: lininteg.hpp:25
DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
Definition: lininteg.hpp:295
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:57
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:22
LinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:27
VectorBoundaryFluxLFIntegrator(Coefficient &f, double s=1.0, const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:200
Class for boundary integration .
Definition: lininteg.hpp:93
DomainLFIntegrator(Coefficient &QF, int a=2, int b=0)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:54
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:200
VectorFEBoundaryFluxLFIntegrator(Coefficient &f)
Definition: lininteg.hpp:221
BoundaryLFIntegrator(Coefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:80
MatrixCoefficient * MQ
Definition: lininteg.hpp:287
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:118
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:431
DomainLFIntegrator(Coefficient &QF, const IntegrationRule *ir)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:60
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)=0
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
VectorDomainLFIntegrator(VectorCoefficient &QF)
Constructs a domain integrator with a given VectorCoefficient.
Definition: lininteg.hpp:138
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:270
VectorFEDomainLFIntegrator(VectorCoefficient &F)
Definition: lininteg.hpp:179
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:333
BoundaryTangentialLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:118
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:86
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:171
Vector data type.
Definition: vector.hpp:33
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:235
Class for boundary integration .
Definition: lininteg.hpp:231
DGDirichletLFIntegrator(Coefficient &u, MatrixCoefficient &q, const double s, const double k)
Definition: lininteg.hpp:300
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:26
Class for boundary integration in 2D.
Definition: lininteg.hpp:111
DGDirichletLFIntegrator(Coefficient &u, Coefficient &q, const double s, const double k)
Definition: lininteg.hpp:297
BoundaryFlowIntegrator(Coefficient &_f, VectorCoefficient &_u, double a, double b)
Definition: lininteg.hpp:261