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