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