MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilinearform.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_BILINEARFORM
13 #define MFEM_BILINEARFORM
14 
15 #include "../config/config.hpp"
16 #include "../linalg/linalg.hpp"
17 #include "fespace.hpp"
18 #include "gridfunc.hpp"
19 #include "linearform.hpp"
20 #include "bilininteg.hpp"
21 #include "bilinearform_ext.hpp"
22 #include "staticcond.hpp"
23 #include "hybridization.hpp"
24 
25 namespace mfem
26 {
27 
28 /// Enumeration defining the assembly level for bilinear and nonlinear form
29 /// classes derived from Operator.
30 enum class AssemblyLevel
31 {
32  /// Fully assembled form, i.e. a global sparse matrix in MFEM, Hypre or PETSC
33  /// format.
34  FULL,
35  /// Form assembled at element level, which computes and stores dense element
36  /// matrices.
37  ELEMENT,
38  /// Partially-assembled form, which computes and stores data only at
39  /// quadrature points.
40  PARTIAL,
41  /// "Matrix-free" form that computes all of its action on-the-fly without any
42  /// substantial storage.
43  NONE,
44 };
45 
46 
47 /** Class for bilinear form - "Matrix" with associated FE space and
48  BLFIntegrators. */
49 class BilinearForm : public Matrix
50 {
51 protected:
52  /// Sparse matrix to be associated with the form. Owned.
54 
55  /// Matrix used to eliminate b.c. Owned.
57 
58  /// FE space on which the form lives. Not owned.
60 
61  /// The form assembly level (full, partial, etc.)
63  /// Element batch size used in the form action (1, 8, num_elems, etc.)
64  int batch;
65  /** Extension for supporting Full Assembly (FA), Element Assembly (EA),
66  Partial Assembly (PA), or Matrix Free assembly (MF). */
68 
69  /// Indicates the Mesh::sequence corresponding to the current state of the
70  /// BilinearForm.
71  long sequence;
72 
73  /** @brief Indicates the BilinearFormIntegrator%s stored in #dbfi, #bbfi,
74  #fbfi, and #bfbfi are owned by another BilinearForm. */
76 
77  /// Set of Domain Integrators to be applied.
79 
80  /// Set of Boundary Integrators to be applied.
82  Array<Array<int>*> bbfi_marker; ///< Entries are not owned.
83 
84  /// Set of interior face Integrators to be applied.
86 
87  /// Set of boundary face Integrators to be applied.
89  Array<Array<int>*> bfbfi_marker; ///< Entries are not owned.
90 
93 
95 
98 
99  /**
100  * This member allows one to specify what should be done
101  * to the diagonal matrix entries and corresponding RHS
102  * values upon elimination of the constrained DoFs.
103  */
105 
107  // Allocate appropriate SparseMatrix and assign it to mat
108  void AllocMat();
109 
110  void ConformingAssemble();
111 
112  // may be used in the construction of derived classes
114  {
115  fes = NULL; sequence = -1;
116  mat = mat_e = NULL; extern_bfs = 0; element_matrices = NULL;
117  static_cond = NULL; hybridization = NULL;
121  batch = 1;
122  ext = NULL;
123  }
124 
125 private:
126  /// Copy construction is not supported; body is undefined.
127  BilinearForm(const BilinearForm &);
128 
129  /// Copy assignment is not supported; body is undefined.
130  BilinearForm &operator=(const BilinearForm &);
131 
132 public:
133  /// Creates bilinear form associated with FE space @a *f.
134  /** The pointer @a f is not owned by the newly constructed object. */
136 
137  /** @brief Create a BilinearForm on the FiniteElementSpace @a f, using the
138  same integrators as the BilinearForm @a bf.
139 
140  The pointer @a f is not owned by the newly constructed object.
141 
142  The integrators in @a bf are copied as pointers and they are not owned by
143  the newly constructed BilinearForm.
144 
145  The optional parameter @a ps is used to initialize the internal flag
146  #precompute_sparsity, see UsePrecomputedSparsity() for details. */
147  BilinearForm(FiniteElementSpace *f, BilinearForm *bf, int ps = 0);
148 
149  /// Get the size of the BilinearForm as a square matrix.
150  int Size() const { return height; }
151 
152  /// Set the desired assembly level. The default is AssemblyLevel::FULL.
153  /** This method must be called before assembly. */
154  void SetAssemblyLevel(AssemblyLevel assembly_level);
155 
156  /** Enable the use of static condensation. For details see the description
157  for class StaticCondensation in fem/staticcond.hpp This method should be
158  called before assembly. If the number of unknowns after static
159  condensation is not reduced, it is not enabled. */
161 
162  /** Check if static condensation was actually enabled by a previous call to
163  EnableStaticCondensation(). */
164  bool StaticCondensationIsEnabled() const { return static_cond; }
165 
166  /// Return the trace FE space associated with static condensation.
168  { return static_cond ? static_cond->GetTraceFESpace() : NULL; }
169 
170  /** Enable hybridization; for details see the description for class
171  Hybridization in fem/hybridization.hpp. This method should be called
172  before assembly. */
173  void EnableHybridization(FiniteElementSpace *constr_space,
174  BilinearFormIntegrator *constr_integ,
175  const Array<int> &ess_tdof_list);
176 
177  /** For scalar FE spaces, precompute the sparsity pattern of the matrix
178  (assuming dense element matrices) based on the types of integrators
179  present in the bilinear form. */
180  void UsePrecomputedSparsity(int ps = 1) { precompute_sparsity = ps; }
181 
182  /** @brief Use the given CSR sparsity pattern to allocate the internal
183  SparseMatrix.
184 
185  - The @a I and @a J arrays must define a square graph with size equal to
186  GetVSize() of the associated FiniteElementSpace.
187  - This method should be called after enabling static condensation or
188  hybridization, if used.
189  - In the case of static condensation, @a I and @a J are not used.
190  - The ownership of the arrays @a I and @a J remains with the caller. */
191  void UseSparsity(int *I, int *J, bool isSorted);
192 
193  /// Use the sparsity of @a A to allocate the internal SparseMatrix.
194  void UseSparsity(SparseMatrix &A);
195 
196  /** Pre-allocate the internal SparseMatrix before assembly. If the flag
197  'precompute sparsity' is set, the matrix is allocated in CSR format (i.e.
198  finalized) and the entries are initialized with zeros. */
199  void AllocateMatrix() { if (mat == NULL) { AllocMat(); } }
200 
201  /// Access all integrators added with AddDomainIntegrator().
203 
204  /// Access all integrators added with AddBoundaryIntegrator().
206  /** @brief Access all boundary markers added with AddBoundaryIntegrator().
207  If no marker was specified when the integrator was added, the
208  corresponding pointer (to Array<int>) will be NULL. */
210 
211  /// Access all integrators added with AddInteriorFaceIntegrator().
213 
214  /// Access all integrators added with AddBdrFaceIntegrator().
216  /** @brief Access all boundary markers added with AddBdrFaceIntegrator().
217  If no marker was specified when the integrator was added, the
218  corresponding pointer (to Array<int>) will be NULL. */
220 
221  const double &operator()(int i, int j) { return (*mat)(i,j); }
222 
223  /// Returns reference to a_{ij}.
224  virtual double &Elem(int i, int j);
225 
226  /// Returns constant reference to a_{ij}.
227  virtual const double &Elem(int i, int j) const;
228 
229  /// Matrix vector multiplication.
230  virtual void Mult(const Vector &x, Vector &y) const { mat->Mult(x, y); }
231 
232  void FullMult(const Vector &x, Vector &y) const
233  { mat->Mult(x, y); mat_e->AddMult(x, y); }
234 
235  virtual void AddMult(const Vector &x, Vector &y, const double a = 1.0) const
236  { mat -> AddMult (x, y, a); }
237 
238  void FullAddMult(const Vector &x, Vector &y) const
239  { mat->AddMult(x, y); mat_e->AddMult(x, y); }
240 
241  virtual void AddMultTranspose(const Vector & x, Vector & y,
242  const double a = 1.0) const
243  { mat->AddMultTranspose(x, y, a); }
244 
245  void FullAddMultTranspose(const Vector & x, Vector & y) const
246  { mat->AddMultTranspose(x, y); mat_e->AddMultTranspose(x, y); }
247 
248  virtual void MultTranspose(const Vector & x, Vector & y) const
249  { y = 0.0; AddMultTranspose (x, y); }
250 
251  double InnerProduct(const Vector &x, const Vector &y) const
252  { return mat->InnerProduct (x, y); }
253 
254  /// Returns a pointer to (approximation) of the matrix inverse.
255  virtual MatrixInverse *Inverse() const;
256 
257  /// Finalizes the matrix initialization.
258  virtual void Finalize(int skip_zeros = 1);
259 
260  /// Returns a reference to the sparse matrix
261  const SparseMatrix &SpMat() const
262  {
263  MFEM_VERIFY(mat, "mat is NULL and can't be dereferenced");
264  return *mat;
265  }
267  {
268  MFEM_VERIFY(mat, "mat is NULL and can't be dereferenced");
269  return *mat;
270  }
271  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
272 
273  /// Returns a reference to the sparse matrix of eliminated b.c.
274  const SparseMatrix &SpMatElim() const
275  {
276  MFEM_VERIFY(mat_e, "mat_e is NULL and can't be dereferenced");
277  return *mat_e;
278  }
280  {
281  MFEM_VERIFY(mat_e, "mat_e is NULL and can't be dereferenced");
282  return *mat_e;
283  }
284 
285  /// Adds new Domain Integrator. Assumes ownership of @a bfi.
287 
288  /// Adds new Boundary Integrator. Assumes ownership of @a bfi.
290 
291  /** @brief Adds new Boundary Integrator, restricted to specific boundary
292  attributes.
293 
294  Assumes ownership of @a bfi. The array @a bdr_marker is stored internally
295  as a pointer to the given Array<int> object. */
297  Array<int> &bdr_marker);
298 
299  /// Adds new interior Face Integrator. Assumes ownership of @a bfi.
301 
302  /// Adds new boundary Face Integrator. Assumes ownership of @a bfi.
304 
305  /** @brief Adds new boundary Face Integrator, restricted to specific boundary
306  attributes.
307 
308  Assumes ownership of @a bfi. The array @a bdr_marker is stored internally
309  as a pointer to the given Array<int> object. */
311  Array<int> &bdr_marker);
312 
313  void operator=(const double a)
314  {
315  if (mat != NULL) { *mat = a; }
316  if (mat_e != NULL) { *mat_e = a; }
317  }
318 
319  /// Assembles the form i.e. sums over all domain/bdr integrators.
320  void Assemble(int skip_zeros = 1);
321 
322  /// Get the finite element space prolongation matrix
323  virtual const Operator *GetProlongation() const
324  { return fes->GetConformingProlongation(); }
325  /// Get the finite element space restriction matrix
326  virtual const Operator *GetRestriction() const
327  { return fes->GetConformingRestriction(); }
328 
329  /** @brief Form the linear system A X = B, corresponding to this bilinear
330  form and the linear form @a b(.). */
331  /** This method applies any necessary transformations to the linear system
332  such as: eliminating boundary conditions; applying conforming constraints
333  for non-conforming AMR; parallel assembly; static condensation;
334  hybridization.
335 
336  The GridFunction-size vector @a x must contain the essential b.c. The
337  BilinearForm and the LinearForm-size vector @a b must be assembled.
338 
339  The vector @a X is initialized with a suitable initial guess: when using
340  hybridization, the vector @a X is set to zero; otherwise, the essential
341  entries of @a X are set to the corresponding b.c. and all other entries
342  are set to zero (@a copy_interior == 0) or copied from @a x
343  (@a copy_interior != 0).
344 
345  This method can be called multiple times (with the same @a ess_tdof_list
346  array) to initialize different right-hand sides and boundary condition
347  values.
348 
349  After solving the linear system, the finite element solution @a x can be
350  recovered by calling RecoverFEMSolution() (with the same vectors @a X,
351  @a b, and @a x).
352 
353  NOTE: If there are no transformations, @a X simply reuses the data of
354  @a x. */
355  virtual void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x,
356  Vector &b, OperatorHandle &A, Vector &X,
357  Vector &B, int copy_interior = 0);
358 
359  /** @brief Form the linear system A X = B, corresponding to this bilinear
360  form and the linear form @a b(.). */
361  /** Version of the method FormLinearSystem() where the system matrix is
362  returned in the variable @a A, of type OpType, holding a *reference* to
363  the system matrix (created with the method OpType::MakeRef()). The
364  reference will be invalidated when SetOperatorType(), Update(), or the
365  destructor is called.
366 
367  Currently, this method can be used only with AssemblyLevel::FULL. */
368  template <typename OpType>
369  void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x, Vector &b,
370  OpType &A, Vector &X, Vector &B,
371  int copy_interior = 0)
372  {
373  OperatorHandle Ah;
374  FormLinearSystem(ess_tdof_list, x, b, Ah, X, B, copy_interior);
375  OpType *A_ptr = Ah.Is<OpType>();
376  MFEM_VERIFY(A_ptr, "invalid OpType used");
377  A.MakeRef(*A_ptr);
378  }
379 
380  /// Form the linear system matrix @a A, see FormLinearSystem() for details.
381  virtual void FormSystemMatrix(const Array<int> &ess_tdof_list,
382  OperatorHandle &A);
383 
384  /// Form the linear system matrix A, see FormLinearSystem() for details.
385  /** Version of the method FormSystemMatrix() where the system matrix is
386  returned in the variable @a A, of type OpType, holding a *reference* to
387  the system matrix (created with the method OpType::MakeRef()). The
388  reference will be invalidated when SetOperatorType(), Update(), or the
389  destructor is called.
390 
391  Currently, this method can be used only with AssemblyLevel::FULL. */
392  template <typename OpType>
393  void FormSystemMatrix(const Array<int> &ess_tdof_list, OpType &A)
394  {
395  OperatorHandle Ah;
396  FormSystemMatrix(ess_tdof_list, Ah);
397  OpType *A_ptr = Ah.Is<OpType>();
398  MFEM_VERIFY(A_ptr, "invalid OpType used");
399  A.MakeRef(*A_ptr);
400  }
401 
402  /// Recover the solution of a linear system formed with FormLinearSystem().
403  /** Call this method after solving a linear system constructed using the
404  FormLinearSystem() method to recover the solution as a GridFunction-size
405  vector in @a x. Use the same arguments as in the FormLinearSystem() call.
406  */
407  virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
408 
409  /// Compute and store internally all element matrices.
410  void ComputeElementMatrices();
411 
412  /// Free the memory used by the element matrices.
414  { delete element_matrices; element_matrices = NULL; }
415 
416  void ComputeElementMatrix(int i, DenseMatrix &elmat);
417  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
418  Array<int> &vdofs, int skip_zeros = 1);
419  void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat,
420  Array<int> &vdofs, int skip_zeros = 1);
421 
422  /// Eliminate essential boundary DOFs from the system.
423  /** The array @a bdr_attr_is_ess marks boundary attributes that constitute
424  the essential part of the boundary. By default, the diagonal at the
425  essential DOFs is set to 1.0. This behavior is controlled by the argument
426  @a dpolicy. */
427  void EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
428  const Vector &sol, Vector &rhs,
429  DiagonalPolicy dpolicy = DIAG_ONE);
430 
431  /// Eliminate essential boundary DOFs from the system matrix.
432  void EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
433  DiagonalPolicy dpolicy = DIAG_ONE);
434  /// Perform elimination and set the diagonal entry to the given value
435  void EliminateEssentialBCDiag(const Array<int> &bdr_attr_is_ess,
436  double value);
437 
438  /// Eliminate the given @a vdofs. NOTE: here, @a vdofs is a list of DOFs.
439  void EliminateVDofs(const Array<int> &vdofs, const Vector &sol, Vector &rhs,
440  DiagonalPolicy dpolicy = DIAG_ONE);
441 
442  /// Eliminate the given @a vdofs, storing the eliminated part internally.
443  /** This method works in conjunction with EliminateVDofsInRHS() and allows
444  elimination of boundary conditions in multiple right-hand sides. In this
445  method, @a vdofs is a list of DOFs. */
446  void EliminateVDofs(const Array<int> &vdofs,
447  DiagonalPolicy dpolicy = DIAG_ONE);
448 
449  /** @brief Similar to
450  EliminateVDofs(const Array<int> &, const Vector &, Vector &, DiagonalPolicy)
451  but here @a ess_dofs is a marker (boolean) array on all vector-dofs
452  (@a ess_dofs[i] < 0 is true). */
453  void EliminateEssentialBCFromDofs(const Array<int> &ess_dofs, const Vector &sol,
454  Vector &rhs, DiagonalPolicy dpolicy = DIAG_ONE);
455 
456  /** @brief Similar to EliminateVDofs(const Array<int> &, DiagonalPolicy) but
457  here @a ess_dofs is a marker (boolean) array on all vector-dofs
458  (@a ess_dofs[i] < 0 is true). */
459  void EliminateEssentialBCFromDofs(const Array<int> &ess_dofs,
460  DiagonalPolicy dpolicy = DIAG_ONE);
461  /// Perform elimination and set the diagonal entry to the given value
462  void EliminateEssentialBCFromDofsDiag(const Array<int> &ess_dofs,
463  double value);
464 
465  /** @brief Use the stored eliminated part of the matrix (see
466  EliminateVDofs(const Array<int> &, DiagonalPolicy)) to modify the r.h.s.
467  @a b; @a vdofs is a list of DOFs (non-directional, i.e. >= 0). */
468  void EliminateVDofsInRHS(const Array<int> &vdofs, const Vector &x,
469  Vector &b);
470 
471  double FullInnerProduct(const Vector &x, const Vector &y) const
472  { return mat->InnerProduct(x, y) + mat_e->InnerProduct(x, y); }
473 
474  virtual void Update(FiniteElementSpace *nfes = NULL);
475 
476  /// (DEPRECATED) Return the FE space associated with the BilinearForm.
477  /** @deprecated Use FESpace() instead. */
479 
480  /// Return the FE space associated with the BilinearForm.
482  /// Read-only access to the associated FiniteElementSpace.
483  const FiniteElementSpace *FESpace() const { return fes; }
484 
485  /// Sets diagonal policy used upon construction of the linear system
486  void SetDiagonalPolicy(DiagonalPolicy policy);
487 
488  /// Destroys bilinear form.
489  virtual ~BilinearForm();
490 };
491 
492 
493 /**
494  Class for assembling of bilinear forms `a(u,v)` defined on different
495  trial and test spaces. The assembled matrix `A` is such that
496 
497  a(u,v) = V^t A U
498 
499  where `U` and `V` are the vectors representing the functions `u` and `v`,
500  respectively. The first argument, `u`, of `a(,)` is in the trial space
501  and the second argument, `v`, is in the test space. Thus,
502 
503  # of rows of A = dimension of the test space and
504  # of cols of A = dimension of the trial space.
505 
506  Both trial and test spaces should be defined on the same mesh.
507 */
508 class MixedBilinearForm : public Matrix
509 {
510 protected:
511  SparseMatrix *mat; ///< Owned.
512 
513  FiniteElementSpace *trial_fes, ///< Not owned
514  *test_fes; ///< Not owned
515 
516  /** @brief Indicates the BilinearFormIntegrator%s stored in #dom, #bdr, and
517  #skt are owned by another MixedBilinearForm. */
519 
520  /// Domain integrators.
522  /// Boundary integrators.
524  /// Trace face (skeleton) integrators.
526 
527 private:
528  /// Copy construction is not supported; body is undefined.
530 
531  /// Copy assignment is not supported; body is undefined.
532  MixedBilinearForm &operator=(const MixedBilinearForm &);
533 
534 public:
535  /** @brief Construct a MixedBilinearForm on the given trial, @a tr_fes, and
536  test, @a te_fes, FiniteElementSpace%s. */
537  /** The pointers @a tr_fes and @a te_fes are not owned by the newly
538  constructed object. */
540  FiniteElementSpace *te_fes);
541 
542  /** @brief Create a MixedBilinearForm on the given trial, @a tr_fes, and
543  test, @a te_fes, FiniteElementSpace%s, using the same integrators as the
544  MixedBilinearForm @a mbf.
545 
546  The pointers @a tr_fes and @a te_fes are not owned by the newly
547  constructed object.
548 
549  The integrators in @a mbf are copied as pointers and they are not owned
550  by the newly constructed MixedBilinearForm. */
552  FiniteElementSpace *te_fes,
553  MixedBilinearForm *mbf);
554 
555  virtual double &Elem(int i, int j);
556 
557  virtual const double &Elem(int i, int j) const;
558 
559  virtual void Mult(const Vector & x, Vector & y) const;
560 
561  virtual void AddMult(const Vector & x, Vector & y,
562  const double a = 1.0) const;
563 
564  virtual void AddMultTranspose(const Vector & x, Vector & y,
565  const double a = 1.0) const;
566 
567  virtual void MultTranspose(const Vector & x, Vector & y) const
568  { y = 0.0; AddMultTranspose (x, y); }
569 
570  virtual MatrixInverse *Inverse() const;
571 
572  virtual void Finalize(int skip_zeros = 1);
573 
574  /** Extract the associated matrix as SparseMatrix blocks. The number of
575  block rows and columns is given by the vector dimensions (vdim) of the
576  test and trial spaces, respectively. */
577  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
578 
579  const SparseMatrix &SpMat() const { return *mat; }
580  SparseMatrix &SpMat() { return *mat; }
581  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
582 
583  /// Adds a domain integrator. Assumes ownership of @a bfi.
585 
586  /// Adds a boundary integrator. Assumes ownership of @a bfi.
588 
589  /** @brief Add a trace face integrator. Assumes ownership of @a bfi.
590 
591  This type of integrator assembles terms over all faces of the mesh using
592  the face FE from the trial space and the two adjacent volume FEs from the
593  test space. */
595 
596  /// Access all integrators added with AddDomainIntegrator().
598 
599  /// Access all integrators added with AddBoundaryIntegrator().
601 
602  /// Access all integrators added with AddTraceFaceIntegrator().
604 
605  void operator=(const double a) { *mat = a; }
606 
607  void Assemble(int skip_zeros = 1);
608 
609  /** For partially conforming trial and/or test FE spaces, complete the
610  assembly process by performing A := P2^t A P1 where A is the internal
611  sparse matrix; P1 and P2 are the conforming prolongation matrices of the
612  trial and test FE spaces, respectively. After this call the
613  MixedBilinearForm becomes an operator on the conforming FE spaces. */
614  void ConformingAssemble();
615 
616  void EliminateTrialDofs(Array<int> &bdr_attr_is_ess,
617  const Vector &sol, Vector &rhs);
618 
620  const Vector &sol, Vector &rhs);
621 
622  virtual void EliminateTestDofs(Array<int> &bdr_attr_is_ess);
623 
624  void Update();
625 
626  virtual ~MixedBilinearForm();
627 };
628 
629 
630 /**
631  Class for constructing the matrix representation of a linear operator,
632  `v = L u`, from one FiniteElementSpace (domain) to another FiniteElementSpace
633  (range). The constructed matrix `A` is such that
634 
635  V = A U
636 
637  where `U` and `V` are the vectors of degrees of freedom representing the
638  functions `u` and `v`, respectively. The dimensions of `A` are
639 
640  number of rows of A = dimension of the range space and
641  number of cols of A = dimension of the domain space.
642 
643  This class is very similar to MixedBilinearForm. One difference is that
644  the linear operator `L` is defined using a special kind of
645  BilinearFormIntegrator (we reuse its functionality instead of defining a
646  new class). The other difference with the MixedBilinearForm class is that
647  the "assembly" process overwrites the global matrix entries using the
648  local element matrices instead of adding them.
649 
650  Note that if we define the bilinear form `b(u,v) := (Lu,v)` using an inner
651  product in the range space, then its matrix representation, `B`, is
652 
653  B = M A, (since V^t B U = b(u,v) = (Lu,v) = V^t M A U)
654 
655  where `M` denotes the mass matrix for the inner product in the range space:
656  `V1^t M V2 = (v1,v2)`. Similarly, if `c(u,w) := (Lu,Lw)` then
657 
658  C = A^t M A.
659 */
661 {
662 private:
663  /// Copy construction is not supported; body is undefined.
665 
666  /// Copy assignment is not supported; body is undefined.
668 
669 public:
670  /** @brief Construct a DiscreteLinearOperator on the given
671  FiniteElementSpace%s @a domain_fes and @a range_fes. */
672  /** The pointers @a domain_fes and @a range_fes are not owned by the newly
673  constructed object. */
675  FiniteElementSpace *range_fes)
676  : MixedBilinearForm(domain_fes, range_fes) { }
677 
678  /// Adds a domain interpolator. Assumes ownership of @a di.
680  { AddDomainIntegrator(di); }
681 
682  /// Adds a trace face interpolator. Assumes ownership of @a di.
684  { AddTraceFaceIntegrator(di); }
685 
686  /// Access all interpolators added with AddDomainInterpolator().
688 
689  /** @brief Construct the internal matrix representation of the discrete
690  linear operator. */
691  virtual void Assemble(int skip_zeros = 1);
692 };
693 
694 }
695 
696 #endif
void FormSystemMatrix(const Array< int > &ess_tdof_list, OpType &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void EliminateEssentialBCFromDofs(const Array< int > &ess_dofs, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Similar to EliminateVDofs(const Array&lt;int&gt; &amp;, const Vector &amp;, Vector &amp;, DiagonalPolicy) but here ess_...
void FullAddMult(const Vector &x, Vector &y) const
SparseMatrix & SpMat()
Array< BilinearFormIntegrator * > * GetBBFI()
Access all integrators added with AddBoundaryIntegrator().
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Array< BilinearFormIntegrator * > skt
Trace face (skeleton) integrators.
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:769
virtual MatrixInverse * Inverse() const
Returns a pointer to (approximation) of the matrix inverse.
Array< BilinearFormIntegrator * > fbfi
Set of interior face Integrators to be applied.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
BilinearFormExtension * ext
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
Array< BilinearFormIntegrator * > dbfi
Set of Domain Integrators to be applied.
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
void EliminateEssentialBCFromTrialDofs(Array< int > &marked_vdofs, const Vector &sol, Vector &rhs)
Array< Array< int > * > bbfi_marker
Entries are not owned.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
void EliminateTrialDofs(Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs)
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
Definition: sparsemat.cpp:552
int batch
Element batch size used in the form action (1, 8, num_elems, etc.)
Array< Array< int > * > * GetBFBFI_Marker()
Access all boundary markers added with AddBdrFaceIntegrator(). If no marker was specified when the in...
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Definition: sparsemat.cpp:623
Abstract data type for matrix inverse.
Definition: matrix.hpp:66
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::FULL.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation matrix.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OpType &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
double FullInnerProduct(const Vector &x, const Vector &y) const
Array< Array< int > * > * GetBBFI_Marker()
Access all boundary markers added with AddBoundaryIntegrator(). If no marker was specified when the i...
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
void SetDiagonalPolicy(DiagonalPolicy policy)
Sets diagonal policy used upon construction of the linear system.
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
void EliminateEssentialBCFromDofsDiag(const Array< int > &ess_dofs, double value)
Perform elimination and set the diagonal entry to the given value.
Array< BilinearFormIntegrator * > bdr
Boundary integrators.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
SparseMatrix & SpMat()
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:776
void FullAddMultTranspose(const Vector &x, Vector &y) const
OpType * Is() const
Return the Operator pointer dynamically cast to a specified OpType.
Definition: handle.hpp:100
Keep the diagonal value.
Definition: matrix.hpp:36
DiscreteLinearOperator(FiniteElementSpace *domain_fes, FiniteElementSpace *range_fes)
Construct a DiscreteLinearOperator on the given FiniteElementSpaces domain_fes and range_fes...
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
Data type sparse matrix.
Definition: sparsemat.hpp:40
StaticCondensation * static_cond
Owned.
const FiniteElementSpace * FESpace() const
Read-only access to the associated FiniteElementSpace.
SparseMatrix * mat
Sparse matrix to be associated with the form. Owned.
void EliminateVDofs(const Array< int > &vdofs, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs.
int Size() const
Get the size of the BilinearForm as a square matrix.
Abstract data type matrix.
Definition: matrix.hpp:27
int extern_bfs
Indicates the BilinearFormIntegrators stored in dbfi, bbfi, fbfi, and bfbfi are owned by another Bili...
void EliminateEssentialBC(const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate essential boundary DOFs from the system.
void Assemble(int skip_zeros=1)
void FullMult(const Vector &x, Vector &y) const
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void EliminateEssentialBCDiag(const Array< int > &bdr_attr_is_ess, double value)
Perform elimination and set the diagonal entry to the given value.
Class extending the BilinearForm class to support the different AssemblyLevels.
virtual void Update(FiniteElementSpace *nfes=NULL)
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
SparseMatrix * mat
Owned.
FiniteElementSpace * GetFES()
(DEPRECATED) Return the FE space associated with the BilinearForm.
void AddTraceFaceIntegrator(BilinearFormIntegrator *bfi)
Add a trace face integrator. Assumes ownership of bfi.
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
Array< BilinearFormIntegrator * > * GetTFBFI()
Access all integrators added with AddTraceFaceIntegrator().
Set the diagonal value to one.
Definition: matrix.hpp:35
Dynamic 2D array using row-major layout.
Definition: array.hpp:304
Array< BilinearFormIntegrator * > * GetBBFI()
Access all integrators added with AddBoundaryIntegrator().
void ComputeElementMatrix(int i, DenseMatrix &elmat)
SparseMatrix & SpMatElim()
SparseMatrix * mat_e
Matrix used to eliminate b.c. Owned.
bool StaticCondensationIsEnabled() const
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:23
virtual const Operator * GetRestriction() const
Get the finite element space restriction matrix.
void AssembleElementMatrix(int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
void UsePrecomputedSparsity(int ps=1)
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
SparseMatrix * LoseMat()
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:545
Array< BilinearFormIntegrator * > bfbfi
Set of boundary face Integrators to be applied.
Array< BilinearFormIntegrator * > bbfi
Set of Boundary Integrators to be applied.
DenseTensor * element_matrices
Owned.
SparseMatrix * LoseMat()
void operator=(const double a)
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
void AddTraceFaceInterpolator(DiscreteInterpolator *di)
Adds a trace face interpolator. Assumes ownership of di.
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Array< BilinearFormIntegrator * > dom
Domain integrators.
double InnerProduct(const Vector &x, const Vector &y) const
DiagonalPolicy diag_policy
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:24
const SparseMatrix & SpMat() const
FiniteElementSpace * test_fes
Not owned.
const double & operator()(int i, int j)
FiniteElementSpace * GetTraceFESpace()
Return a pointer to the reduced/trace FE space.
Definition: staticcond.hpp:109
virtual void EliminateTestDofs(Array< int > &bdr_attr_is_ess)
int extern_bfs
Indicates the BilinearFormIntegrators stored in dom, bdr, and skt are owned by another MixedBilinearF...
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void EliminateVDofsInRHS(const Array< int > &vdofs, const Vector &x, Vector &b)
Use the stored eliminated part of the matrix (see EliminateVDofs(const Array&lt;int&gt; &amp;...
const SparseMatrix & SpMatElim() const
Returns a reference to the sparse matrix of eliminated b.c.
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
void AddInteriorFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new interior Face Integrator. Assumes ownership of bfi.
FiniteElementSpace * fes
FE space on which the form lives. Not owned.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void ComputeElementMatrices()
Compute and store internally all element matrices.
FiniteElementSpace * SCFESpace() const
Return the trace FE space associated with static condensation.
void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1)
DenseMatrix elemmat
FiniteElementSpace * trial_fes
Not owned.
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
virtual ~BilinearForm()
Destroys bilinear form.
Vector data type.
Definition: vector.hpp:48
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds a boundary integrator. Assumes ownership of bfi.
virtual MatrixInverse * Inverse() const
Returns a pointer to (an approximation) of the matrix inverse.
Hybridization * hybridization
Owned.
void operator=(const double a)
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
const SparseMatrix & SpMat() const
Returns a reference to the sparse matrix.
Array< Array< int > * > bfbfi_marker
Entries are not owned.
Array< int > vdofs
void FreeElementMatrices()
Free the memory used by the element matrices.
Abstract operator.
Definition: operator.hpp:21
AssemblyLevel assembly
The form assembly level (full, partial, etc.)
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:656
Array< BilinearFormIntegrator * > * GetDI()
Access all interpolators added with AddDomainInterpolator().
void UseSparsity(int *I, int *J, bool isSorted)
Use the given CSR sparsity pattern to allocate the internal SparseMatrix.
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:779
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void EnableStaticCondensation()
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.