MFEM
v3.3
Finite element discretization library
|
#include <bilinearform.hpp>
Public Member Functions | |
MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes) | |
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 |
Operator application: y=A(x) . More... | |
virtual void | AddMult (const Vector &x, Vector &y, const double a=1.0) const |
virtual void | AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) 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... | |
virtual MatrixInverse * | Inverse () const |
Returns a pointer to (an approximation) of the matrix inverse. More... | |
virtual void | Finalize (int skip_zeros=1) |
Finalizes the matrix initialization. More... | |
void | GetBlocks (Array2D< SparseMatrix * > &blocks) const |
const SparseMatrix & | SpMat () const |
SparseMatrix & | SpMat () |
SparseMatrix * | LoseMat () |
void | AddDomainIntegrator (BilinearFormIntegrator *bfi) |
void | AddBoundaryIntegrator (BilinearFormIntegrator *bfi) |
void | AddTraceFaceIntegrator (BilinearFormIntegrator *bfi) |
Array< BilinearFormIntegrator * > * | GetDBFI () |
Array< BilinearFormIntegrator * > * | GetBBFI () |
Array< BilinearFormIntegrator * > * | GetTFBFI () |
void | operator= (const double a) |
void | Assemble (int skip_zeros=1) |
void | ConformingAssemble () |
void | EliminateTrialDofs (Array< int > &bdr_attr_is_ess, Vector &sol, Vector &rhs) |
void | EliminateEssentialBCFromTrialDofs (Array< int > &marked_vdofs, Vector &sol, Vector &rhs) |
virtual void | EliminateTestDofs (Array< int > &bdr_attr_is_ess) |
void | Update () |
virtual | ~MixedBilinearForm () |
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) | |
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... | |
virtual const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. 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... | |
virtual void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem(). 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 Attributes | |
SparseMatrix * | mat |
FiniteElementSpace * | trial_fes |
FiniteElementSpace * | test_fes |
Array< BilinearFormIntegrator * > | dom |
Array< BilinearFormIntegrator * > | bdr |
Array< BilinearFormIntegrator * > | skt |
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 } |
Enumeration defining IDs for some classes derived from Operator. More... | |
Class for assembling of bilinear forms a(u,v)
defined on different trial and test spaces. The assembled matrix A
is such that
a(u,v) = V^t A U
where U
and V
are the vectors representing the functions u
and v
, respectively. The first argument, u
, of a(,)
is in the trial space and the second argument, v
, is in the test space. Thus,
# of rows of A = dimension of the test space and # of cols of A = dimension of the trial space.
Both trial and test spaces should be defined on the same mesh.
Definition at line 345 of file bilinearform.hpp.
mfem::MixedBilinearForm::MixedBilinearForm | ( | FiniteElementSpace * | tr_fes, |
FiniteElementSpace * | te_fes | ||
) |
Definition at line 898 of file bilinearform.cpp.
|
virtual |
Definition at line 1139 of file bilinearform.cpp.
void mfem::MixedBilinearForm::AddBoundaryIntegrator | ( | BilinearFormIntegrator * | bfi | ) |
Definition at line 961 of file bilinearform.cpp.
void mfem::MixedBilinearForm::AddDomainIntegrator | ( | BilinearFormIntegrator * | bfi | ) |
Definition at line 956 of file bilinearform.cpp.
|
virtual |
Definition at line 922 of file bilinearform.cpp.
|
virtual |
Definition at line 928 of file bilinearform.cpp.
void mfem::MixedBilinearForm::AddTraceFaceIntegrator | ( | BilinearFormIntegrator * | bfi | ) |
Add a trace face integrator. This type of integrator assembles terms over all faces of the mesh using the face FE from the trial space and the two adjacent volume FEs from the test space.
Definition at line 966 of file bilinearform.cpp.
void mfem::MixedBilinearForm::Assemble | ( | int | skip_zeros = 1 | ) |
Definition at line 971 of file bilinearform.cpp.
void mfem::MixedBilinearForm::ConformingAssemble | ( | ) |
For partially conforming trial and/or test FE spaces, complete the assembly process by performing A := P2^t A P1 where A is the internal sparse matrix; P1 and P2 are the conforming prolongation matrices of the trial and test FE spaces, respectively. After this call the MixedBilinearForm becomes an operator on the conforming FE spaces.
Definition at line 1056 of file bilinearform.cpp.
|
virtual |
Returns reference to a_{ij}.
Implements mfem::Matrix.
Definition at line 907 of file bilinearform.cpp.
|
virtual |
Returns constant reference to a_{ij}.
Implements mfem::Matrix.
Definition at line 912 of file bilinearform.cpp.
void mfem::MixedBilinearForm::EliminateEssentialBCFromTrialDofs | ( | Array< int > & | marked_vdofs, |
Vector & | sol, | ||
Vector & | rhs | ||
) |
Definition at line 1105 of file bilinearform.cpp.
|
virtual |
Definition at line 1111 of file bilinearform.cpp.
void mfem::MixedBilinearForm::EliminateTrialDofs | ( | Array< int > & | bdr_attr_is_ess, |
Vector & | sol, | ||
Vector & | rhs | ||
) |
Definition at line 1082 of file bilinearform.cpp.
|
virtual |
Finalizes the matrix initialization.
Reimplemented from mfem::Matrix.
Definition at line 939 of file bilinearform.cpp.
|
inline |
Definition at line 399 of file bilinearform.hpp.
void mfem::MixedBilinearForm::GetBlocks | ( | Array2D< SparseMatrix * > & | blocks | ) | const |
Extract the associated matrix as SparseMatrix blocks. The number of block rows and columns is given by the vector dimensions (vdim) of the test and trial spaces, respectively.
Definition at line 944 of file bilinearform.cpp.
|
inline |
Definition at line 397 of file bilinearform.hpp.
|
inline |
Definition at line 401 of file bilinearform.hpp.
|
virtual |
Returns a pointer to (an approximation) of the matrix inverse.
Implements mfem::Matrix.
Definition at line 934 of file bilinearform.cpp.
|
inline |
Definition at line 386 of file bilinearform.hpp.
Operator application: y=A(x)
.
Implements mfem::Operator.
Definition at line 917 of file bilinearform.cpp.
|
inlinevirtual |
Action of the transpose operator: y=A^t(x)
. The default behavior in class Operator is to generate an error.
Reimplemented from mfem::Operator.
Definition at line 372 of file bilinearform.hpp.
|
inline |
Definition at line 403 of file bilinearform.hpp.
|
inline |
Definition at line 384 of file bilinearform.hpp.
|
inline |
Definition at line 385 of file bilinearform.hpp.
void mfem::MixedBilinearForm::Update | ( | ) |
Definition at line 1131 of file bilinearform.cpp.
|
protected |
Definition at line 353 of file bilinearform.hpp.
|
protected |
Definition at line 352 of file bilinearform.hpp.
|
protected |
Definition at line 348 of file bilinearform.hpp.
|
protected |
Definition at line 354 of file bilinearform.hpp.
|
protected |
Definition at line 350 of file bilinearform.hpp.
|
protected |
Definition at line 350 of file bilinearform.hpp.