MFEM
v3.3.2
Finite element discretization library
|
Class for parallel bilinear form. More...
#include <pbilinearform.hpp>
Public Member Functions | |
ParBilinearForm (ParFiniteElementSpace *pf) | |
ParBilinearForm (ParFiniteElementSpace *pf, ParBilinearForm *bf) | |
void | KeepNbrBlock (bool knb=true) |
void | SetOperatorType (Operator::Type tid) |
Set the operator type id for the parallel matrix/operator. More... | |
void | Assemble (int skip_zeros=1) |
Assemble the local matrix. More... | |
HypreParMatrix * | ParallelAssemble () |
Returns the matrix assembled on the true dofs, i.e. P^t A P. More... | |
HypreParMatrix * | ParallelAssembleElim () |
Returns the eliminated matrix assembled on the true dofs, i.e. P^t A_e P. More... | |
HypreParMatrix * | ParallelAssemble (SparseMatrix *m) |
Return the matrix m assembled on the true dofs, i.e. P^t A P. More... | |
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) specified by A. More... | |
void | ParallelAssembleElim (OperatorHandle &A_elim) |
void | ParallelAssemble (OperatorHandle &A, SparseMatrix *A_local) |
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. More... | |
HypreParMatrix * | ParallelEliminateEssentialBC (const Array< int > &bdr_attr_is_ess, HypreParMatrix &A) const |
Eliminate essential boundary DOFs from a parallel assembled matrix A. More... | |
HypreParMatrix * | ParallelEliminateTDofs (const Array< int > &tdofs_list, HypreParMatrix &A) const |
Eliminate essential true DOFs from a parallel assembled matrix A. More... | |
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. More... | |
ParFiniteElementSpace * | ParFESpace () const |
Return the parallel FE space associated with the ParBilinearForm. More... | |
ParFiniteElementSpace * | SCParFESpace () const |
Return the parallel trace FE space associated with static condensation. More... | |
virtual const Operator * | GetProlongation () const |
Get the parallel finite element space prolongation matrix. More... | |
virtual const Operator * | GetRestriction () const |
Get the parallel finite element space restriction matrix. More... | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0) |
template<typename OpType > | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, OpType &A, Vector &X, Vector &B, int copy_interior=0) |
void | FormSystemMatrix (const Array< int > &ess_tdof_list, OperatorHandle &A) |
Form the linear system matrix A, see FormLinearSystem() for details. More... | |
template<typename OpType > | |
void | FormSystemMatrix (const Array< int > &ess_tdof_list, OpType &A) |
virtual void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
virtual void | Update (FiniteElementSpace *nfes=NULL) |
virtual | ~ParBilinearForm () |
Public Member Functions inherited from mfem::BilinearForm | |
BilinearForm (FiniteElementSpace *f) | |
Creates bilinear form associated with FE space *f. More... | |
BilinearForm (FiniteElementSpace *f, BilinearForm *bf, int ps=0) | |
int | Size () const |
Get the size of the BilinearForm as a square matrix. More... | |
void | EnableStaticCondensation () |
bool | StaticCondensationIsEnabled () const |
FiniteElementSpace * | SCFESpace () const |
Return the trace FE space associated with static condensation. More... | |
void | EnableHybridization (FiniteElementSpace *constr_space, BilinearFormIntegrator *constr_integ, const Array< int > &ess_tdof_list) |
void | UsePrecomputedSparsity (int ps=1) |
void | UseSparsity (int *I, int *J, bool isSorted) |
Use the given CSR sparsity pattern to allocate the internal SparseMatrix. More... | |
void | UseSparsity (SparseMatrix &A) |
Use the sparsity of A to allocate the internal SparseMatrix. More... | |
void | AllocateMatrix () |
Array< BilinearFormIntegrator * > * | GetDBFI () |
Array< BilinearFormIntegrator * > * | GetBBFI () |
Array< BilinearFormIntegrator * > * | GetFBFI () |
Array< BilinearFormIntegrator * > * | GetBFBFI () |
const double & | operator() (int i, int j) |
virtual double & | Elem (int i, int j) |
Returns reference to a_{ij}. More... | |
virtual const double & | Elem (int i, int j) const |
Returns constant reference to a_{ij}. More... | |
virtual void | Mult (const Vector &x, Vector &y) const |
Matrix vector multiplication. More... | |
void | FullMult (const Vector &x, Vector &y) const |
virtual void | AddMult (const Vector &x, Vector &y, const double a=1.0) const |
void | FullAddMult (const Vector &x, Vector &y) const |
virtual void | AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const |
void | FullAddMultTranspose (const Vector &x, Vector &y) const |
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 error. More... | |
double | InnerProduct (const Vector &x, const Vector &y) const |
virtual MatrixInverse * | Inverse () const |
Returns a pointer to (approximation) of the matrix inverse. More... | |
virtual void | Finalize (int skip_zeros=1) |
Finalizes the matrix initialization. More... | |
const SparseMatrix & | SpMat () const |
Returns a reference to the sparse matrix. More... | |
SparseMatrix & | SpMat () |
SparseMatrix * | LoseMat () |
const SparseMatrix & | SpMatElim () const |
Returns a reference to the sparse matrix of eliminated b.c. More... | |
SparseMatrix & | SpMatElim () |
void | AddDomainIntegrator (BilinearFormIntegrator *bfi) |
Adds new Domain Integrator. More... | |
void | AddBoundaryIntegrator (BilinearFormIntegrator *bfi) |
Adds new Boundary Integrator. More... | |
void | AddInteriorFaceIntegrator (BilinearFormIntegrator *bfi) |
Adds new interior Face Integrator. More... | |
void | AddBdrFaceIntegrator (BilinearFormIntegrator *bfi) |
Adds new boundary Face Integrator. More... | |
void | AddBdrFaceIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker) |
Adds new boundary Face Integrator, restricted to specific boundary attributes. More... | |
void | operator= (const double a) |
void | Assemble (int skip_zeros=1) |
Assembles the form i.e. sums over all domain/bdr integrators. More... | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0) |
void | FormSystemMatrix (const Array< int > &ess_tdof_list, SparseMatrix &A) |
Form the linear system matrix A, see FormLinearSystem for details. More... | |
void | ComputeElementMatrices () |
Compute and store internally all element matrices. More... | |
void | FreeElementMatrices () |
Free the memory used by the element matrices. More... | |
void | ComputeElementMatrix (int i, DenseMatrix &elmat) |
void | AssembleElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1) |
void | AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1) |
void | EliminateEssentialBC (const Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0) |
void | EliminateEssentialBC (const Array< int > &bdr_attr_is_ess, int d=0) |
void | EliminateEssentialBCDiag (const Array< int > &bdr_attr_is_ess, double value) |
Perform elimination and set the diagonal entry to the given value. More... | |
void | EliminateVDofs (const Array< int > &vdofs, Vector &sol, Vector &rhs, int d=0) |
Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs. More... | |
void | EliminateVDofs (const Array< int > &vdofs, int d=0) |
void | EliminateEssentialBCFromDofs (const Array< int > &ess_dofs, Vector &sol, Vector &rhs, int d=0) |
void | EliminateEssentialBCFromDofs (const Array< int > &ess_dofs, int d=0) |
void | EliminateEssentialBCFromDofsDiag (const Array< int > &ess_dofs, double value) |
Perform elimination and set the diagonal entry to the given value. More... | |
void | EliminateVDofsInRHS (const Array< int > &vdofs, const Vector &x, Vector &b) |
double | FullInnerProduct (const Vector &x, const Vector &y) const |
FiniteElementSpace * | GetFES () |
(DEPRECATED) Return the FE space associated with the BilinearForm. More... | |
FiniteElementSpace * | FESpace () |
Return the FE space associated with the BilinearForm. More... | |
const FiniteElementSpace * | FESpace () const |
Read-only access to the associated FiniteElementSpace. More... | |
virtual | ~BilinearForm () |
Destroys bilinear form. More... | |
Public Member Functions inherited from mfem::Matrix | |
Matrix (int s) | |
Creates a square matrix of size s. More... | |
Matrix (int h, int w) | |
Creates a matrix of the given height and width. More... | |
virtual void | Print (std::ostream &out=mfem::out, int width_=4) const |
Prints matrix to stream out. More... | |
virtual | ~Matrix () |
Destroys matrix. More... | |
Public Member Functions inherited from mfem::Operator | |
Operator (int s=0) | |
Construct a square Operator with given size s (default 0). More... | |
Operator (int h, int w) | |
Construct an Operator with the given height (output size) and width (input size). More... | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows(). More... | |
int | NumRows () const |
Get the number of rows (size of output) of the Operator. Synonym with Height(). More... | |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols(). More... | |
int | NumCols () const |
Get the number of columns (size of input) of the Operator. Synonym with Width(). More... | |
virtual Operator & | GetGradient (const Vector &x) const |
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More... | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0) |
Form a constrained linear system using a matrix-free approach. More... | |
void | PrintMatlab (std::ostream &out, int n=0, int m=0) const |
Prints operator with input size n and output size m in Matlab format. More... | |
virtual | ~Operator () |
Virtual destructor. More... | |
Protected Member Functions | |
void | pAllocMat () |
void | AssembleSharedFaces (int skip_zeros=1) |
Protected Member Functions inherited from mfem::BilinearForm | |
void | AllocMat () |
void | ConformingAssemble () |
BilinearForm () | |
Protected Attributes | |
ParFiniteElementSpace * | pfes |
ParGridFunction | X |
ParGridFunction | Y |
OperatorHandle | p_mat |
OperatorHandle | p_mat_e |
bool | keep_nbr_block |
Protected Attributes inherited from mfem::BilinearForm | |
SparseMatrix * | mat |
Sparse matrix to be associated with the form. More... | |
SparseMatrix * | mat_e |
Matrix used to eliminate b.c. More... | |
FiniteElementSpace * | fes |
FE space on which the form lives. More... | |
long | sequence |
int | extern_bfs |
Array< BilinearFormIntegrator * > | dbfi |
Set of Domain Integrators to be applied. More... | |
Array< BilinearFormIntegrator * > | bbfi |
Set of Boundary Integrators to be applied. More... | |
Array< BilinearFormIntegrator * > | fbfi |
Set of interior face Integrators to be applied. More... | |
Array< BilinearFormIntegrator * > | bfbfi |
Set of boundary face Integrators to be applied. More... | |
Array< Array< int > * > | bfbfi_marker |
DenseMatrix | elemmat |
Array< int > | vdofs |
DenseTensor * | element_matrices |
StaticCondensation * | static_cond |
Hybridization * | hybridization |
int | precompute_sparsity |
Protected Attributes inherited from mfem::Operator | |
int | height |
Dimension of the output / number of rows in the matrix. More... | |
int | width |
Dimension of the input / number of columns in the matrix. More... | |
Additional Inherited Members | |
Public Types inherited from mfem::Operator | |
enum | Type { MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ, PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE, PETSC_MATGENERIC } |
Enumeration defining IDs for some classes derived from Operator. More... | |
Class for parallel bilinear form.
Definition at line 28 of file pbilinearform.hpp.
|
inline |
Definition at line 44 of file pbilinearform.hpp.
|
inline |
Definition at line 49 of file pbilinearform.hpp.
|
inlinevirtual |
Definition at line 215 of file pbilinearform.hpp.
void mfem::ParBilinearForm::Assemble | ( | int | skip_zeros = 1 | ) |
Assemble the local matrix.
Definition at line 222 of file pbilinearform.cpp.
|
protected |
Definition at line 186 of file pbilinearform.cpp.
void mfem::ParBilinearForm::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 the current bilinear form and b(.), by applying any necessary transformations such as: eliminating boundary conditions; applying conforming constraints for non-conforming AMR; parallel assembly; static condensation; hybridization.
The ParGridFunction-size vector x must contain the essential b.c. The ParBilinearForm and the ParLinearForm-size vector b must be assembled.
The vector X is initialized with a suitable initial guess: when using hybridization, the vector X is set to zero; otherwise, the essential entries of X are set to the corresponding b.c. and all other entries are set to zero (copy_interior == 0) or copied from x (copy_interior != 0).
This method can be called multiple times (with the same ess_tdof_list array) to initialize different right-hand sides and boundary condition values.
After solving the linear system, the finite element solution x can be recovered by calling RecoverFEMSolution (with the same vectors X, b, and x).
Definition at line 279 of file pbilinearform.cpp.
|
inline |
Version of the method FormLinearSystem() where the system matrix is returned in the variable A, of type OpType, holding a reference to the system matrix (created with the method OpType::MakeRef()). The reference will be invalidated when SetOperatorType(), Update(), or the destructor is called.
Definition at line 179 of file pbilinearform.hpp.
void mfem::ParBilinearForm::FormSystemMatrix | ( | const Array< int > & | ess_tdof_list, |
OperatorHandle & | A | ||
) |
Form the linear system matrix A, see FormLinearSystem() for details.
Definition at line 323 of file pbilinearform.cpp.
|
inline |
Version of the method FormSystemMatrix() where the system matrix is returned in the variable A, of type OpType, holding a reference to the system matrix (created with the method OpType::MakeRef()). The reference will be invalidated when SetOperatorType(), Update(), or the destructor is called.
Definition at line 199 of file pbilinearform.hpp.
|
inlinevirtual |
Get the parallel finite element space prolongation matrix.
Reimplemented from mfem::BilinearForm.
Definition at line 142 of file pbilinearform.hpp.
|
inlinevirtual |
Get the parallel finite element space restriction matrix.
Reimplemented from mfem::BilinearForm.
Definition at line 145 of file pbilinearform.hpp.
|
inline |
When set to true and the ParBilinearForm has interior face integrators, the local SparseMatrix will include the rows (in addition to the columns) corresponding to face-neighbor dofs. The default behavior is to disregard those rows. Must be called before the first Assemble call.
Definition at line 58 of file pbilinearform.hpp.
|
protected |
Definition at line 22 of file pbilinearform.cpp.
|
inline |
Returns the matrix assembled on the true dofs, i.e. P^t A P.
The returned matrix has to be deleted by the caller.
Definition at line 75 of file pbilinearform.hpp.
HypreParMatrix * mfem::ParBilinearForm::ParallelAssemble | ( | SparseMatrix * | m | ) |
Return the matrix m assembled on the true dofs, i.e. P^t A P.
The returned matrix has to be deleted by the caller.
Definition at line 178 of file pbilinearform.cpp.
|
inline |
Returns the matrix assembled on the true dofs, i.e. A = P^t A_local P, in the format (type id) specified by A.
Definition at line 87 of file pbilinearform.hpp.
void mfem::ParBilinearForm::ParallelAssemble | ( | OperatorHandle & | A, |
SparseMatrix * | A_local | ||
) |
Returns the matrix A_local assembled on the true dofs, i.e. A = P^t A_local P in the format (type id) specified by A.
Definition at line 124 of file pbilinearform.cpp.
|
inline |
Returns the eliminated matrix assembled on the true dofs, i.e. P^t A_e P.
The returned matrix has to be deleted by the caller.
Definition at line 79 of file pbilinearform.hpp.
|
inline |
Returns the eliminated matrix assembled on the true dofs, i.e. A_elim = P^t A_elim_local P in the format (type id) specified by A.
Definition at line 92 of file pbilinearform.hpp.
void mfem::ParBilinearForm::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.
The array bdr_attr_is_ess marks boundary attributes that constitute the essential part of the boundary.
Definition at line 239 of file pbilinearform.cpp.
HypreParMatrix * mfem::ParBilinearForm::ParallelEliminateEssentialBC | ( | const Array< int > & | bdr_attr_is_ess, |
HypreParMatrix & | A | ||
) | const |
Eliminate essential boundary DOFs from a parallel assembled matrix A.
The array bdr_attr_is_ess marks boundary attributes that constitute the essential part of the boundary. The eliminated part is stored in a matrix A_elim such that A_original = A_new + A_elim. Returns a pointer to the newly allocated matrix A_elim which should be deleted by the caller. The matrices A and A_elim can be used to eliminate boundary conditions in multiple right-hand sides, by calling the function EliminateBC() (from hypre.hpp).
Definition at line 252 of file pbilinearform.cpp.
|
inline |
Eliminate essential true DOFs from a parallel assembled matrix A.
Given a list of essential true dofs and the parallel assembled matrix A, eliminate the true dofs from the matrix, storing the eliminated part in a matrix A_elim such that A_original = A_new + A_elim. Returns a pointer to the newly allocated matrix A_elim which should be deleted by the caller. The matrices A and A_elim can be used to eliminate boundary conditions in multiple right-hand sides, by calling the function EliminateBC() (from hypre.hpp).
Definition at line 126 of file pbilinearform.hpp.
|
inline |
Return the parallel FE space associated with the ParBilinearForm.
Definition at line 135 of file pbilinearform.hpp.
|
virtual |
Call this method after solving a linear system constructed using the FormLinearSystem method to recover the solution as a ParGridFunction-size vector in x. Use the same arguments as in the FormLinearSystem call.
Reimplemented from mfem::BilinearForm.
Definition at line 365 of file pbilinearform.cpp.
|
inline |
Return the parallel trace FE space associated with static condensation.
Definition at line 138 of file pbilinearform.hpp.
|
inline |
Set the operator type id for the parallel matrix/operator.
If using static condensation or hybridization, call this method after enabling it.
Definition at line 63 of file pbilinearform.hpp.
void mfem::ParBilinearForm::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.
Definition at line 262 of file pbilinearform.cpp.
|
virtual |
Reimplemented from mfem::BilinearForm.
Definition at line 394 of file pbilinearform.cpp.
|
protected |
Definition at line 36 of file pbilinearform.hpp.
|
protected |
Definition at line 34 of file pbilinearform.hpp.
|
protected |
Definition at line 34 of file pbilinearform.hpp.
|
protected |
Definition at line 31 of file pbilinearform.hpp.
|
mutableprotected |
Definition at line 32 of file pbilinearform.hpp.
|
mutableprotected |
Definition at line 32 of file pbilinearform.hpp.