MFEM  v4.6.0
Finite element discretization library
lininteg.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 #include "bilininteg.hpp"
18 #include <random>
19 
20 namespace mfem
21 {
22 
23 /// Abstract base class LinearFormIntegrator
25 {
26 protected:
28 
29  LinearFormIntegrator(const IntegrationRule *ir = NULL) { IntRule = ir; }
30 
31 public:
32 
33  /// Method probing for assembly on device
34  virtual bool SupportsDevice() const { return false; }
35 
36  /// Method defining assembly on device
37  virtual void AssembleDevice(const FiniteElementSpace &fes,
38  const Array<int> &markers,
39  Vector &b);
40 
41  /** Given a particular Finite Element and a transformation (Tr)
42  computes the element vector, elvect. */
43  virtual void AssembleRHSElementVect(const FiniteElement &el,
45  Vector &elvect) = 0;
46  virtual void AssembleRHSElementVect(const FiniteElement &el,
48  Vector &elvect);
49  virtual void AssembleRHSElementVect(const FiniteElement &el1,
50  const FiniteElement &el2,
52  Vector &elvect);
53 
54  virtual void SetIntRule(const IntegrationRule *ir) { IntRule = ir; }
55  const IntegrationRule* GetIntRule() { return IntRule; }
56 
57  virtual ~LinearFormIntegrator() { }
58 };
59 
60 
61 /// Abstract class for integrators that support delta coefficients
63 {
64 protected:
67 
68  /** @brief This constructor should be used by derived classes that use a
69  scalar DeltaCoefficient. */
72  delta(dynamic_cast<DeltaCoefficient*>(&q)),
73  vec_delta(NULL) { }
74 
75  /** @brief This constructor should be used by derived classes that use a
76  VectorDeltaCoefficient. */
78  const IntegrationRule *ir = NULL)
80  delta(NULL),
81  vec_delta(dynamic_cast<VectorDeltaCoefficient*>(&vq)) { }
82 
83 public:
84  /// Returns true if the derived class instance uses a delta coefficient.
85  bool IsDelta() const { return (delta || vec_delta); }
86 
87  /// Returns the center of the delta coefficient.
88  void GetDeltaCenter(Vector &center)
89  {
90  if (delta) { delta->GetDeltaCenter(center); return; }
91  if (vec_delta) { vec_delta->GetDeltaCenter(center); return; }
92  center.SetSize(0);
93  }
94 
95  /** @brief Assemble the delta coefficient at the IntegrationPoint set in
96  @a Trans which is assumed to map to the delta coefficient center.
97 
98  @note This method should be called for one mesh element only, including
99  in parallel, even when the center of the delta coefficient is shared by
100  multiple elements. */
101  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
102  ElementTransformation &Trans,
103  Vector &elvect) = 0;
104 };
105 
106 
107 /// Class for domain integration L(v) := (f, v)
109 {
110  Vector shape;
111  Coefficient &Q;
112  int oa, ob;
113 public:
114  /// Constructs a domain integrator with a given Coefficient
115  DomainLFIntegrator(Coefficient &QF, int a = 2, int b = 0)
116  // the old default was a = 1, b = 1
117  // for simple elliptic problems a = 2, b = -2 is OK
118  : DeltaLFIntegrator(QF), Q(QF), oa(a), ob(b) { }
119 
120  /// Constructs a domain integrator with a given Coefficient
122  : DeltaLFIntegrator(QF, ir), Q(QF), oa(1), ob(1) { }
123 
124  virtual bool SupportsDevice() const { return true; }
125 
126  /// Method defining assembly on device
127  virtual void AssembleDevice(const FiniteElementSpace &fes,
128  const Array<int> &markers,
129  Vector &b);
130 
131  /** Given a particular Finite Element and a transformation (Tr)
132  computes the element right hand side element vector, elvect. */
133  virtual void AssembleRHSElementVect(const FiniteElement &el,
135  Vector &elvect);
136 
137  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
138  ElementTransformation &Trans,
139  Vector &elvect);
140 
142 };
143 
144 /// Class for domain integrator L(v) := (f, grad v)
146 {
147 private:
148  Vector shape, Qvec;
150  DenseMatrix dshape;
151 
152 public:
153  /// Constructs the domain integrator (Q, grad v)
155  : DeltaLFIntegrator(QF), Q(QF) { }
156 
157  virtual bool SupportsDevice() const { return true; }
158 
159  /// Method defining assembly on device
160  virtual void AssembleDevice(const FiniteElementSpace &fes,
161  const Array<int> &markers,
162  Vector &b);
163 
164  /** Given a particular Finite Element and a transformation (Tr)
165  computes the element right hand side element vector, elvect. */
166  virtual void AssembleRHSElementVect(const FiniteElement &el,
168  Vector &elvect);
169 
170  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
171  ElementTransformation &Trans,
172  Vector &elvect);
173 
175 };
176 
177 
178 /// Class for boundary integration L(v) := (g, v)
180 {
181  Vector shape;
182  Coefficient &Q;
183  int oa, ob;
184 public:
185  /** @brief Constructs a boundary integrator with a given Coefficient @a QG.
186  Integration order will be @a a * basis_order + @a b. */
187  BoundaryLFIntegrator(Coefficient &QG, int a = 1, int b = 1)
188  : Q(QG), oa(a), ob(b) { }
189 
190  virtual bool SupportsDevice() const { return true; }
191 
192  /// Method defining assembly on device
193  virtual void AssembleDevice(const FiniteElementSpace &fes,
194  const Array<int> &markers,
195  Vector &b);
196 
197  /** Given a particular boundary Finite Element and a transformation (Tr)
198  computes the element boundary vector, elvect. */
199  virtual void AssembleRHSElementVect(const FiniteElement &el,
201  Vector &elvect);
202  virtual void AssembleRHSElementVect(const FiniteElement &el,
204  Vector &elvect);
205 
207 };
208 
209 /// Class for boundary integration \f$ L(v) = (g \cdot n, v) \f$
211 {
212  Vector shape;
214  int oa, ob;
215 public:
216  /// Constructs a boundary integrator with a given Coefficient QG
218  : Q(QG), oa(a), ob(b) { }
219 
220  virtual bool SupportsDevice() const { return true; }
221 
222  /// Method defining assembly on device
223  virtual void AssembleDevice(const FiniteElementSpace &fes,
224  const Array<int> &markers,
225  Vector &b);
226 
227  virtual void AssembleRHSElementVect(const FiniteElement &el,
229  Vector &elvect);
230 
232 };
233 
234 /// Class for boundary integration \f$ L(v) = (g \cdot \tau, v) \f$ in 2D
236 {
237  Vector shape;
239  int oa, ob;
240 public:
241  /// Constructs a boundary integrator with a given Coefficient QG
243  : Q(QG), oa(a), ob(b) { }
244 
245  virtual void AssembleRHSElementVect(const FiniteElement &el,
247  Vector &elvect);
248 
250 };
251 
252 /** Class for domain integration of L(v) := (f, v), where
253  f=(f1,...,fn) and v=(v1,...,vn). */
255 {
256 private:
257  Vector shape, Qvec;
259 
260 public:
261  /// Constructs a domain integrator with a given VectorCoefficient
263  : DeltaLFIntegrator(QF), Q(QF) { }
264 
265  virtual bool SupportsDevice() const { return true; }
266 
267  /// Method defining assembly on device
268  virtual void AssembleDevice(const FiniteElementSpace &fes,
269  const Array<int> &markers,
270  Vector &b);
271 
272  /** Given a particular Finite Element and a transformation (Tr)
273  computes the element right hand side element vector, elvect. */
274  virtual void AssembleRHSElementVect(const FiniteElement &el,
276  Vector &elvect);
277 
278  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
279  ElementTransformation &Trans,
280  Vector &elvect);
281 
283 };
284 
285 /** Class for domain integrator L(v) := (f, grad v), where
286  f=(f1x,f1y,f1z,...,fnx,fny,fnz) and v=(v1,...,vn). */
288 {
289 private:
290  Vector shape, Qvec;
292  DenseMatrix dshape;
293 
294 public:
295  /// Constructs the domain integrator (Q, grad v)
297  : DeltaLFIntegrator(QF), Q(QF) { }
298 
299  virtual bool SupportsDevice() const override { return true; }
300 
301  /// Method defining assembly on device
302  virtual void AssembleDevice(const FiniteElementSpace &fes,
303  const Array<int> &markers,
304  Vector &b) override;
305 
306  /** Given a particular Finite Element and a transformation (Tr)
307  computes the element right hand side element vector, elvect. */
308  virtual void AssembleRHSElementVect(const FiniteElement &el,
310  Vector &elvect) override;
311 
312  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
313  ElementTransformation &Trans,
314  Vector &elvect) override;
315 
317 };
318 
319 /** Class for boundary integration of L(v) := (g, v), where
320  f=(f1,...,fn) and v=(v1,...,vn). */
322 {
323 private:
324  Vector shape, vec;
326 
327 public:
328  /// Constructs a boundary integrator with a given VectorCoefficient QG
330 
331  /** Given a particular boundary Finite Element and a transformation (Tr)
332  computes the element boundary vector, elvect. */
333  virtual void AssembleRHSElementVect(const FiniteElement &el,
335  Vector &elvect);
336 
337  // For DG spaces
338  virtual void AssembleRHSElementVect(const FiniteElement &el,
340  Vector &elvect);
341 
343 };
344 
345 /// \f$ (f, v)_{\Omega} \f$ for VectorFiniteElements (Nedelec, Raviart-Thomas)
347 {
348 private:
349  VectorCoefficient &QF;
350  DenseMatrix vshape;
351  Vector vec;
352 
353 public:
355  : DeltaLFIntegrator(F), QF(F) { }
356 
357  virtual void AssembleRHSElementVect(const FiniteElement &el,
359  Vector &elvect);
360 
361  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
362  ElementTransformation &Trans,
363  Vector &elvect);
364 
365  virtual bool SupportsDevice() const { return true; }
366 
367  virtual void AssembleDevice(const FiniteElementSpace &fes,
368  const Array<int> &markers,
369  Vector &b);
370 
372 };
373 
374 /// \f$ (Q, curl v)_{\Omega} \f$ for Nedelec Elements)
376 {
377 private:
378  VectorCoefficient *QF=nullptr;
379  DenseMatrix curlshape;
380  Vector vec;
381 
382 public:
383  /// Constructs the domain integrator (Q, curl v)
385  : DeltaLFIntegrator(F), QF(&F) { }
386 
387  virtual void AssembleRHSElementVect(const FiniteElement &el,
389  Vector &elvect);
390 
391  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
392  ElementTransformation &Trans,
393  Vector &elvect);
394 
396 };
397 
398 /// \f$ (Q, div v)_{\Omega} \f$ for RT Elements)
400 {
401 private:
402  Vector divshape;
403  Coefficient &Q;
404 public:
405  /// Constructs the domain integrator (Q, div v)
407  : DeltaLFIntegrator(QF), Q(QF) { }
408 
409  /** Given a particular Finite Element and a transformation (Tr)
410  computes the element right hand side element vector, elvect. */
411  virtual void AssembleRHSElementVect(const FiniteElement &el,
413  Vector &elvect);
414 
415  virtual void AssembleDeltaElementVect(const FiniteElement &fe,
416  ElementTransformation &Trans,
417  Vector &elvect);
418 
420 };
421 
422 /** \f$ (f, v \cdot n)_{\partial\Omega} \f$ for vector test function
423  v=(v1,...,vn) where all vi are in the same scalar FE space and f is a
424  scalar function. */
426 {
427 private:
428  double Sign;
429  Coefficient *F;
430  Vector shape, nor;
431 
432 public:
434  const IntegrationRule *ir = NULL)
435  : LinearFormIntegrator(ir), Sign(s), F(&f) { }
436 
437  virtual void AssembleRHSElementVect(const FiniteElement &el,
439  Vector &elvect);
440 
442 };
443 
444 /** Class for boundary integration of (f, v.n) for scalar coefficient f and
445  RT vector test function v. This integrator works with RT spaces defined
446  using the RT_FECollection class. */
448 {
449 private:
450  Coefficient *F;
451  Vector shape;
452  int oa, ob; // these control the quadrature order, see DomainLFIntegrator
453 
454 public:
456  : F(NULL), oa(a), ob(b) { }
458  : F(&f), oa(a), ob(b) { }
459 
460  virtual void AssembleRHSElementVect(const FiniteElement &el,
462  Vector &elvect);
463 
465 
466  virtual bool SupportsDevice() const { return true; }
467 
468  virtual void AssembleDevice(const FiniteElementSpace &fes,
469  const Array<int> &markers,
470  Vector &b);
471 };
472 
473 /// Class for boundary integration \f$ L(v) = (n \times f, v) \f$
475 {
476 private:
478  int oa, ob;
479 
480 public:
482  int a = 2, int b = 0)
483  : f(QG), oa(a), ob(b) { }
484 
485  virtual void AssembleRHSElementVect(const FiniteElement &el,
487  Vector &elvect);
488 
490 };
491 
492 
493 /** Class for boundary integration of the linear form:
494  (alpha/2) < (u.n) f, w > - beta < |u.n| f, w >,
495  where f and u are given scalar and vector coefficients, respectively,
496  and w is the scalar test function. */
498 {
499 private:
500  Coefficient *f;
502  double alpha, beta;
503 
504  Vector shape;
505 
506 public:
508  double a)
509  { f = &f_; u = &u_; alpha = a; beta = 0.5*a; }
510 
512  double a, double b)
513  { f = &f_; u = &u_; alpha = a; beta = b; }
514 
515  virtual void AssembleRHSElementVect(const FiniteElement &el,
517  Vector &elvect);
518  virtual void AssembleRHSElementVect(const FiniteElement &el,
520  Vector &elvect);
521 
523 };
524 
525 
526 /** Boundary linear integrator for imposing non-zero Dirichlet boundary
527  conditions, to be used in conjunction with DGDiffusionIntegrator.
528  Specifically, given the Dirichlet data u_D, the linear form assembles the
529  following integrals on the boundary:
530 
531  sigma < u_D, (Q grad(v)).n > + kappa < {h^{-1} Q} u_D, v >,
532 
533  where Q is a scalar or matrix diffusion coefficient and v is the test
534  function. The parameters sigma and kappa should be the same as the ones
535  used in the DGDiffusionIntegrator. */
537 {
538 protected:
541  double sigma, kappa;
542 
543  // these are not thread-safe!
546 
547 public:
548  DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
549  : uD(&u), Q(NULL), MQ(NULL), sigma(s), kappa(k) { }
551  const double s, const double k)
552  : uD(&u), Q(&q), MQ(NULL), sigma(s), kappa(k) { }
554  const double s, const double k)
555  : uD(&u), Q(NULL), MQ(&q), sigma(s), kappa(k) { }
556 
557  virtual void AssembleRHSElementVect(const FiniteElement &el,
559  Vector &elvect);
560  virtual void AssembleRHSElementVect(const FiniteElement &el,
562  Vector &elvect);
563 
565 };
566 
567 
568 /** Boundary linear form integrator for imposing non-zero Dirichlet boundary
569  conditions, in a DG elasticity formulation. Specifically, the linear form is
570  given by
571 
572  alpha < u_D, (lambda div(v) I + mu (grad(v) + grad(v)^T)) . n > +
573  + kappa < h^{-1} (lambda + 2 mu) u_D, v >,
574 
575  where u_D is the given Dirichlet data. The parameters alpha, kappa, lambda
576  and mu, should match the parameters with the same names used in the bilinear
577  form integrator, DGElasticityIntegrator. */
579 {
580 protected:
583  double alpha, kappa;
584 
585 #ifndef MFEM_THREAD_SAFE
594 #endif
595 
596 public:
598  Coefficient &lambda_, Coefficient &mu_,
599  double alpha_, double kappa_)
600  : uD(uD_), lambda(&lambda_), mu(&mu_), alpha(alpha_), kappa(kappa_) { }
601 
602  virtual void AssembleRHSElementVect(const FiniteElement &el,
604  Vector &elvect);
605  virtual void AssembleRHSElementVect(const FiniteElement &el,
607  Vector &elvect);
608 
610 };
611 
612 
613 /** Class for spatial white Gaussian noise integration.
614 
615  The target problem is the linear SPDE a(u,v) = F(v) with F(v) := <Ẇ,v>,
616  where Ẇ is spatial white Gaussian noise. When the Galerkin method is used to
617  discretize this problem into a linear system of equations Ax = b, the RHS is
618  a Gaussian random vector b~N(0,M) whose covariance matrix is the same as the
619  mass matrix M_ij = (v_i,v_j). This property can be ensured if b = H w, where
620  HHáµ€ = M and each component w_i~N(0,1).
621 
622  There is much flexibility in how we may wish to define H. In this PR, we
623  define H = Páµ€ diag(L_e), where P is the local-to-global dof assembly matrix
624  and diag(L_e) is a block-diagonal matrix with L_e L_eáµ€ = M_e, where M_e is
625  the element mass matrix for element e. A straightforward computation shows
626  that HHáµ€ = Páµ€ diag(M_e) P = M, as necessary. */
628 {
629 #ifdef MFEM_USE_MPI
630  MPI_Comm comm;
631 #endif
632  MassIntegrator massinteg;
634 
635  // Define random generator with Gaussian distribution
636  std::default_random_engine generator;
637  std::normal_distribution<double> dist;
638 
639  bool save_factors = false;
640 public:
641 
642 #ifdef MFEM_USE_MPI
643  /** @brief Sets the @a seed_ of the random number generator. A fixed seed
644  allows for a reproducible sequence of white noise vectors. */
646  : LinearFormIntegrator(), comm(MPI_COMM_NULL)
647  {
648  if (seed_ > 0) { SetSeed(seed_); }
649  }
650 
651  /** @brief Sets the MPI communicator @a comm_ and the @a seed_ of the random
652  number generator. A fixed seed allows for a reproducible sequence of
653  white noise vectors. */
654  WhiteGaussianNoiseDomainLFIntegrator(MPI_Comm comm_, int seed_)
655  : LinearFormIntegrator(), comm(comm_)
656  {
657  int myid;
658  MPI_Comm_rank(comm, &myid);
659 
660  int seed = (seed_ > 0) ? seed_ + myid : time(0) + myid;
661  SetSeed(seed);
662  }
663 #else
664  /** @brief Sets the @a seed_ of the random number generator. A fixed seed
665  allows for a reproducible sequence of white noise vectors. */
668  {
669  if (seed_ > 0) { SetSeed(seed_); }
670  }
671 #endif
672  /// @brief Sets/resets the @a seed of the random number generator.
673  void SetSeed(int seed)
674  {
675  generator.seed(seed);
676  }
677 
679  virtual void AssembleRHSElementVect(const FiniteElement &el,
681  Vector &elvect);
682 
683  /** @brief Saves the lower triangular matrices in the element-wise Cholesky
684  decomposition. The parameter @a NE should be the number of elements in
685  the mesh. */
686  void SaveFactors(int NE)
687  {
688  save_factors = true;
689  ResetFactors(NE);
690  }
691 
692  /** @brief Resets the array of saved lower triangular Cholesky decomposition
693  matrices. The parameter @a NE should be the number of elements in the
694  mesh. */
695  void ResetFactors(int NE)
696  {
697  for (int i = 0; i<L.Size(); i++)
698  {
699  delete L[i];
700  }
701  L.DeleteAll();
702 
703  L.SetSize(NE);
704  for (int i = 0; i<NE; i++)
705  {
706  L[i] = nullptr;
707  }
708  }
709 
711  {
712  for (int i = 0; i<L.Size(); i++)
713  {
714  delete L[i];
715  }
716  L.DeleteAll();
717  }
718 };
719 
720 
721 /** Class for domain integration of L(v) := (f, v), where
722  f=(f1,...,fn) and v=(v1,...,vn). that makes use of
723  VectorQuadratureFunctionCoefficient*/
725 {
726 private:
728 
729 public:
731  const IntegrationRule *ir)
732  : LinearFormIntegrator(ir), vqfc(vqfc)
733  {
734  if (ir)
735  {
736  MFEM_WARNING("Integration rule not used in this class. "
737  "The QuadratureFunction integration rules are used instead");
738  }
739  }
740 
742  virtual void AssembleRHSElementVect(const FiniteElement &fe,
744  Vector &elvect);
745 
746  virtual void SetIntRule(const IntegrationRule *ir)
747  {
748  MFEM_WARNING("Integration rule not used in this class. "
749  "The QuadratureFunction integration rules are used instead");
750  }
751 };
752 
753 
754 /** Class for domain integration L(v) := (f, v) that makes use
755  of QuadratureFunctionCoefficient. */
757 {
758 private:
760 
761 public:
763  const IntegrationRule *ir)
764  : LinearFormIntegrator(ir), qfc(qfc)
765  {
766  if (ir)
767  {
768  MFEM_WARNING("Integration rule not used in this class. "
769  "The QuadratureFunction integration rules are used instead");
770  }
771  }
772 
774  virtual void AssembleRHSElementVect(const FiniteElement &fe,
776  Vector &elvect);
777 
778  virtual void SetIntRule(const IntegrationRule *ir)
779  {
780  MFEM_WARNING("Integration rule not used in this class. "
781  "The QuadratureFunction integration rules are used instead");
782  }
783 };
784 
785 }
786 
787 
788 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:233
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:110
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:622
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
void GetDeltaCenter(Vector &center)
VectorFEDomainLFCurlIntegrator(VectorCoefficient &F)
Constructs the domain integrator (Q, curl v)
Definition: lininteg.hpp:384
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:220
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
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:580
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:1006
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:552
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:710
VectorBoundaryLFIntegrator(VectorCoefficient &QG)
Constructs a boundary integrator with a given VectorCoefficient QG.
Definition: lininteg.hpp:329
Class for domain integrator L(v) := (f, grad v)
Definition: lininteg.hpp:145
BoundaryNormalLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:217
DGElasticityDirichletLFIntegrator(VectorCoefficient &uD_, Coefficient &lambda_, Coefficient &mu_, double alpha_, double kappa_)
Definition: lininteg.hpp:597
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:867
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:78
WhiteGaussianNoiseDomainLFIntegrator(int seed_=0)
Sets the seed_ of the random number generator. A fixed seed allows for a reproducible sequence of whi...
Definition: lininteg.hpp:645
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:69
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:77
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
virtual ~LinearFormIntegrator()
Definition: lininteg.hpp:57
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:55
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:268
Class for boundary integration L(v) := (g, v)
Definition: lininteg.hpp:179
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:54
Vector beta
const IntegrationRule * IntRule
Definition: lininteg.hpp:27
BoundaryFlowIntegrator(Coefficient &f_, VectorCoefficient &u_, double a, double b)
Definition: lininteg.hpp:511
DGDirichletLFIntegrator(Coefficient &u, const double s, const double k)
Definition: lininteg.hpp:548
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:126
Abstract base class LinearFormIntegrator.
Definition: lininteg.hpp:24
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) override
Method defining assembly on device.
void DeleteAll()
Delete the whole array.
Definition: array.hpp:854
WhiteGaussianNoiseDomainLFIntegrator(MPI_Comm comm_, int seed_)
Sets the MPI communicator comm_ and the seed_ of the random number generator. A fixed seed allows for...
Definition: lininteg.hpp:654
void GetDeltaCenter(Vector &center)
Returns the center of the delta coefficient.
Definition: lininteg.hpp:88
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:190
LinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:29
VectorBoundaryFluxLFIntegrator(Coefficient &f, double s=1.0, const IntegrationRule *ir=NULL)
Definition: lininteg.hpp:433
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
virtual bool SupportsDevice() const override
Method probing for assembly on device.
Definition: lininteg.hpp:299
DomainLFGradIntegrator(VectorCoefficient &QF)
Constructs the domain integrator (Q, grad v)
Definition: lininteg.hpp:154
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
virtual void AssembleDeltaElementVect(const FiniteElement &fe, ElementTransformation &Trans, Vector &elvect) override
Assemble the delta coefficient at the IntegrationPoint set in Trans which is assumed to map to the de...
Definition: lininteg.cpp:369
BoundaryFlowIntegrator(Coefficient &f_, VectorCoefficient &u_, double a)
Definition: lininteg.hpp:507
Class for boundary integration .
Definition: lininteg.hpp:210
DomainLFIntegrator(Coefficient &QF, int a=2, int b=0)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:115
double b
Definition: lissajous.cpp:42
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:375
VectorFEBoundaryFluxLFIntegrator(int a=1, int b=-1)
Definition: lininteg.hpp:455
VectorFEBoundaryTangentLFIntegrator(VectorCoefficient &QG, int a=2, int b=0)
Definition: lininteg.hpp:481
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:187
VectorDeltaCoefficient * vec_delta
Definition: lininteg.hpp: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:535
MatrixCoefficient * MQ
Definition: lininteg.hpp:540
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
Definition: lininteg.cpp:18
VectorFEBoundaryFluxLFIntegrator(Coefficient &f, int a=2, int b=0)
Definition: lininteg.hpp:457
void SetSeed(int seed)
Sets/resets the seed of the random number generator.
Definition: lininteg.hpp:673
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:746
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:228
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:776
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:466
DomainLFIntegrator(Coefficient &QF, const IntegrationRule *ir)
Constructs a domain integrator with a given Coefficient.
Definition: lininteg.hpp:121
void GetDeltaCenter(Vector &center)
Write the center of the delta function into center.
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:265
VectorFEDomainLFDivIntegrator(Coefficient &QF)
Constructs the domain integrator (Q, div v)
Definition: lininteg.hpp:406
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)=0
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:157
DeltaLFIntegrator(Coefficient &q, const IntegrationRule *ir=NULL)
This constructor should be used by derived classes that use a scalar DeltaCoefficient.
Definition: lininteg.hpp:70
QuadratureLFIntegrator(QuadratureFunctionCoefficient &qfc, const IntegrationRule *ir)
Definition: lininteg.hpp:762
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
DeltaCoefficient * delta
Definition: lininteg.hpp:65
Base class for Matrix Coefficients that optionally depend on time and space.
Abstract class for integrators that support delta coefficients.
Definition: lininteg.hpp:62
VectorDomainLFIntegrator(VectorCoefficient &QF)
Constructs a domain integrator with a given VectorCoefficient.
Definition: lininteg.hpp:262
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:589
VectorFEDomainLFIntegrator(VectorCoefficient &F)
Definition: lininteg.hpp:354
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:654
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:124
double a
Definition: lissajous.cpp:41
for Nedelec Elements)
Definition: lininteg.hpp:375
void ResetFactors(int NE)
Resets the array of saved lower triangular Cholesky decomposition matrices. The parameter NE should b...
Definition: lininteg.hpp:695
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:487
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
BoundaryTangentialLFIntegrator(VectorCoefficient &QG, int a=1, int b=1)
Constructs a boundary integrator with a given Coefficient QG.
Definition: lininteg.hpp:242
virtual void SetIntRule(const IntegrationRule *ir)
Definition: lininteg.hpp:778
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:34
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:1044
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:189
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:503
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:346
const double alpha
Definition: ex15.cpp:369
virtual void AssembleRHSElementVect(const FiniteElement &fe, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:1075
Vector data type.
Definition: vector.hpp:58
void SaveFactors(int NE)
Saves the lower triangular matrices in the element-wise Cholesky decomposition. The parameter NE shou...
Definition: lininteg.hpp:686
bool IsDelta() const
Returns true if the derived class instance uses a delta coefficient.
Definition: lininteg.hpp:85
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:453
virtual bool SupportsDevice() const
Method probing for assembly on device.
Definition: lininteg.hpp:365
Class for boundary integration .
Definition: lininteg.hpp:474
VectorDomainLFGradIntegrator(VectorCoefficient &QF)
Constructs the domain integrator (Q, grad v)
Definition: lininteg.hpp:296
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:310
virtual void AssembleDevice(const FiniteElementSpace &fes, const Array< int > &markers, Vector &b)
Method defining assembly on device.
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
DGDirichletLFIntegrator(Coefficient &u, MatrixCoefficient &q, const double s, const double k)
Definition: lininteg.hpp:553
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect)
Definition: lininteg.cpp:38
Class for boundary integration in 2D.
Definition: lininteg.hpp:235
VectorQuadratureLFIntegrator(VectorQuadratureFunctionCoefficient &vqfc, const IntegrationRule *ir)
Definition: lininteg.hpp:730
virtual void AssembleRHSElementVect(const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
Definition: lininteg.cpp:327
DGDirichletLFIntegrator(Coefficient &u, Coefficient &q, const double s, const double k)
Definition: lininteg.hpp:550
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...
double f(const Vector &p)