MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pbilinearform.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_PBILINEARFORM
13 #define MFEM_PBILINEARFORM
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include <mpi.h>
20 #include "pfespace.hpp"
21 #include "pgridfunc.hpp"
22 #include "bilinearform.hpp"
23 
24 namespace mfem
25 {
26 
27 /// Class for parallel bilinear form
29 {
30 protected:
31  ParFiniteElementSpace *pfes; ///< Points to the same object as #fes
32 
33  /// Auxiliary objects used in TrueAddMult().
34  mutable ParGridFunction X, Y;
35 
37 
39 
40  // Allocate mat - called when (mat == NULL && fbfi.Size() > 0)
41  void pAllocMat();
42 
43  void AssembleSharedFaces(int skip_zeros = 1);
44 
45 private:
46  /// Copy construction is not supported; body is undefined.
48 
49  /// Copy assignment is not supported; body is undefined.
50  ParBilinearForm &operator=(const ParBilinearForm &);
51 
52 public:
53  /// Creates parallel bilinear form associated with the FE space @a *pf.
54  /** The pointer @a pf is not owned by the newly constructed object. */
56  : BilinearForm(pf), pfes(pf),
58  { keep_nbr_block = false; }
59 
60  /** @brief Create a ParBilinearForm on the ParFiniteElementSpace @a *pf,
61  using the same integrators as the ParBilinearForm @a *bf.
62 
63  The pointer @a pf is not owned by the newly constructed object.
64 
65  The integrators in @a bf are copied as pointers and they are not owned by
66  the newly constructed ParBilinearForm. */
68  : BilinearForm(pf, bf), pfes(pf),
70  { keep_nbr_block = false; }
71 
72  /** When set to true and the ParBilinearForm has interior face integrators,
73  the local SparseMatrix will include the rows (in addition to the columns)
74  corresponding to face-neighbor dofs. The default behavior is to disregard
75  those rows. Must be called before the first Assemble call. */
76  void KeepNbrBlock(bool knb = true) { keep_nbr_block = knb; }
77 
78  /** @brief Set the operator type id for the parallel matrix/operator when
79  using AssemblyLevel::LEGACYFULL. */
80  /** If using static condensation or hybridization, call this method *after*
81  enabling it. */
83  {
84  p_mat.SetType(tid); p_mat_e.SetType(tid);
87  }
88 
89  /// Assemble the local matrix
90  void Assemble(int skip_zeros = 1);
91 
92  /// Returns the matrix assembled on the true dofs, i.e. P^t A P.
93  /** The returned matrix has to be deleted by the caller. */
95 
96  /// Returns the eliminated matrix assembled on the true dofs, i.e. P^t A_e P.
97  /** The returned matrix has to be deleted by the caller. */
99 
100  /// Return the matrix @a m assembled on the true dofs, i.e. P^t A P.
101  /** The returned matrix has to be deleted by the caller. */
103 
104  /** @brief Returns the matrix assembled on the true dofs, i.e.
105  @a A = P^t A_local P, in the format (type id) specified by @a A. */
107 
108  /** Returns the eliminated matrix assembled on the true dofs, i.e.
109  @a A_elim = P^t A_elim_local P in the format (type id) specified by @a A.
110  */
112  { ParallelAssemble(A_elim, mat_e); }
113 
114  /** Returns the matrix @a A_local assembled on the true dofs, i.e.
115  @a A = P^t A_local P in the format (type id) specified by @a A. */
116  void ParallelAssemble(OperatorHandle &A, SparseMatrix *A_local);
117 
118  /// Eliminate essential boundary DOFs from a parallel assembled system.
119  /** The array @a bdr_attr_is_ess marks boundary attributes that constitute
120  the essential part of the boundary. */
121  void ParallelEliminateEssentialBC(const Array<int> &bdr_attr_is_ess,
122  HypreParMatrix &A,
123  const HypreParVector &X,
124  HypreParVector &B) const;
125 
126  /// Eliminate essential boundary DOFs from a parallel assembled matrix @a A.
127  /** The array @a bdr_attr_is_ess marks boundary attributes that constitute
128  the essential part of the boundary. The eliminated part is stored in a
129  matrix A_elim such that A_original = A_new + A_elim. Returns a pointer to
130  the newly allocated matrix A_elim which should be deleted by the caller.
131  The matrices @a A and A_elim can be used to eliminate boundary conditions
132  in multiple right-hand sides, by calling the function EliminateBC() (from
133  hypre.hpp). */
135  HypreParMatrix &A) const;
136 
137  /// Eliminate essential true DOFs from a parallel assembled matrix @a A.
138  /** Given a list of essential true dofs and the parallel assembled matrix
139  @a A, eliminate the true dofs from the matrix, storing the eliminated
140  part in a matrix A_elim such that A_original = A_new + A_elim. Returns a
141  pointer to the newly allocated matrix A_elim which should be deleted by
142  the caller. The matrices @a A and A_elim can be used to eliminate
143  boundary conditions in multiple right-hand sides, by calling the function
144  EliminateBC() (from hypre.hpp). */
146  HypreParMatrix &A) const
147  { return A.EliminateRowsCols(tdofs_list); }
148 
149  /** @brief Compute @a y += @a a (P^t A P) @a x, where @a x and @a y are
150  vectors on the true dofs. */
151  void TrueAddMult(const Vector &x, Vector &y, const double a = 1.0) const;
152 
153  /// Return the parallel FE space associated with the ParBilinearForm.
154  ParFiniteElementSpace *ParFESpace() const { return pfes; }
155 
156  /// Return the parallel trace FE space associated with static condensation.
158  { return static_cond ? static_cond->GetParTraceFESpace() : NULL; }
159 
160  /// Get the parallel finite element space prolongation matrix
161  virtual const Operator *GetProlongation() const
162  { return pfes->GetProlongationMatrix(); }
163  /// Get the parallel finite element space restriction matrix
164  virtual const Operator *GetRestriction() const
165  { return pfes->GetRestrictionMatrix(); }
166 
169 
170  virtual void FormLinearSystem(const Array<int> &ess_tdof_list, Vector &x,
171  Vector &b, OperatorHandle &A, Vector &X,
172  Vector &B, int copy_interior = 0);
173 
174  virtual void FormSystemMatrix(const Array<int> &ess_tdof_list,
175  OperatorHandle &A);
176 
177  /** Call this method after solving a linear system constructed using the
178  FormLinearSystem method to recover the solution as a ParGridFunction-size
179  vector in x. Use the same arguments as in the FormLinearSystem call. */
180  virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x);
181 
182  virtual void Update(FiniteElementSpace *nfes = NULL);
183 
184  virtual ~ParBilinearForm() { }
185 };
186 
187 /// Class for parallel bilinear form using different test and trial FE spaces.
189 {
190 protected:
191  /// Points to the same object as #trial_fes
193  /// Points to the same object as #test_fes
195  /// Auxiliary objects used in TrueAddMult().
196  mutable ParGridFunction X, Y;
197 
198  /// Matrix and eliminated matrix
200 
201 private:
202  /// Copy construction is not supported; body is undefined.
204 
205  /// Copy assignment is not supported; body is undefined.
206  ParMixedBilinearForm &operator=(const ParMixedBilinearForm &);
207 
208 public:
209  /** @brief Construct a ParMixedBilinearForm on the given FiniteElementSpace%s
210  @a trial_fes and @a test_fes. */
211  /** The pointers @a trial_fes and @a test_fes are not owned by the newly
212  constructed object. */
215  : MixedBilinearForm(trial_fes, test_fes),
217  {
220  }
221 
222  /** @brief Create a ParMixedBilinearForm on the given FiniteElementSpace%s
223  @a trial_fes and @a test_fes, using the same integrators as the
224  ParMixedBilinearForm @a mbf.
225 
226  The pointers @a trial_fes and @a test_fes are not owned by the newly
227  constructed object.
228 
229  The integrators in @a mbf are copied as pointers and they are not owned
230  by the newly constructed ParMixedBilinearForm. */
233  ParMixedBilinearForm * mbf)
234  : MixedBilinearForm(trial_fes, test_fes, mbf),
236  {
239  }
240 
241  /// Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
243 
244  /** @brief Returns the matrix assembled on the true dofs, i.e.
245  @a A = P_test^t A_local P_trial, in the format (type id) specified by
246  @a A. */
248 
249  /** @brief Return in @a A a parallel (on truedofs) version of this operator.
250 
251  This returns the same operator as FormRectangularLinearSystem(), but does
252  without the transformations of the right-hand side. */
253  virtual void FormRectangularSystemMatrix(const Array<int> &trial_tdof_list,
254  const Array<int> &test_tdof_list,
255  OperatorHandle &A);
256 
257  /** @brief Form the parallel linear system A X = B, corresponding to this mixed
258  bilinear form and the linear form @a b(.).
259 
260  Return in @a A a *reference* to the system matrix that is column-constrained.
261  The reference will be invalidated when SetOperatorType(), Update(), or the
262  destructor is called. */
263  virtual void FormRectangularLinearSystem(const Array<int> &trial_tdof_list,
264  const Array<int> &test_tdof_list, Vector &x,
265  Vector &b, OperatorHandle &A, Vector &X,
266  Vector &B);
267 
268  /// Compute y += a (P^t A P) x, where x and y are vectors on the true dofs
269  void TrueAddMult(const Vector &x, Vector &y, const double a = 1.0) const;
270 
271  virtual ~ParMixedBilinearForm() { }
272 };
273 
274 /** The parallel matrix representation a linear operator between parallel finite
275  element spaces */
277 {
278 protected:
279  /// Points to the same object as #trial_fes
281  /// Points to the same object as #test_fes
283 
284 private:
285  /// Copy construction is not supported; body is undefined.
287 
288  /// Copy assignment is not supported; body is undefined.
290 
291 public:
292  /** @brief Construct a ParDiscreteLinearOperator on the given
293  FiniteElementSpace%s @a dfes (domain FE space) and @a rfes (range FE
294  space). */
295  /** The pointers @a dfes and @a rfes are not owned by the newly constructed
296  object. */
298  ParFiniteElementSpace *rfes)
299  : DiscreteLinearOperator(dfes, rfes) { domain_fes=dfes; range_fes=rfes; }
300 
301  /// Returns the matrix "assembled" on the true dofs
303 
304  /** Extract the parallel blocks corresponding to the vector dimensions of the
305  domain and range parallel finite element spaces */
306  void GetParBlocks(Array2D<HypreParMatrix *> &blocks) const;
307 
309 };
310 
311 }
312 
313 #endif // MFEM_USE_MPI
314 
315 #endif
ParDiscreteLinearOperator(ParFiniteElementSpace *dfes, ParFiniteElementSpace *rfes)
Construct a ParDiscreteLinearOperator on the given FiniteElementSpaces dfes (domain FE space) and rfe...
void SetOperatorType(Operator::Type tid)
Set the operator type id for the parallel reduced matrix/operator.
Definition: staticcond.hpp:174
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1368
ParFiniteElementSpace * SCParFESpace() const
Return the parallel trace FE space associated with static condensation.
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(...
ParGridFunction X
Auxiliary objects used in TrueAddMult().
void KeepNbrBlock(bool knb=true)
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:895
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
ParFiniteElementSpace * GetParTraceFESpace()
Return a pointer to the parallel reduced/trace FE space.
Definition: staticcond.hpp:113
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(...
void TrueAddMult(const Vector &x, Vector &y, const double a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
OperatorHandle p_mat_e
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void AssembleSharedFaces(int skip_zeros=1)
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
void SetOperatorType(Operator::Type tid)
Set the operator type id for the parallel hybridized matrix/operator.
void SetOperatorType(Operator::Type tid)
Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACYFULL.
ParFiniteElementSpace * pfes
Points to the same object as fes.
void ParallelAssembleElim(OperatorHandle &A_elim)
Data type sparse matrix.
Definition: sparsemat.hpp:46
StaticCondensation * static_cond
Owned.
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
SparseMatrix * mat
Sparse matrix to be associated with the form. Owned.
double b
Definition: lissajous.cpp:42
HypreParMatrix * ParallelEliminateTDofs(const Array< int > &tdofs_list, HypreParMatrix &A) const
Eliminate essential true DOFs from a parallel assembled matrix A.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Return in A a parallel (on truedofs) version of this operator.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
ParFiniteElementSpace * test_pfes
Points to the same object as test_fes.
Dynamic 2D array using row-major layout.
Definition: array.hpp:333
HypreParMatrix * ParallelAssemble() const
Returns the matrix &quot;assembled&quot; on the true dofs.
SparseMatrix * mat_e
Sparse Matrix used to store the eliminations from the b.c. Owned. .
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:70
ParFiniteElementSpace * domain_fes
Points to the same object as trial_fes.
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:237
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
HypreParMatrix * ParallelAssembleElim()
Returns the eliminated matrix assembled on the true dofs, i.e. P^t A_e P.
ParGridFunction X
Auxiliary objects used in TrueAddMult().
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
ParBilinearForm(ParFiniteElementSpace *pf)
Creates parallel bilinear form associated with the FE space *pf.
virtual const Operator * GetProlongation() const
Get the parallel finite element space prolongation matrix.
FiniteElementSpace * test_fes
Not owned.
double a
Definition: lissajous.cpp:41
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Class for parallel bilinear form using different test and trial FE spaces.
ParBilinearForm(ParFiniteElementSpace *pf, ParBilinearForm *bf)
Create a ParBilinearForm on the ParFiniteElementSpace *pf, using the same integrators as the ParBilin...
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:339
ParFiniteElementSpace * range_fes
Points to the same object as test_fes.
void ParallelEliminateEssentialBC(const Array< int > &bdr_attr_is_ess, HypreParMatrix &A, const HypreParVector &X, HypreParVector &B) const
Eliminate essential boundary DOFs from a parallel assembled system.
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 parallel linear system A X = B, corresponding to this mixed bilinear form and the linear for...
ParMixedBilinearForm(ParFiniteElementSpace *trial_fes, ParFiniteElementSpace *test_fes, ParMixedBilinearForm *mbf)
Create a ParMixedBilinearForm on the given FiniteElementSpaces trial_fes and test_fes, using the same integrators as the ParMixedBilinearForm mbf.
Class for parallel bilinear form.
virtual const Operator * GetRestriction() const
Get the parallel finite element space restriction matrix.
ParFiniteElementSpace * trial_pfes
Points to the same object as trial_fes.
FiniteElementSpace * trial_fes
Not owned.
OperatorHandle p_mat
Vector data type.
Definition: vector.hpp:51
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
ID for class HypreParMatrix.
Definition: operator.hpp:241
Hybridization * hybridization
Owned.
ParMixedBilinearForm(ParFiniteElementSpace *trial_fes, ParFiniteElementSpace *test_fes)
Construct a ParMixedBilinearForm on the given FiniteElementSpaces trial_fes and test_fes.
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
void TrueAddMult(const Vector &x, Vector &y, const double a=1.0) const
Compute y += a (P^t A P) x, where x and y are vectors on the true dofs.
OperatorHandle p_mat
Matrix and eliminated matrix.
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:124
void ParallelAssemble(OperatorHandle &A)
Returns the matrix assembled on the true dofs, i.e. A = P^t A_local P, in the format (type id) specif...