MFEM  v4.2.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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_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 /** @brief Enumeration defining the assembly level for bilinear and nonlinear
29  form classes derived from Operator. */
30 enum class AssemblyLevel
31 {
32  /// Legacy fully assembled form, i.e. a global sparse matrix in MFEM, Hypre
33  /// or PETSC format. This assembly is ALWAYS performed on the host.
34  LEGACYFULL = 0,
35  /// Fully assembled form, i.e. a global sparse matrix in MFEM format. This
36  /// assembly is compatible with device execution.
37  FULL,
38  /// Form assembled at element level, which computes and stores dense element
39  /// matrices.
40  ELEMENT,
41  /// Partially-assembled form, which computes and stores data only at
42  /// quadrature points.
43  PARTIAL,
44  /// "Matrix-free" form that computes all of its action on-the-fly without any
45  /// substantial storage.
46  NONE,
47 };
48 
49 
50 /** @brief A "square matrix" operator for the associated FE space and
51  BLFIntegrators The sum of all the BLFIntegrators can be used form the matrix
52  M. This class also supports other assembly levels specified via the
53  SetAssemblyLevel() function. */
54 class BilinearForm : public Matrix
55 {
56 protected:
57  /// Sparse matrix \f$ M \f$ to be associated with the form. Owned.
59 
60  /** @brief Sparse Matrix \f$ M_e \f$ used to store the eliminations
61  from the b.c. Owned.
62  \f$ M + M_e = M_{original} \f$ */
64 
65  /// FE space on which the form lives. Not owned.
67 
68  /// The assembly level of the form (full, partial, etc.)
70  /// Element batch size used in the form action (1, 8, num_elems, etc.)
71  int batch;
72  /** @brief Extension for supporting Full Assembly (FA), Element Assembly (EA),
73  Partial Assembly (PA), or Matrix Free assembly (MF). */
75 
76  /** @brief Indicates the Mesh::sequence corresponding to the current state of
77  the BilinearForm. */
78  long sequence;
79 
80  /** @brief Indicates the BilinearFormIntegrator%s stored in #dbfi, #bbfi,
81  #fbfi, and #bfbfi are owned by another BilinearForm. */
83 
84  /// Set of Domain Integrators to be applied.
86 
87  /// Set of Boundary Integrators to be applied.
89  Array<Array<int>*> bbfi_marker; ///< Entries are not owned.
90 
91  /// Set of interior face Integrators to be applied.
93 
94  /// Set of boundary face Integrators to be applied.
96  Array<Array<int>*> bfbfi_marker; ///< Entries are not owned.
97 
100 
102 
105 
106  /** This data member allows one to specify what should be done to the
107  diagonal matrix entries and corresponding RHS values upon elimination of
108  the constrained DoFs. */
110 
112  // Allocate appropriate SparseMatrix and assign it to mat
113  void AllocMat();
114 
115  void ConformingAssemble();
116 
117  // may be used in the construction of derived classes
119  {
120  fes = NULL; sequence = -1;
121  mat = mat_e = NULL; extern_bfs = 0; element_matrices = NULL;
122  static_cond = NULL; hybridization = NULL;
126  batch = 1;
127  ext = NULL;
128  }
129 
130 private:
131  /// Copy construction is not supported; body is undefined.
132  BilinearForm(const BilinearForm &);
133 
134  /// Copy assignment is not supported; body is undefined.
135  BilinearForm &operator=(const BilinearForm &);
136 
137 public:
138  /// Creates bilinear form associated with FE space @a *f.
139  /** The pointer @a f is not owned by the newly constructed object. */
141 
142  /** @brief Create a BilinearForm on the FiniteElementSpace @a f, using the
143  same integrators as the BilinearForm @a bf.
144 
145  The pointer @a f is not owned by the newly constructed object.
146 
147  The integrators in @a bf are copied as pointers and they are not owned by
148  the newly constructed BilinearForm.
149 
150  The optional parameter @a ps is used to initialize the internal flag
151  #precompute_sparsity, see UsePrecomputedSparsity() for details. */
152  BilinearForm(FiniteElementSpace *f, BilinearForm *bf, int ps = 0);
153 
154  /// Get the size of the BilinearForm as a square matrix.
155  int Size() const { return height; }
156 
157  /// Set the desired assembly level.
158  /** Valid choices are:
159 
160  - AssemblyLevel::LEGACYFULL (default)
161  - AssemblyLevel::FULL
162  - AssemblyLevel::PARTIAL
163  - AssemblyLevel::ELEMENT
164  - AssemblyLevel::NONE
165 
166  This method must be called before assembly. */
167  void SetAssemblyLevel(AssemblyLevel assembly_level);
168 
169  /// Returns the assembly level
171 
172  /** @brief Enable the use of static condensation. For details see the
173  description for class StaticCondensation in fem/staticcond.hpp This method
174  should be called before assembly. If the number of unknowns after static
175  condensation is not reduced, it is not enabled. */
177 
178  /** @brief Check if static condensation was actually enabled by a previous
179  call to EnableStaticCondensation(). */
180  bool StaticCondensationIsEnabled() const { return static_cond; }
181 
182  /// Return the trace FE space associated with static condensation.
184  { return static_cond ? static_cond->GetTraceFESpace() : NULL; }
185 
186  /// Enable hybridization.
187  /** For details see the description for class
188  Hybridization in fem/hybridization.hpp. This method should be called
189  before assembly. */
190  void EnableHybridization(FiniteElementSpace *constr_space,
191  BilinearFormIntegrator *constr_integ,
192  const Array<int> &ess_tdof_list);
193 
194  /** @brief For scalar FE spaces, precompute the sparsity pattern of the matrix
195  (assuming dense element matrices) based on the types of integrators
196  present in the bilinear form. */
197  void UsePrecomputedSparsity(int ps = 1) { precompute_sparsity = ps; }
198 
199  /** @brief Use the given CSR sparsity pattern to allocate the internal
200  SparseMatrix.
201 
202  - The @a I and @a J arrays must define a square graph with size equal to
203  GetVSize() of the associated FiniteElementSpace.
204  - This method should be called after enabling static condensation or
205  hybridization, if used.
206  - In the case of static condensation, @a I and @a J are not used.
207  - The ownership of the arrays @a I and @a J remains with the caller. */
208  void UseSparsity(int *I, int *J, bool isSorted);
209 
210  /// Use the sparsity of @a A to allocate the internal SparseMatrix.
211  void UseSparsity(SparseMatrix &A);
212 
213  /// Pre-allocate the internal SparseMatrix before assembly.
214  /** If the flag 'precompute sparsity'
215  is set, the matrix is allocated in CSR format (i.e.
216  finalized) and the entries are initialized with zeros. */
217  void AllocateMatrix() { if (mat == NULL) { AllocMat(); } }
218 
219  /// Access all the integrators added with AddDomainIntegrator().
221 
222  /// Access all the integrators added with AddBoundaryIntegrator().
224  /** @brief Access all boundary markers added with AddBoundaryIntegrator().
225  If no marker was specified when the integrator was added, the
226  corresponding pointer (to Array<int>) will be NULL. */
228 
229  /// Access all integrators added with AddInteriorFaceIntegrator().
231 
232  /// Access all integrators added with AddBdrFaceIntegrator().
234  /** @brief Access all boundary markers added with AddBdrFaceIntegrator().
235  If no marker was specified when the integrator was added, the
236  corresponding pointer (to Array<int>) will be NULL. */
238 
239  /// Returns a reference to: \f$ M_{ij} \f$
240  const double &operator()(int i, int j) { return (*mat)(i,j); }
241 
242  /// Returns a reference to: \f$ M_{ij} \f$
243  virtual double &Elem(int i, int j);
244 
245  /// Returns constant reference to: \f$ M_{ij} \f$
246  virtual const double &Elem(int i, int j) const;
247 
248  /// Matrix vector multiplication: \f$ y = M x \f$
249  virtual void Mult(const Vector &x, Vector &y) const;
250 
251  /** @brief Matrix vector multiplication with the original uneliminated
252  matrix. The original matrix is \f$ M + M_e \f$ so we have:
253  \f$ y = M x + M_e x \f$ */
254  void FullMult(const Vector &x, Vector &y) const
255  { mat->Mult(x, y); mat_e->AddMult(x, y); }
256 
257  /// Add the matrix vector multiple to a vector: \f$ y += a M x \f$
258  virtual void AddMult(const Vector &x, Vector &y, const double a = 1.0) const
259  { mat -> AddMult (x, y, a); }
260 
261  /** @brief Add the original uneliminated matrix vector multiple to a vector.
262  The original matrix is \f$ M + Me \f$ so we have:
263  \f$ y += M x + M_e x \f$ */
264  void FullAddMult(const Vector &x, Vector &y) const
265  { mat->AddMult(x, y); mat_e->AddMult(x, y); }
266 
267  /// Add the matrix transpose vector multiplication: \f$ y += a M^T x \f$
268  virtual void AddMultTranspose(const Vector & x, Vector & y,
269  const double a = 1.0) const
270  { mat->AddMultTranspose(x, y, a); }
271 
272  /** @brief Add the original uneliminated matrix transpose vector
273  multiple to a vector. The original matrix is \f$ M + M_e \f$
274  so we have: \f$ y += M^T x + {M_e}^T x \f$ */
275  void FullAddMultTranspose(const Vector & x, Vector & y) const
276  { mat->AddMultTranspose(x, y); mat_e->AddMultTranspose(x, y); }
277 
278  /// Matrix transpose vector multiplication: \f$ y = M^T x \f$
279  virtual void MultTranspose(const Vector & x, Vector & y) const
280  { y = 0.0; AddMultTranspose (x, y); }
281 
282  /// Compute \f$ y^T M x \f$
283  double InnerProduct(const Vector &x, const Vector &y) const
284  { return mat->InnerProduct (x, y); }
285 
286  /// Returns a pointer to (approximation) of the matrix inverse: \f$ M^{-1} \f$
287  virtual MatrixInverse *Inverse() const;
288 
289  /// Finalizes the matrix initialization.
290  virtual void Finalize(int skip_zeros = 1);
291 
292  /// Returns a const reference to the sparse matrix.
293  const SparseMatrix &SpMat() const
294  {
295  MFEM_VERIFY(mat, "mat is NULL and can't be dereferenced");
296  return *mat;
297  }
298 
299  /// Returns a reference to the sparse matrix: \f$ M \f$
301  {
302  MFEM_VERIFY(mat, "mat is NULL and can't be dereferenced");
303  return *mat;
304  }
305 
306  /** @brief Nullifies the internal matrix \f$ M \f$ and returns a pointer
307  to it. Used for transfering ownership. */
308  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
309 
310  /// Returns a const reference to the sparse matrix of eliminated b.c.: \f$ M_e \f$
311  const SparseMatrix &SpMatElim() const
312  {
313  MFEM_VERIFY(mat_e, "mat_e is NULL and can't be dereferenced");
314  return *mat_e;
315  }
316 
317  /// Returns a reference to the sparse matrix of eliminated b.c.: \f$ M_e \f$
319  {
320  MFEM_VERIFY(mat_e, "mat_e is NULL and can't be dereferenced");
321  return *mat_e;
322  }
323 
324  /// Adds new Domain Integrator. Assumes ownership of @a bfi.
326 
327  /// Adds new Boundary Integrator. Assumes ownership of @a bfi.
329 
330  /** @brief Adds new Boundary Integrator, restricted to specific boundary
331  attributes.
332 
333  Assumes ownership of @a bfi. The array @a bdr_marker is stored internally
334  as a pointer to the given Array<int> object. */
336  Array<int> &bdr_marker);
337 
338  /// Adds new interior Face Integrator. Assumes ownership of @a bfi.
340 
341  /// Adds new boundary Face Integrator. Assumes ownership of @a bfi.
343 
344  /** @brief Adds new boundary Face Integrator, restricted to specific boundary
345  attributes.
346 
347  Assumes ownership of @a bfi. The array @a bdr_marker is stored internally
348  as a pointer to the given Array<int> object. */
350  Array<int> &bdr_marker);
351 
352  /// Sets all sparse values of \f$ M \f$ and \f$ M_e \f$ to 'a'.
353  void operator=(const double a)
354  {
355  if (mat != NULL) { *mat = a; }
356  if (mat_e != NULL) { *mat_e = a; }
357  }
358 
359  /// Assembles the form i.e. sums over all domain/bdr integrators.
360  void Assemble(int skip_zeros = 1);
361 
362  /** @brief Assemble the diagonal of the bilinear form into diag
363 
364  For adaptively refined meshes, this returns P^T d_e, where d_e is the
365  locally assembled diagonal on each element and P^T is the transpose of
366  the conforming prolongation. In general this is not the correct diagonal
367  for an AMR mesh. */
368  void AssembleDiagonal(Vector &diag) const;
369 
370  /// Get the finite element space prolongation operator.
371  virtual const Operator *GetProlongation() const
372  { return fes->GetConformingProlongation(); }
373  /// Get the finite element space restriction operator
374  virtual const Operator *GetRestriction() const
375  { return fes->GetConformingRestriction(); }
376  /// Get the output finite element space prolongation matrix
377  virtual const Operator *GetOutputProlongation() const
378  { return GetProlongation(); }
379  /// Get the output finite element space restriction matrix
380  virtual const Operator *GetOutputRestriction() const
381  { return GetRestriction(); }
382 
383  /** @brief Form the linear system A X = B, corresponding to this bilinear
384  form and the linear form @a b(.). */
385  /** This method applies any necessary transformations to the linear system
386  such as: eliminating boundary conditions; applying conforming constraints
387  for non-conforming AMR; parallel assembly; static condensation;
388  hybridization.
389 
390  The GridFunction-size vector @a x must contain the essential b.c. The
391  BilinearForm and the LinearForm-size vector @a b must be assembled.
392 
393  The vector @a X is initialized with a suitable initial guess: when using
394  hybridization, the vector @a X is set to zero; otherwise, the essential
395  entries of @a X are set to the corresponding b.c. and all other entries
396  are set to zero (@a copy_interior == 0) or copied from @a x
397  (@a copy_interior != 0).
398 
399  This method can be called multiple times (with the same @a ess_tdof_list
400  array) to initialize different right-hand sides and boundary condition
401  values.
402 
403  After solving the linear system, the finite element solution @a x can be
404  recovered by calling RecoverFEMSolution() (with the same vectors @a X,
405  @a b, and @a x).
406 
407  NOTE: If there are no transformations, @a X simply reuses the data of
408  @a x. */
409  virtual void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x,
410  Vector &b, OperatorHandle &A, Vector &X,
411  Vector &B, int copy_interior = 0);
412 
413  /** @brief Form the linear system A X = B, corresponding to this bilinear
414  form and the linear form @a b(.). */
415  /** Version of the method FormLinearSystem() where the system matrix is
416  returned in the variable @a A, of type OpType, holding a *reference* to
417  the system matrix (created with the method OpType::MakeRef()). The
418  reference will be invalidated when SetOperatorType(), Update(), or the
419  destructor is called. */
420  template <typename OpType>
421  void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x, Vector &b,
422  OpType &A, Vector &X, Vector &B,
423  int copy_interior = 0)
424  {
425  OperatorHandle Ah;
426  FormLinearSystem(ess_tdof_list, x, b, Ah, X, B, copy_interior);
427  OpType *A_ptr = Ah.Is<OpType>();
428  MFEM_VERIFY(A_ptr, "invalid OpType used");
429  A.MakeRef(*A_ptr);
430  }
431 
432  /// Form the linear system matrix @a A, see FormLinearSystem() for details.
433  virtual void FormSystemMatrix(const Array<int> &ess_tdof_list,
434  OperatorHandle &A);
435 
436  /// Form the linear system matrix A, see FormLinearSystem() for details.
437  /** Version of the method FormSystemMatrix() where the system matrix is
438  returned in the variable @a A, of type OpType, holding a *reference* to
439  the system matrix (created with the method OpType::MakeRef()). The
440  reference will be invalidated when SetOperatorType(), Update(), or the
441  destructor is called. */
442  template <typename OpType>
443  void FormSystemMatrix(const Array<int> &ess_tdof_list, OpType &A)
444  {
445  OperatorHandle Ah;
446  FormSystemMatrix(ess_tdof_list, Ah);
447  OpType *A_ptr = Ah.Is<OpType>();
448  MFEM_VERIFY(A_ptr, "invalid OpType used");
449  A.MakeRef(*A_ptr);
450  }
451 
452  /// Recover the solution of a linear system formed with FormLinearSystem().
453  /** Call this method after solving a linear system constructed using the
454  FormLinearSystem() method to recover the solution as a GridFunction-size
455  vector in @a x. Use the same arguments as in the FormLinearSystem() call.
456  */
457  virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
458 
459  /// Compute and store internally all element matrices.
460  void ComputeElementMatrices();
461 
462  /// Free the memory used by the element matrices.
464  { delete element_matrices; element_matrices = NULL; }
465 
466  /// Compute the element matrix of the given element
467  /** The element matrix is computed by calling the domain integrators
468  or the one stored internally by a prior call of ComputeElementMatrices()
469  is returned when available.
470  */
471  void ComputeElementMatrix(int i, DenseMatrix &elmat);
472 
473  /// Compute the boundary element matrix of the given boundary element
474  void ComputeBdrElementMatrix(int i, DenseMatrix &elmat);
475 
476  /// Assemble the given element matrix
477  /** The element matrix @a elmat is assembled for the element @a i, i.e.
478  added to the system matrix. The flag @a skip_zeros skips the zero
479  elements of the matrix, unless they are breaking the symmetry of
480  the system matrix.
481  */
482  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
483  int skip_zeros = 1);
484 
485  /// Assemble the given element matrix
486  /** The element matrix @a elmat is assembled for the element @a i, i.e.
487  added to the system matrix. The vdofs of the element are returned
488  in @a vdofs. The flag @a skip_zeros skips the zero elements of the
489  matrix, unless they are breaking the symmetry of the system matrix.
490  */
491  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
492  Array<int> &vdofs, int skip_zeros = 1);
493 
494  /// Assemble the given boundary element matrix
495  /** The boundary element matrix @a elmat is assembled for the boundary
496  element @a i, i.e. added to the system matrix. The flag @a skip_zeros
497  skips the zero elements of the matrix, unless they are breaking the
498  symmetry of the system matrix.
499  */
500  void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat,
501  int skip_zeros = 1);
502 
503  /// Assemble the given boundary element matrix
504  /** The boundary element matrix @a elmat is assembled for the boundary
505  element @a i, i.e. added to the system matrix. The vdofs of the element
506  are returned in @a vdofs. The flag @a skip_zeros skips the zero elements
507  of the matrix, unless they are breaking the symmetry of the system matrix.
508  */
509  void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat,
510  Array<int> &vdofs, int skip_zeros = 1);
511 
512  /// Eliminate essential boundary DOFs from the system.
513  /** The array @a bdr_attr_is_ess marks boundary attributes that constitute
514  the essential part of the boundary. By default, the diagonal at the
515  essential DOFs is set to 1.0. This behavior is controlled by the argument
516  @a dpolicy. */
517  void EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
518  const Vector &sol, Vector &rhs,
519  DiagonalPolicy dpolicy = DIAG_ONE);
520 
521  /// Eliminate essential boundary DOFs from the system matrix.
522  void EliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
523  DiagonalPolicy dpolicy = DIAG_ONE);
524  /// Perform elimination and set the diagonal entry to the given value
525  void EliminateEssentialBCDiag(const Array<int> &bdr_attr_is_ess,
526  double value);
527 
528  /// Eliminate the given @a vdofs. NOTE: here, @a vdofs is a list of DOFs.
529  /** In this case the eliminations are applied to the internal \f$ M \f$
530  and @a rhs without storing the elimination matrix \f$ M_e \f$. */
531  void EliminateVDofs(const Array<int> &vdofs, const Vector &sol, Vector &rhs,
532  DiagonalPolicy dpolicy = DIAG_ONE);
533 
534  /// Eliminate the given @a vdofs, storing the eliminated part internally in \f$ M_e \f$.
535  /** This method works in conjunction with EliminateVDofsInRHS() and allows
536  elimination of boundary conditions in multiple right-hand sides. In this
537  method, @a vdofs is a list of DOFs. */
538  void EliminateVDofs(const Array<int> &vdofs,
539  DiagonalPolicy dpolicy = DIAG_ONE);
540 
541  /** @brief Similar to
542  EliminateVDofs(const Array<int> &, const Vector &, Vector &, DiagonalPolicy)
543  but here @a ess_dofs is a marker (boolean) array on all vector-dofs
544  (@a ess_dofs[i] < 0 is true). */
545  void EliminateEssentialBCFromDofs(const Array<int> &ess_dofs, const Vector &sol,
546  Vector &rhs, DiagonalPolicy dpolicy = DIAG_ONE);
547 
548  /** @brief Similar to EliminateVDofs(const Array<int> &, DiagonalPolicy) but
549  here @a ess_dofs is a marker (boolean) array on all vector-dofs
550  (@a ess_dofs[i] < 0 is true). */
551  void EliminateEssentialBCFromDofs(const Array<int> &ess_dofs,
552  DiagonalPolicy dpolicy = DIAG_ONE);
553  /// Perform elimination and set the diagonal entry to the given value
554  void EliminateEssentialBCFromDofsDiag(const Array<int> &ess_dofs,
555  double value);
556 
557  /** @brief Use the stored eliminated part of the matrix (see
558  EliminateVDofs(const Array<int> &, DiagonalPolicy)) to modify the r.h.s.
559  @a b; @a vdofs is a list of DOFs (non-directional, i.e. >= 0). */
560  void EliminateVDofsInRHS(const Array<int> &vdofs, const Vector &x,
561  Vector &b);
562 
563  /// Compute inner product for full uneliminated matrix \f$ y^T M x + y^T M_e x \f$
564  double FullInnerProduct(const Vector &x, const Vector &y) const
565  { return mat->InnerProduct(x, y) + mat_e->InnerProduct(x, y); }
566 
567  /// Update the @a FiniteElementSpace and delete all data associated with the old one.
568  virtual void Update(FiniteElementSpace *nfes = NULL);
569 
570  /// (DEPRECATED) Return the FE space associated with the BilinearForm.
571  /** @deprecated Use FESpace() instead. */
572  MFEM_DEPRECATED FiniteElementSpace *GetFES() { return fes; }
573 
574  /// Return the FE space associated with the BilinearForm.
576  /// Read-only access to the associated FiniteElementSpace.
577  const FiniteElementSpace *FESpace() const { return fes; }
578 
579  /// Sets diagonal policy used upon construction of the linear system.
580  /** Policies include:
581 
582  - DIAG_ZERO (Set the diagonal values to zero)
583  - DIAG_ONE (Set the diagonal values to one)
584  - DIAG_KEEP (Keep the diagonal values)
585  */
586  void SetDiagonalPolicy(DiagonalPolicy policy);
587 
588  /// Indicate that integrators are not owned by the BilinearForm
590 
591  /// Destroys bilinear form.
592  virtual ~BilinearForm();
593 };
594 
595 
596 /**
597  Class for assembling of bilinear forms `a(u,v)` defined on different
598  trial and test spaces. The assembled matrix `M` is such that
599 
600  a(u,v) = V^t M U
601 
602  where `U` and `V` are the vectors representing the functions `u` and `v`,
603  respectively. The first argument, `u`, of `a(,)` is in the trial space
604  and the second argument, `v`, is in the test space. Thus,
605 
606  # of rows of M = dimension of the test space and
607  # of cols of M = dimension of the trial space.
608 
609  Both trial and test spaces should be defined on the same mesh.
610 */
611 class MixedBilinearForm : public Matrix
612 {
613 protected:
614  SparseMatrix *mat; ///< Owned.
615  SparseMatrix *mat_e; ///< Owned.
616 
617  FiniteElementSpace *trial_fes, ///< Not owned
618  *test_fes; ///< Not owned
619 
620  /// The form assembly level (full, partial, etc.)
622 
623  /** Extension for supporting Full Assembly (FA), Element Assembly (EA),
624  Partial Assembly (PA), or Matrix Free assembly (MF). */
626 
627  /** @brief Indicates the BilinearFormIntegrator%s stored in #dbfi, #bbfi,
628  #tfbfi and #btfbfi are owned by another MixedBilinearForm. */
630 
631  /// Domain integrators.
633 
634  /// Boundary integrators.
636  Array<Array<int>*> bbfi_marker;///< Entries are not owned.
637 
638  /// Trace face (skeleton) integrators.
640 
641  /// Boundary trace face (skeleton) integrators.
643  Array<Array<int>*> btfbfi_marker;///< Entries are not owned.
644 
647 
648 private:
649  /// Copy construction is not supported; body is undefined.
651 
652  /// Copy assignment is not supported; body is undefined.
653  MixedBilinearForm &operator=(const MixedBilinearForm &);
654 
655 public:
656  /** @brief Construct a MixedBilinearForm on the given trial, @a tr_fes, and
657  test, @a te_fes, FiniteElementSpace%s. */
658  /** The pointers @a tr_fes and @a te_fes are not owned by the newly
659  constructed object. */
661  FiniteElementSpace *te_fes);
662 
663  /** @brief Create a MixedBilinearForm on the given trial, @a tr_fes, and
664  test, @a te_fes, FiniteElementSpace%s, using the same integrators as the
665  MixedBilinearForm @a mbf.
666 
667  The pointers @a tr_fes and @a te_fes are not owned by the newly
668  constructed object.
669 
670  The integrators in @a mbf are copied as pointers and they are not owned
671  by the newly constructed MixedBilinearForm. */
673  FiniteElementSpace *te_fes,
674  MixedBilinearForm *mbf);
675 
676  /// Returns a reference to: \f$ M_{ij} \f$
677  virtual double &Elem(int i, int j);
678 
679  /// Returns a reference to: \f$ M_{ij} \f$
680  virtual const double &Elem(int i, int j) const;
681 
682  /// Matrix multiplication: \f$ y = M x \f$
683  virtual void Mult(const Vector & x, Vector & y) const;
684 
685  virtual void AddMult(const Vector & x, Vector & y,
686  const double a = 1.0) const;
687 
688  virtual void MultTranspose(const Vector & x, Vector & y) const;
689  virtual void AddMultTranspose(const Vector & x, Vector & y,
690  const double a = 1.0) const;
691 
692  virtual MatrixInverse *Inverse() const;
693 
694  /// Finalizes the matrix initialization.
695  virtual void Finalize(int skip_zeros = 1);
696 
697  /** Extract the associated matrix as SparseMatrix blocks. The number of
698  block rows and columns is given by the vector dimensions (vdim) of the
699  test and trial spaces, respectively. */
700  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
701 
702  /// Returns a const reference to the sparse matrix: \f$ M \f$
703  const SparseMatrix &SpMat() const { return *mat; }
704 
705  /// Returns a reference to the sparse matrix: \f$ M \f$
706  SparseMatrix &SpMat() { return *mat; }
707 
708  /** @brief Nullifies the internal matrix \f$ M \f$ and returns a pointer
709  to it. Used for transfering ownership. */
710  SparseMatrix *LoseMat() { SparseMatrix *tmp = mat; mat = NULL; return tmp; }
711 
712  /// Adds a domain integrator. Assumes ownership of @a bfi.
714 
715  /// Adds a boundary integrator. Assumes ownership of @a bfi.
717 
718  /// Adds a boundary integrator. Assumes ownership of @a bfi.
720  Array<int> &bdr_marker);
721 
722  /** @brief Add a trace face integrator. Assumes ownership of @a bfi.
723 
724  This type of integrator assembles terms over all faces of the mesh using
725  the face FE from the trial space and the two adjacent volume FEs from the
726  test space. */
728 
729  /// Adds a boundary trace face integrator. Assumes ownership of @a bfi.
731 
732  /// Adds a boundary trace face integrator. Assumes ownership of @a bfi.
734  Array<int> &bdr_marker);
735 
736  /// Access all integrators added with AddDomainIntegrator().
738 
739  /// Access all integrators added with AddBoundaryIntegrator().
741  /** @brief Access all boundary markers added with AddBoundaryIntegrator().
742  If no marker was specified when the integrator was added, the
743  corresponding pointer (to Array<int>) will be NULL. */
745 
746  /// Access all integrators added with AddTraceFaceIntegrator().
748 
749  /// Access all integrators added with AddBdrTraceFaceIntegrator().
751  /** @brief Access all boundary markers added with AddBdrTraceFaceIntegrator().
752  If no marker was specified when the integrator was added, the
753  corresponding pointer (to Array<int>) will be NULL. */
755 
756  /// Sets all sparse values of \f$ M \f$ to @a a.
757  void operator=(const double a) { *mat = a; }
758 
759  /// Set the desired assembly level. The default is AssemblyLevel::LEGACYFULL.
760  /** This method must be called before assembly. */
761  void SetAssemblyLevel(AssemblyLevel assembly_level);
762 
763  void Assemble(int skip_zeros = 1);
764 
765  /** @brief Assemble the diagonal of ADA^T into diag, where A is this mixed
766  bilinear form and D is a diagonal. */
767  void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const;
768 
769  /// Get the input finite element space prolongation matrix
770  virtual const Operator *GetProlongation() const
771  { return trial_fes->GetProlongationMatrix(); }
772 
773  /// Get the input finite element space restriction matrix
774  virtual const Operator *GetRestriction() const
775  { return trial_fes->GetRestrictionMatrix(); }
776 
777  /// Get the test finite element space prolongation matrix
778  virtual const Operator *GetOutputProlongation() const
779  { return test_fes->GetProlongationMatrix(); }
780 
781  /// Get the test finite element space restriction matrix
782  virtual const Operator *GetOutputRestriction() const
783  { return test_fes->GetRestrictionMatrix(); }
784 
785  /** For partially conforming trial and/or test FE spaces, complete the
786  assembly process by performing A := P2^t A P1 where A is the internal
787  sparse matrix; P1 and P2 are the conforming prolongation matrices of the
788  trial and test FE spaces, respectively. After this call the
789  MixedBilinearForm becomes an operator on the conforming FE spaces. */
790  void ConformingAssemble();
791 
792  /// Compute the element matrix of the given element
793  void ComputeElementMatrix(int i, DenseMatrix &elmat);
794 
795  /// Compute the boundary element matrix of the given boundary element
796  void ComputeBdrElementMatrix(int i, DenseMatrix &elmat);
797 
798  /// Assemble the given element matrix
799  /** The element matrix @a elmat is assembled for the element @a i, i.e.
800  added to the system matrix. The flag @a skip_zeros skips the zero
801  elements of the matrix, unless they are breaking the symmetry of
802  the system matrix.
803  */
804  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
805  int skip_zeros = 1);
806 
807  /// Assemble the given element matrix
808  /** The element matrix @a elmat is assembled for the element @a i, i.e.
809  added to the system matrix. The vdofs of the element are returned
810  in @a trial_vdofs and @a test_vdofs. The flag @a skip_zeros skips
811  the zero elements of the matrix, unless they are breaking the symmetry
812  of the system matrix.
813  */
814  void AssembleElementMatrix(int i, const DenseMatrix &elmat,
816  int skip_zeros = 1);
817 
818  /// Assemble the given boundary element matrix
819  /** The boundary element matrix @a elmat is assembled for the boundary
820  element @a i, i.e. added to the system matrix. The flag @a skip_zeros
821  skips the zero elements of the matrix, unless they are breaking the
822  symmetry of the system matrix.
823  */
824  void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat,
825  int skip_zeros = 1);
826 
827  /// Assemble the given boundary element matrix
828  /** The boundary element matrix @a elmat is assembled for the boundary
829  element @a i, i.e. added to the system matrix. The vdofs of the element
830  are returned in @a trial_vdofs and @a test_vdofs. The flag @a skip_zeros
831  skips the zero elements of the matrix, unless they are breaking the
832  symmetry of the system matrix.
833  */
834  void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat,
836  int skip_zeros = 1);
837 
838  void EliminateTrialDofs(const Array<int> &bdr_attr_is_ess,
839  const Vector &sol, Vector &rhs);
840 
841  void EliminateEssentialBCFromTrialDofs(const Array<int> &marked_vdofs,
842  const Vector &sol, Vector &rhs);
843 
844  virtual void EliminateTestDofs(const Array<int> &bdr_attr_is_ess);
845 
846  /** @brief Return in @a A that is column-constrained.
847 
848  This returns the same operator as FormRectangularLinearSystem(), but does
849  without the transformations of the right-hand side. */
850  void FormRectangularSystemMatrix(const Array<int> &trial_tdof_list,
851  const Array<int> &test_tdof_list,
852  OperatorHandle &A);
853 
854  /** @brief Form the column-constrained linear system matrix A.
855  See FormRectangularSystemMatrix() for details.
856 
857  Version of the method FormRectangularSystemMatrix() where the system matrix is
858  returned in the variable @a A, of type OpType, holding a *reference* to
859  the system matrix (created with the method OpType::MakeRef()). The
860  reference will be invalidated when SetOperatorType(), Update(), or the
861  destructor is called. */
862  template <typename OpType>
863  void FormRectangularSystemMatrix(const Array<int> &trial_tdof_list,
864  const Array<int> &test_tdof_list, OpType &A)
865  {
866  OperatorHandle Ah;
867  FormRectangularSystemMatrix(trial_tdof_list, test_tdof_list, Ah);
868  OpType *A_ptr = Ah.Is<OpType>();
869  MFEM_VERIFY(A_ptr, "invalid OpType used");
870  A.MakeRef(*A_ptr);
871  }
872 
873  /** @brief Form the linear system A X = B, corresponding to this mixed bilinear
874  form and the linear form @a b(.).
875 
876  Return in @a A a *reference* to the system matrix that is column-constrained.
877  The reference will be invalidated when SetOperatorType(), Update(), or the
878  destructor is called. */
879  void FormRectangularLinearSystem(const Array<int> &trial_tdof_list,
880  const Array<int> &test_tdof_list,
881  Vector &x, Vector &b,
882  OperatorHandle &A, Vector &X, Vector &B);
883 
884  /** @brief Form the linear system A X = B, corresponding to this bilinear
885  form and the linear form @a b(.).
886 
887  Version of the method FormRectangularLinearSystem() where the system matrix is
888  returned in the variable @a A, of type OpType, holding a *reference* to
889  the system matrix (created with the method OpType::MakeRef()). The
890  reference will be invalidated when SetOperatorType(), Update(), or the
891  destructor is called. */
892  template <typename OpType>
893  void FormRectangularLinearSystem(const Array<int> &trial_tdof_list,
894  const Array<int> &test_tdof_list,
895  Vector &x, Vector &b,
896  OpType &A, Vector &X, Vector &B)
897  {
898  OperatorHandle Ah;
899  FormRectangularLinearSystem(trial_tdof_list, test_tdof_list, x, b, Ah, X, B);
900  OpType *A_ptr = Ah.Is<OpType>();
901  MFEM_VERIFY(A_ptr, "invalid OpType used");
902  A.MakeRef(*A_ptr);
903  }
904 
905  void Update();
906 
907  /// Return the trial FE space associated with the BilinearForm.
909  /// Read-only access to the associated trial FiniteElementSpace.
910  const FiniteElementSpace *TrialFESpace() const { return trial_fes; }
911 
912  /// Return the test FE space associated with the BilinearForm.
914  /// Read-only access to the associated test FiniteElementSpace.
915  const FiniteElementSpace *TestFESpace() const { return test_fes; }
916 
917  virtual ~MixedBilinearForm();
918 };
919 
920 
921 /**
922  Class for constructing the matrix representation of a linear operator,
923  `v = L u`, from one FiniteElementSpace (domain) to another FiniteElementSpace
924  (range). The constructed matrix `A` is such that
925 
926  V = A U
927 
928  where `U` and `V` are the vectors of degrees of freedom representing the
929  functions `u` and `v`, respectively. The dimensions of `A` are
930 
931  number of rows of A = dimension of the range space and
932  number of cols of A = dimension of the domain space.
933 
934  This class is very similar to MixedBilinearForm. One difference is that
935  the linear operator `L` is defined using a special kind of
936  BilinearFormIntegrator (we reuse its functionality instead of defining a
937  new class). The other difference with the MixedBilinearForm class is that
938  the "assembly" process overwrites the global matrix entries using the
939  local element matrices instead of adding them.
940 
941  Note that if we define the bilinear form `b(u,v) := (Lu,v)` using an inner
942  product in the range space, then its matrix representation, `B`, is
943 
944  B = M A, (since V^t B U = b(u,v) = (Lu,v) = V^t M A U)
945 
946  where `M` denotes the mass matrix for the inner product in the range space:
947  `V1^t M V2 = (v1,v2)`. Similarly, if `c(u,w) := (Lu,Lw)` then
948 
949  C = A^t M A.
950 */
952 {
953 private:
954  /// Copy construction is not supported; body is undefined.
956 
957  /// Copy assignment is not supported; body is undefined.
959 
960 public:
961  /** @brief Construct a DiscreteLinearOperator on the given
962  FiniteElementSpace%s @a domain_fes and @a range_fes. */
963  /** The pointers @a domain_fes and @a range_fes are not owned by the newly
964  constructed object. */
966  FiniteElementSpace *range_fes)
967  : MixedBilinearForm(domain_fes, range_fes) { }
968 
969  /// Adds a domain interpolator. Assumes ownership of @a di.
971  { AddDomainIntegrator(di); }
972 
973  /// Adds a trace face interpolator. Assumes ownership of @a di.
975  { AddTraceFaceIntegrator(di); }
976 
977  /// Access all interpolators added with AddDomainInterpolator().
979 
980  /** @brief Construct the internal matrix representation of the discrete
981  linear operator. */
982  virtual void Assemble(int skip_zeros = 1);
983 };
984 
985 }
986 
987 #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
Add the original uneliminated matrix vector multiple to a vector. The original matrix is so we have:...
SparseMatrix & SpMat()
Returns a reference to the sparse matrix: .
Array< BilinearFormIntegrator * > * GetBBFI()
Access all the integrators added with AddBoundaryIntegrator().
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication: .
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
virtual const Operator * GetOutputRestriction() const
Get the output finite element space restriction matrix.
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:874
virtual MatrixInverse * Inverse() const
Returns a pointer to (approximation) of the matrix inverse: .
Array< BilinearFormIntegrator * > fbfi
Set of interior face Integrators to be applied.
MFEM_DEPRECATED FiniteElementSpace * GetFES()
(DEPRECATED) Return the FE space associated with the BilinearForm.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
virtual const Operator * GetProlongation() const
Get the input finite element space prolongation matrix.
virtual void EliminateTestDofs(const Array< int > &bdr_attr_is_ess)
void EliminateTrialDofs(const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs)
void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Return in A that is column-constrained.
FiniteElementSpace * FESpace()
Return the FE space associated with the BilinearForm.
BilinearFormExtension * ext
Extension for supporting Full Assembly (FA), Element Assembly (EA), Partial Assembly (PA)...
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
Array< BilinearFormIntegrator * > dbfi
Set of Domain Integrators to be applied.
virtual const Operator * GetOutputProlongation() const
Get the output finite element space prolongation matrix.
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
Add the matrix vector multiple to a vector: .
Array< BilinearFormIntegrator * > * GetBFBFI()
Access all integrators added with AddBdrFaceIntegrator().
AssemblyLevel GetAssemblyLevel() const
Returns the assembly level.
Array< Array< int > * > bbfi_marker
Entries are not owned.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
Array< Array< int > * > * GetBTFBFI_Marker()
Access all boundary markers added with AddBdrTraceFaceIntegrator(). If no marker was specified when t...
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().
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:601
int batch
Element batch size used in the form action (1, 8, num_elems, etc.)
AssemblyLevel assembly
The form assembly level (full, partial, 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:728
void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag.
Abstract data type for matrix inverse.
Definition: matrix.hpp:62
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
Array< BilinearFormIntegrator * > btfbfi
Boundary trace face (skeleton) integrators.
Array< BilinearFormIntegrator * > bbfi
Boundary integrators.
virtual const Operator * GetProlongation() const
Get the finite element space prolongation operator.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
Array< BilinearFormIntegrator * > * GetFBFI()
Access all integrators added with AddInteriorFaceIntegrator().
virtual const Operator * GetRestriction() const
Get the input finite element space restriction matrix.
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
Compute inner product for full uneliminated matrix .
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.
Class extending the MixedBilinearForm class to support different AssemblyLevels.
void ComputeBdrElementMatrix(int i, DenseMatrix &elmat)
Compute the boundary element matrix of the given boundary element.
void EnableHybridization(FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list)
Enable hybridization.
void AddBdrTraceFaceIntegrator(BilinearFormIntegrator *bfi)
Adds a boundary trace face integrator. Assumes ownership of bfi.
void AssembleElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given element matrix.
void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OpType &A)
Form the column-constrained linear system matrix A. See FormRectangularSystemMatrix() for details...
void EliminateEssentialBCFromDofsDiag(const Array< int > &ess_dofs, double value)
Perform elimination and set the diagonal entry to the given value.
long sequence
Indicates the Mesh::sequence corresponding to the current state of the BilinearForm.
void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given boundary element matrix.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
SparseMatrix & SpMat()
Returns a reference to the sparse matrix: .
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:881
void FullAddMultTranspose(const Vector &x, Vector &y) const
Add the original uneliminated matrix transpose vector multiple to a vector. The original matrix is s...
OpType * Is() const
Return the Operator pointer dynamically cast to a specified OpType.
Definition: handle.hpp:100
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:46
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 FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
Form the linear system A X = B, corresponding to this mixed bilinear form and the linear form b(...
double b
Definition: lissajous.cpp:42
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 AssembleElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given element matrix.
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.
FiniteElementSpace * TestFESpace()
Return the test FE space associated with the BilinearForm.
void Assemble(int skip_zeros=1)
virtual const Operator * GetOutputRestriction() const
Get the test finite element space restriction matrix.
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:330
void FullMult(const Vector &x, Vector &y) const
Matrix vector multiplication with the original uneliminated matrix. The original matrix is so we hav...
void EliminateEssentialBCDiag(const Array< int > &bdr_attr_is_ess, double value)
Perform elimination and set the diagonal entry to the given value.
void AllocateMatrix()
Pre-allocate the internal SparseMatrix before assembly.
Array< BilinearFormIntegrator * > * GetBTFBFI()
Access all integrators added with AddBdrTraceFaceIntegrator().
Class extending the BilinearForm class to support different AssemblyLevels.
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
void FormRectangularLinearSystem(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OpType &A, Vector &X, Vector &B)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
const FiniteElementSpace * TestFESpace() const
Read-only access to the associated test FiniteElementSpace.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T into diag, where A is this mixed bilinear form and D is a diagonal...
void ComputeElementMatrix(int i, DenseMatrix &elmat)
Compute the element matrix of the given element.
SparseMatrix * mat
Owned.
Set the diagonal value to one.
Definition: operator.hpp:48
void AssembleBdrElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given boundary element matrix.
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().
Dynamic 2D array using row-major layout.
Definition: array.hpp:333
Array< BilinearFormIntegrator * > * GetBBFI()
Access all integrators added with AddBoundaryIntegrator().
void ComputeElementMatrix(int i, DenseMatrix &elmat)
Compute the element matrix of the given element.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACYFULL.
SparseMatrix & SpMatElim()
Returns a reference to the sparse matrix of eliminated b.c.: .
SparseMatrix * mat_e
Sparse Matrix used to store the eliminations from the b.c. Owned. .
void UseExternalIntegrators()
Indicate that integrators are not owned by the BilinearForm.
bool StaticCondensationIsEnabled() const
Check if static condensation was actually enabled by a previous call to EnableStaticCondensation().
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:31
virtual const Operator * GetRestriction() const
Get the finite element space restriction operator.
void UsePrecomputedSparsity(int ps=1)
For scalar FE spaces, precompute the sparsity pattern of the matrix (assuming dense element matrices)...
Array< Array< int > * > btfbfi_marker
Entries are not owned.
virtual void MultTranspose(const Vector &x, Vector &y) const
Matrix transpose vector multiplication: .
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
MixedBilinearFormExtension * ext
SparseMatrix * LoseMat()
Nullifies the internal matrix and returns a pointer to it. Used for transfering ownership.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:594
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()
Nullifies the internal matrix and returns a pointer to it. Used for transfering ownership.
void operator=(const double a)
Sets all sparse values of and to &#39;a&#39;.
const FiniteElementSpace * TrialFESpace() const
Read-only access to the associated trial FiniteElementSpace.
virtual void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
Add the matrix transpose vector multiplication: .
void AddTraceFaceInterpolator(DiscreteInterpolator *di)
Adds a trace face interpolator. Assumes ownership of di.
virtual double & Elem(int i, int j)
Returns a reference to: .
double InnerProduct(const Vector &x, const Vector &y) const
Compute .
DiagonalPolicy diag_policy
Array< Array< int > * > bbfi_marker
Entries are not owned.
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
void EliminateEssentialBCFromTrialDofs(const Array< int > &marked_vdofs, const Vector &sol, Vector &rhs)
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
FiniteElementSpace * test_fes
Not owned.
double a
Definition: lissajous.cpp:41
Array< Array< int > * > * GetBBFI_Marker()
Access all boundary markers added with AddBoundaryIntegrator(). If no marker was specified when the i...
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
const double & operator()(int i, int j)
Returns a reference to: .
FiniteElementSpace * GetTraceFESpace()
Return a pointer to the reduced/trace FE space.
Definition: staticcond.hpp:109
int extern_bfs
Indicates the BilinearFormIntegrators stored in dbfi, bbfi, tfbfi and btfbfi are owned by another Mix...
virtual const Operator * GetOutputProlongation() const
Get the test finite element space prolongation matrix.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
SparseMatrix * mat_e
Owned.
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;...
Keep the diagonal value.
Definition: operator.hpp:49
const SparseMatrix & SpMatElim() const
Returns a const reference to the sparse matrix of eliminated b.c.: .
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:45
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.
void ComputeElementMatrices()
Compute and store internally all element matrices.
FiniteElementSpace * SCFESpace() const
Return the trace FE space associated with static condensation.
DenseMatrix elemmat
FiniteElementSpace * trial_fes
Not owned.
virtual double & Elem(int i, int j)
Returns a reference to: .
virtual ~BilinearForm()
Destroys bilinear form.
Vector data type.
Definition: vector.hpp:51
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds a boundary 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 ...
virtual MatrixInverse * Inverse() const
Returns a pointer to (an approximation) of the matrix inverse.
Hybridization * hybridization
Owned.
void operator=(const double a)
Sets all sparse values of to a.
Array< BilinearFormIntegrator * > * GetDBFI()
Access all integrators added with AddDomainIntegrator().
void ComputeBdrElementMatrix(int i, DenseMatrix &elmat)
Compute the boundary element matrix of the given boundary element.
const SparseMatrix & SpMat() const
Returns a const 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:24
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:334
Array< BilinearFormIntegrator * > dbfi
Domain integrators.
AssemblyLevel assembly
The assembly level of the form (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:728
FiniteElementSpace * TrialFESpace()
Return the trial FE space associated with the BilinearForm.
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.
Array< BilinearFormIntegrator * > tfbfi
Trace face (skeleton) integrators.
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:979
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
double f(const Vector &p)