|
MFEM
v3.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 | 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 | ParallelEliminateEssentialBC (const Array< int > &bdr_attr_is_ess, HypreParMatrix &A, const HypreParVector &X, HypreParVector &B) const |
| HypreParMatrix * | ParallelEliminateEssentialBC (const Array< int > &bdr_attr_is_ess, HypreParMatrix &A) const |
| HypreParMatrix * | ParallelEliminateTDofs (const Array< int > &tdofs_list, HypreParMatrix &A) const |
| 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... | |
| void | FormLinearSystem (Array< int > &ess_tdof_list, Vector &x, Vector &b, HypreParMatrix &A, Vector &X, Vector &B, int copy_interior=0) |
| 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 | 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 |
| 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 | operator= (const double a) |
| void | Assemble (int skip_zeros=1) |
| Assembles the form i.e. sums over all domain/bdr integrators. More... | |
| void | FormLinearSystem (Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0) |
| void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
| 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 (Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs, int d=0) |
| void | EliminateEssentialBC (Array< int > &bdr_attr_is_ess, int d=0) |
| void | EliminateEssentialBCDiag (Array< int > &bdr_attr_is_ess, double value) |
| Perform elimination and set the diagonal entry to the given value. More... | |
| void | EliminateVDofs (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 (Array< int > &vdofs, int d=0) |
| void | EliminateEssentialBCFromDofs (Array< int > &ess_dofs, Vector &sol, Vector &rhs, int d=0) |
| void | EliminateEssentialBCFromDofs (Array< int > &ess_dofs, int d=0) |
| void | EliminateEssentialBCFromDofsDiag (Array< int > &ess_dofs, double value) |
| Perform elimination and set the diagonal entry to the given value. More... | |
| void | EliminateVDofsInRHS (Array< int > &vdofs, const Vector &x, Vector &b) |
| double | FullInnerProduct (const Vector &x, const Vector &y) const |
| FiniteElementSpace * | GetFES () |
| Return the FE space associated with 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... | |
| virtual void | Print (std::ostream &out=std::cout, 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) | |
| int | Height () const |
| Get the height (size of output) of the Operator. Synonym with NumRows. More... | |
| int | NumRows () const |
| int | Width () const |
| Get the width (size of input) of the Operator. Synonym with NumCols. More... | |
| int | NumCols () const |
| virtual void | MultTranspose (const Vector &x, Vector &y) const |
| Action of the transpose operator. More... | |
| virtual Operator & | GetGradient (const Vector &x) const |
| Evaluate the gradient operator at the point x. 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 () |
Protected Member Functions | |
| void | pAllocMat () |
| void | AssembleSharedFaces (int skip_zeros=1) |
Protected Member Functions inherited from mfem::BilinearForm | |
| void | AllocMat () |
| void | ConformingAssemble () |
| BilinearForm () | |
Class for parallel bilinear form.
Definition at line 29 of file pbilinearform.hpp.
|
inline |
Definition at line 45 of file pbilinearform.hpp.
|
inline |
Definition at line 49 of file pbilinearform.hpp.
|
inlinevirtual |
Definition at line 143 of file pbilinearform.hpp.
| void mfem::ParBilinearForm::Assemble | ( | int | skip_zeros = 1 | ) |
Assemble the local matrix.
Definition at line 208 of file pbilinearform.cpp.
|
protected |
Definition at line 172 of file pbilinearform.cpp.
| void mfem::ParBilinearForm::FormLinearSystem | ( | Array< int > & | ess_tdof_list, |
| Vector & | x, | ||
| Vector & | b, | ||
| HypreParMatrix & | 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 265 of file pbilinearform.cpp.
|
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 57 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.
Definition at line 63 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.
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.
Definition at line 66 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 225 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_new = A_original + 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 238 of file pbilinearform.cpp.
|
inline |
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_new = A_original + 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 97 of file pbilinearform.hpp.
|
inline |
Return the parallel FE space associated with the ParBilinearForm.
Definition at line 105 of file pbilinearform.hpp.
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.
Definition at line 339 of file pbilinearform.cpp.
|
inline |
Return the parallel trace FE space associated with static condensation.
Definition at line 108 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 248 of file pbilinearform.cpp.
|
virtual |
Reimplemented from mfem::BilinearForm.
Definition at line 368 of file pbilinearform.cpp.
|
protected |
Definition at line 37 of file pbilinearform.hpp.
|
protected |
Definition at line 35 of file pbilinearform.hpp.
|
protected |
Definition at line 35 of file pbilinearform.hpp.
|
protected |
Definition at line 32 of file pbilinearform.hpp.
|
mutableprotected |
Definition at line 33 of file pbilinearform.hpp.
|
mutableprotected |
Definition at line 33 of file pbilinearform.hpp.
1.8.5