MFEM
v4.5.2
Finite element discretization library
|
Class for parallel bilinear form. More...
#include <pbilinearform.hpp>
Public Member Functions | |
ParBilinearForm (ParFiniteElementSpace *pf) | |
Creates parallel bilinear form associated with the FE space *pf. More... | |
ParBilinearForm (ParFiniteElementSpace *pf, ParBilinearForm *bf) | |
Create a ParBilinearForm on the ParFiniteElementSpace *pf, using the same integrators as the ParBilinearForm *bf. More... | |
void | KeepNbrBlock (bool knb=true) |
void | SetOperatorType (Operator::Type tid) |
Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACY. More... | |
void | Assemble (int skip_zeros=1) |
Assemble the local matrix. More... | |
virtual void | AssembleDiagonal (Vector &diag) const |
Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector. 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 | ParallelRAP (SparseMatrix &loc_A, OperatorHandle &A, bool steal_loc_A=false) |
Compute parallel RAP operator and store it in A as a HypreParMatrix. 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 * | GetRestrictionTranspose () const |
Get the transpose of GetRestriction, useful for matrix-free RAP. More... | |
virtual const Operator * | GetRestriction () const |
Get the parallel finite element space restriction matrix. More... | |
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(.). More... | |
virtual void | FormSystemMatrix (const Array< int > &ess_tdof_list, OperatorHandle &A) |
Form the linear system matrix A, see FormLinearSystem() for details. More... | |
virtual void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
virtual void | Update (FiniteElementSpace *nfes=NULL) |
Update the FiniteElementSpace and delete all data associated with the old one. More... | |
void | EliminateVDofsInRHS (const Array< int > &vdofs, const Vector &x, Vector &b) |
virtual | ~ParBilinearForm () |
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(.). More... | |
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) |
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.). More... | |
virtual 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) |
Form the linear system matrix A, see FormLinearSystem() for details. More... | |
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) | |
Create a BilinearForm on the FiniteElementSpace f, using the same integrators as the BilinearForm bf. More... | |
int | Size () const |
Get the size of the BilinearForm as a square matrix. More... | |
void | SetAssemblyLevel (AssemblyLevel assembly_level) |
Set the desired assembly level. More... | |
void | EnableSparseMatrixSorting (bool enable_it) |
Force the sparse matrix column indices to be sorted when using AssemblyLevel::FULL. More... | |
AssemblyLevel | GetAssemblyLevel () const |
Returns the assembly level. More... | |
Hybridization * | GetHybridization () const |
void | EnableStaticCondensation () |
Enable the use of static condensation. For details see the description for class StaticCondensation in fem/staticcond.hpp This method should be called before assembly. If the number of unknowns after static condensation is not reduced, it is not enabled. More... | |
bool | StaticCondensationIsEnabled () const |
Check if static condensation was actually enabled by a previous call to EnableStaticCondensation(). More... | |
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) |
Enable hybridization. More... | |
void | UsePrecomputedSparsity (int ps=1) |
For scalar FE spaces, precompute the sparsity pattern of the matrix (assuming dense element matrices) based on the types of integrators present in the bilinear form. More... | |
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 () |
Pre-allocate the internal SparseMatrix before assembly. More... | |
Array< BilinearFormIntegrator * > * | GetDBFI () |
Access all the integrators added with AddDomainIntegrator(). More... | |
Array< BilinearFormIntegrator * > * | GetBBFI () |
Access all the integrators added with AddBoundaryIntegrator(). More... | |
Array< Array< int > * > * | GetBBFI_Marker () |
Access all boundary markers added with AddBoundaryIntegrator(). If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL. More... | |
Array< BilinearFormIntegrator * > * | GetFBFI () |
Access all integrators added with AddInteriorFaceIntegrator(). More... | |
Array< BilinearFormIntegrator * > * | GetBFBFI () |
Access all integrators added with AddBdrFaceIntegrator(). More... | |
Array< Array< int > * > * | GetBFBFI_Marker () |
Access all boundary markers added with AddBdrFaceIntegrator(). If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL. More... | |
const double & | operator() (int i, int j) |
Returns a reference to: \( M_{ij} \). More... | |
virtual double & | Elem (int i, int j) |
Returns a reference to: \( M_{ij} \). More... | |
virtual const double & | Elem (int i, int j) const |
Returns constant reference to: \( M_{ij} \). More... | |
virtual void | Mult (const Vector &x, Vector &y) const |
Matrix vector multiplication: \( y = M x \). More... | |
void | FullMult (const Vector &x, Vector &y) const |
Matrix vector multiplication with the original uneliminated matrix. The original matrix is \( M + M_e \) so we have: \( y = M x + M_e x \). More... | |
virtual void | AddMult (const Vector &x, Vector &y, const double a=1.0) const |
Add the matrix vector multiple to a vector: \( y += a M x \). More... | |
void | FullAddMult (const Vector &x, Vector &y) const |
Add the original uneliminated matrix vector multiple to a vector. The original matrix is \( M + Me \) so we have: \( y += M x + M_e x \). More... | |
virtual void | AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const |
Add the matrix transpose vector multiplication: \( y += a M^T x \). More... | |
void | FullAddMultTranspose (const Vector &x, Vector &y) const |
Add the original uneliminated matrix transpose vector multiple to a vector. The original matrix is \( M + M_e \) so we have: \( y += M^T x + {M_e}^T x \). More... | |
virtual void | MultTranspose (const Vector &x, Vector &y) const |
Matrix transpose vector multiplication: \( y = M^T x \). More... | |
double | InnerProduct (const Vector &x, const Vector &y) const |
Compute \( y^T M x \). More... | |
virtual MatrixInverse * | Inverse () const |
Returns a pointer to (approximation) of the matrix inverse: \( M^{-1} \). More... | |
virtual void | Finalize (int skip_zeros=1) |
Finalizes the matrix initialization. More... | |
const SparseMatrix & | SpMat () const |
Returns a const reference to the sparse matrix: \( M \). More... | |
SparseMatrix & | SpMat () |
Returns a reference to the sparse matrix: \( M \). More... | |
bool | HasSpMat () |
Returns true if the sparse matrix is not null, false otherwise. More... | |
SparseMatrix * | LoseMat () |
Nullifies the internal matrix \( M \) and returns a pointer to it. Used for transferring ownership. More... | |
const SparseMatrix & | SpMatElim () const |
Returns a const reference to the sparse matrix of eliminated b.c.: \( M_e \). More... | |
SparseMatrix & | SpMatElim () |
Returns a reference to the sparse matrix of eliminated b.c.: \( M_e \). More... | |
bool | HasSpMatElim () |
Returns true if the sparse matrix of eliminated b.c.s is not null, false otherwise. More... | |
void | AddDomainIntegrator (BilinearFormIntegrator *bfi) |
Adds new Domain Integrator. Assumes ownership of bfi. More... | |
void | AddDomainIntegrator (BilinearFormIntegrator *bfi, Array< int > &elem_marker) |
void | AddBoundaryIntegrator (BilinearFormIntegrator *bfi) |
Adds new Boundary Integrator. Assumes ownership of bfi. More... | |
void | AddBoundaryIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker) |
Adds new Boundary Integrator, restricted to specific boundary attributes. More... | |
void | AddInteriorFaceIntegrator (BilinearFormIntegrator *bfi) |
Adds new interior Face Integrator. Assumes ownership of bfi. More... | |
void | AddBdrFaceIntegrator (BilinearFormIntegrator *bfi) |
Adds new boundary Face Integrator. Assumes ownership of bfi. 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) |
Sets all sparse values of \( M \) and \( M_e \) to 'a'. More... | |
void | Assemble (int skip_zeros=1) |
Assembles the form i.e. sums over all domain/bdr integrators. More... | |
virtual const Operator * | GetOutputProlongation () const |
Get the output finite element space prolongation matrix. More... | |
virtual const Operator * | GetOutputRestrictionTranspose () const |
Returns the output fe space restriction matrix, transposed. More... | |
virtual const Operator * | GetOutputRestriction () const |
Get the output finite element space restriction matrix. More... | |
void | SerialRAP (OperatorHandle &A) |
Compute serial RAP operator and store it in A as a SparseMatrix. More... | |
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) |
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.). More... | |
template<typename OpType > | |
void | FormSystemMatrix (const Array< int > &ess_tdof_list, OpType &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) |
Compute the element matrix of the given element. More... | |
void | ComputeBdrElementMatrix (int i, DenseMatrix &elmat) |
Compute the boundary element matrix of the given boundary element. More... | |
void | AssembleElementMatrix (int i, const DenseMatrix &elmat, int skip_zeros=1) |
Assemble the given element matrix. More... | |
void | AssembleElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1) |
Assemble the given element matrix. More... | |
void | AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, int skip_zeros=1) |
Assemble the given boundary element matrix. More... | |
void | AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, Array< int > &vdofs, int skip_zeros=1) |
Assemble the given boundary element matrix. More... | |
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. More... | |
void | EliminateEssentialBC (const Array< int > &bdr_attr_is_ess, DiagonalPolicy dpolicy=DIAG_ONE) |
Eliminate essential boundary DOFs from the system matrix. More... | |
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, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE) |
Eliminate the given vdofs. NOTE: here, vdofs is a list of DOFs. More... | |
void | EliminateVDofs (const Array< int > &vdofs, DiagonalPolicy dpolicy=DIAG_ONE) |
Eliminate the given vdofs, storing the eliminated part internally in \( M_e \). More... | |
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_dofs is a marker (boolean) array on all vector-dofs (ess_dofs[i] < 0 is true). More... | |
void | EliminateEssentialBCFromDofs (const Array< int > &ess_dofs, DiagonalPolicy dpolicy=DIAG_ONE) |
Similar to EliminateVDofs(const Array<int> &, DiagonalPolicy) but here ess_dofs is a marker (boolean) array on all vector-dofs (ess_dofs[i] < 0 is true). More... | |
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) |
Use the stored eliminated part of the matrix (see EliminateVDofs(const Array<int> &, DiagonalPolicy)) to modify the r.h.s. b; vdofs is a list of DOFs (non-directional, i.e. >= 0). More... | |
double | FullInnerProduct (const Vector &x, const Vector &y) const |
Compute inner product for full uneliminated matrix \( y^T M x + y^T M_e x \). More... | |
MFEM_DEPRECATED 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... | |
void | SetDiagonalPolicy (DiagonalPolicy policy) |
Sets diagonal policy used upon construction of the linear system. More... | |
void | UseExternalIntegrators () |
Indicate that integrators are not owned by the BilinearForm. 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... | |
bool | IsSquare () const |
Returns whether the matrix is a square matrix. 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 | |
void | InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const |
Initializes memory for true vectors of linear system. More... | |
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 MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. More... | |
virtual void | ArrayMult (const Array< const Vector *> &X, Array< Vector *> &Y) const |
Operator application on a matrix: Y=A(X) . More... | |
virtual void | ArrayMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y) const |
Action of the transpose operator on a matrix: Y=A^t(X) . More... | |
virtual void | ArrayAddMult (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const |
Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X) . More... | |
virtual void | ArrayAddMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const |
Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X) . 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 | FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B) |
Form a column-constrained linear system using a matrix-free approach. More... | |
void | FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this square operator. More... | |
void | FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator (including constraints). More... | |
void | FormDiscreteOperator (Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator. More... | |
void | PrintMatlab (std::ostream &out, int n, int m=0) const |
Prints operator with input size n and output size m in Matlab format. More... | |
virtual void | PrintMatlab (std::ostream &out) const |
Prints operator in Matlab format. More... | |
virtual | ~Operator () |
Virtual destructor. More... | |
Type | GetType () const |
Return the type ID of the Operator class. 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 Member Functions inherited from mfem::Operator | |
void | FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout) |
see FormSystemOperator() More... | |
void | FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout) |
see FormRectangularSystemOperator() More... | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". More... | |
Protected Attributes | |
ParFiniteElementSpace * | pfes |
Points to the same object as fes. More... | |
Vector | Xaux |
Auxiliary vectors used in TrueAddMult(): L-, L-, and T-vector, resp. More... | |
Vector | Yaux |
Vector | Ytmp |
OperatorHandle | p_mat |
OperatorHandle | p_mat_e |
bool | keep_nbr_block |
Protected Attributes inherited from mfem::BilinearForm | |
SparseMatrix * | mat |
Sparse matrix \( M \) to be associated with the form. Owned. More... | |
SparseMatrix * | mat_e |
Sparse Matrix \( M_e \) used to store the eliminations from the b.c. Owned. \( M + M_e = M_{original} \). More... | |
FiniteElementSpace * | fes |
FE space on which the form lives. Not owned. More... | |
AssemblyLevel | assembly |
The assembly level of the form (full, partial, etc.) More... | |
int | batch |
Element batch size used in the form action (1, 8, num_elems, etc.) More... | |
BilinearFormExtension * | ext |
Extension for supporting Full Assembly (FA), Element Assembly (EA), Partial Assembly (PA), or Matrix Free assembly (MF). More... | |
bool | sort_sparse_matrix = false |
long | sequence |
Indicates the Mesh::sequence corresponding to the current state of the BilinearForm. More... | |
int | extern_bfs |
Indicates the BilinearFormIntegrators stored in domain_integs, boundary_integs, interior_face_integs, and boundary_face_integs are owned by another BilinearForm. More... | |
Array< BilinearFormIntegrator * > | domain_integs |
Set of Domain Integrators to be applied. More... | |
Array< Array< int > * > | domain_integs_marker |
Array< BilinearFormIntegrator * > | boundary_integs |
Set of Boundary Integrators to be applied. More... | |
Array< Array< int > * > | boundary_integs_marker |
Entries are not owned. More... | |
Array< BilinearFormIntegrator * > | interior_face_integs |
Set of interior face Integrators to be applied. More... | |
Array< BilinearFormIntegrator * > | boundary_face_integs |
Set of boundary face Integrators to be applied. More... | |
Array< Array< int > * > | boundary_face_integs_marker |
Entries are not owned. More... | |
DenseMatrix | elemmat |
Array< int > | vdofs |
DenseTensor * | element_matrices |
Owned. More... | |
StaticCondensation * | static_cond |
Owned. More... | |
Hybridization * | hybridization |
Owned. More... | |
DiagonalPolicy | diag_policy |
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 | DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP } |
Defines operator diagonal policy upon elimination of rows and/or columns. More... | |
enum | Type { ANY_TYPE, MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ, PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE, PETSC_MATGENERIC, Complex_Operator, MFEM_ComplexSparseMat, Complex_Hypre_ParCSR, Complex_DenseMat, MFEM_Block_Matrix, MFEM_Block_Operator } |
Enumeration defining IDs for some classes derived from Operator. More... | |
Class for parallel bilinear form.
Definition at line 28 of file pbilinearform.hpp.
|
inline |
Creates parallel bilinear form associated with the FE space *pf.
The pointer pf is not owned by the newly constructed object.
Definition at line 56 of file pbilinearform.hpp.
|
inline |
Create a ParBilinearForm on the ParFiniteElementSpace *pf, using the same integrators as the ParBilinearForm *bf.
The pointer pf is not owned by the newly constructed object.
The integrators in bf are copied as pointers and they are not owned by the newly constructed ParBilinearForm.
Definition at line 68 of file pbilinearform.hpp.
|
inlinevirtual |
Definition at line 212 of file pbilinearform.hpp.
void mfem::ParBilinearForm::Assemble | ( | int | skip_zeros = 1 | ) |
Assemble the local matrix.
Definition at line 265 of file pbilinearform.cpp.
|
virtual |
Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector.
When the AssemblyLevel is not LEGACY, and the mesh is nonconforming, this method returns |P^T| d_l, where d_l is the local diagonal of the form before applying parallel/conforming assembly, P^T is the transpose of the parallel/conforming prolongation, and |.| denotes the entry-wise absolute value. In general, this is just an approximation of the exact diagonal for this case.
Reimplemented from mfem::BilinearForm.
Definition at line 284 of file pbilinearform.cpp.
|
protected |
Definition at line 220 of file pbilinearform.cpp.
void mfem::ParBilinearForm::EliminateVDofsInRHS | ( | const Array< int > & | vdofs, |
const Vector & | x, | ||
Vector & | b | ||
) |
Definition at line 421 of file pbilinearform.cpp.
void mfem::BilinearForm::FormLinearSystem |
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
This method applies any necessary transformations to the linear system such as: eliminating boundary conditions; applying conforming constraints for non-conforming AMR; parallel assembly; static condensation; hybridization.
The GridFunction-size vector x must contain the essential b.c. The BilinearForm and the LinearForm-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).
NOTE: If there are no transformations, X simply reuses the data of x.
Definition at line 705 of file bilinearform.cpp.
|
inline |
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
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 506 of file bilinearform.hpp.
|
virtual |
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
This method applies any necessary transformations to the linear system such as: eliminating boundary conditions; applying conforming constraints for non-conforming AMR; parallel assembly; static condensation; hybridization.
The GridFunction-size vector x must contain the essential b.c. The BilinearForm and the LinearForm-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).
NOTE: If there are no transformations, X simply reuses the data of x.
Reimplemented from mfem::BilinearForm.
Definition at line 371 of file pbilinearform.cpp.
|
inline |
Form the linear system matrix A, see FormLinearSystem() for details.
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 528 of file bilinearform.hpp.
void mfem::BilinearForm::FormSystemMatrix |
Form the linear system matrix A, see FormLinearSystem() for details.
Definition at line 774 of file bilinearform.cpp.
|
virtual |
Form the linear system matrix A, see FormLinearSystem() for details.
Reimplemented from mfem::BilinearForm.
Definition at line 427 of file pbilinearform.cpp.
|
inlinevirtual |
Get the parallel finite element space prolongation matrix.
Reimplemented from mfem::BilinearForm.
Definition at line 184 of file pbilinearform.hpp.
|
inlinevirtual |
Get the parallel finite element space restriction matrix.
Reimplemented from mfem::BilinearForm.
Definition at line 190 of file pbilinearform.hpp.
|
inlinevirtual |
Get the transpose of GetRestriction, useful for matrix-free RAP.
Definition at line 187 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 77 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 106 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 212 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 129 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 154 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 110 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 134 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 324 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 337 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 168 of file pbilinearform.hpp.
void mfem::ParBilinearForm::ParallelRAP | ( | SparseMatrix & | loc_A, |
OperatorHandle & | A, | ||
bool | steal_loc_A = false |
||
) |
Compute parallel RAP operator and store it in A as a HypreParMatrix.
[in] | loc_A | The rank-local SparseMatrix . |
[out] | A | The OperatorHandle containing the global HypreParMatrix . |
[in] | steal_loc_A | Have the HypreParMatrix in A take ownership of the memory objects in loc_A. |
Definition at line 124 of file pbilinearform.cpp.
|
inline |
Return the parallel FE space associated with the ParBilinearForm.
Definition at line 177 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 475 of file pbilinearform.cpp.
|
inline |
Return the parallel trace FE space associated with static condensation.
Definition at line 180 of file pbilinearform.hpp.
|
inline |
Set the operator type id for the parallel matrix/operator when using AssemblyLevel::LEGACY.
If using static condensation or hybridization, call this method after enabling it.
Definition at line 83 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 347 of file pbilinearform.cpp.
|
virtual |
Update the FiniteElementSpace and delete all data associated with the old one.
Reimplemented from mfem::BilinearForm.
Definition at line 510 of file pbilinearform.cpp.
|
protected |
Definition at line 39 of file pbilinearform.hpp.
|
protected |
Definition at line 37 of file pbilinearform.hpp.
|
protected |
Definition at line 37 of file pbilinearform.hpp.
|
protected |
Points to the same object as fes.
Definition at line 32 of file pbilinearform.hpp.
|
mutableprotected |
Auxiliary vectors used in TrueAddMult(): L-, L-, and T-vector, resp.
Definition at line 35 of file pbilinearform.hpp.
|
mutableprotected |
Definition at line 35 of file pbilinearform.hpp.
|
mutableprotected |
Definition at line 35 of file pbilinearform.hpp.