MFEM
v4.0
Finite element discretization library
|
#include <bilinearform.hpp>
Public Member Functions | |
MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes) | |
Construct a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces. More... | |
MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes, MixedBilinearForm *mbf) | |
Create a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces, using the same integrators as the MixedBilinearForm mbf. More... | |
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) |
Adds a domain integrator. Assumes ownership of bfi. More... | |
void | AddBoundaryIntegrator (BilinearFormIntegrator *bfi) |
Adds a boundary integrator. Assumes ownership of bfi. More... | |
void | AddTraceFaceIntegrator (BilinearFormIntegrator *bfi) |
Add a trace face integrator. Assumes ownership of bfi. More... | |
Array< BilinearFormIntegrator * > * | GetDBFI () |
Access all integrators added with AddDomainIntegrator(). More... | |
Array< BilinearFormIntegrator * > * | GetBBFI () |
Access all integrators added with AddBoundaryIntegrator(). More... | |
Array< BilinearFormIntegrator * > * | GetTFBFI () |
Access all integrators added with AddTraceFaceIntegrator(). More... | |
void | operator= (const double a) |
void | Assemble (int skip_zeros=1) |
void | ConformingAssemble () |
void | EliminateTrialDofs (Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs) |
void | EliminateEssentialBCFromTrialDofs (Array< int > &marked_vdofs, const 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=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 MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. 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... | |
Type | GetType () const |
Return the type ID of the Operator class. More... | |
Protected Attributes | |
SparseMatrix * | mat |
Owned. More... | |
FiniteElementSpace * | trial_fes |
Not owned. More... | |
FiniteElementSpace * | test_fes |
Not owned. More... | |
int | extern_bfs |
Indicates the BilinearFormIntegrators stored in dom, bdr, and skt are owned by another MixedBilinearForm. More... | |
Array< BilinearFormIntegrator * > | dom |
Domain integrators. More... | |
Array< BilinearFormIntegrator * > | bdr |
Boundary integrators. More... | |
Array< BilinearFormIntegrator * > | skt |
Trace face (skeleton) integrators. More... | |
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::Matrix | |
enum | DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP } |
Public Types inherited from mfem::Operator | |
enum | Type { ANY_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 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 508 of file bilinearform.hpp.
mfem::MixedBilinearForm::MixedBilinearForm | ( | FiniteElementSpace * | tr_fes, |
FiniteElementSpace * | te_fes | ||
) |
Construct a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces.
The pointers tr_fes and te_fes are not owned by the newly constructed object.
Definition at line 998 of file bilinearform.cpp.
mfem::MixedBilinearForm::MixedBilinearForm | ( | FiniteElementSpace * | tr_fes, |
FiniteElementSpace * | te_fes, | ||
MixedBilinearForm * | mbf | ||
) |
Create a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces, using the same integrators as the MixedBilinearForm mbf.
The pointers tr_fes and te_fes are not owned by the newly constructed object.
The integrators in mbf are copied as pointers and they are not owned by the newly constructed MixedBilinearForm.
Definition at line 1008 of file bilinearform.cpp.
|
virtual |
Definition at line 1256 of file bilinearform.cpp.
void mfem::MixedBilinearForm::AddBoundaryIntegrator | ( | BilinearFormIntegrator * | bfi | ) |
Adds a boundary integrator. Assumes ownership of bfi.
Definition at line 1078 of file bilinearform.cpp.
void mfem::MixedBilinearForm::AddDomainIntegrator | ( | BilinearFormIntegrator * | bfi | ) |
Adds a domain integrator. Assumes ownership of bfi.
Definition at line 1073 of file bilinearform.cpp.
|
virtual |
Definition at line 1039 of file bilinearform.cpp.
|
virtual |
Definition at line 1045 of file bilinearform.cpp.
void mfem::MixedBilinearForm::AddTraceFaceIntegrator | ( | BilinearFormIntegrator * | bfi | ) |
Add a trace face integrator. Assumes ownership of bfi.
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 1083 of file bilinearform.cpp.
void mfem::MixedBilinearForm::Assemble | ( | int | skip_zeros = 1 | ) |
Definition at line 1088 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 1173 of file bilinearform.cpp.
|
virtual |
Returns reference to a_{ij}.
Implements mfem::Matrix.
Definition at line 1024 of file bilinearform.cpp.
|
virtual |
Returns constant reference to a_{ij}.
Implements mfem::Matrix.
Definition at line 1029 of file bilinearform.cpp.
void mfem::MixedBilinearForm::EliminateEssentialBCFromTrialDofs | ( | Array< int > & | marked_vdofs, |
const Vector & | sol, | ||
Vector & | rhs | ||
) |
Definition at line 1222 of file bilinearform.cpp.
|
virtual |
Definition at line 1228 of file bilinearform.cpp.
void mfem::MixedBilinearForm::EliminateTrialDofs | ( | Array< int > & | bdr_attr_is_ess, |
const Vector & | sol, | ||
Vector & | rhs | ||
) |
Definition at line 1199 of file bilinearform.cpp.
|
virtual |
Finalizes the matrix initialization.
Reimplemented from mfem::Matrix.
Definition at line 1056 of file bilinearform.cpp.
|
inline |
Access all integrators added with AddBoundaryIntegrator().
Definition at line 600 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 1061 of file bilinearform.cpp.
|
inline |
Access all integrators added with AddDomainIntegrator().
Definition at line 597 of file bilinearform.hpp.
|
inline |
Access all integrators added with AddTraceFaceIntegrator().
Definition at line 603 of file bilinearform.hpp.
|
virtual |
Returns a pointer to (an approximation) of the matrix inverse.
Implements mfem::Matrix.
Definition at line 1051 of file bilinearform.cpp.
|
inline |
Definition at line 581 of file bilinearform.hpp.
Operator application: y=A(x)
.
Implements mfem::Operator.
Definition at line 1034 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 567 of file bilinearform.hpp.
|
inline |
Definition at line 605 of file bilinearform.hpp.
|
inline |
Definition at line 579 of file bilinearform.hpp.
|
inline |
Definition at line 580 of file bilinearform.hpp.
void mfem::MixedBilinearForm::Update | ( | ) |
Definition at line 1248 of file bilinearform.cpp.
|
protected |
Boundary integrators.
Definition at line 523 of file bilinearform.hpp.
|
protected |
Domain integrators.
Definition at line 521 of file bilinearform.hpp.
|
protected |
Indicates the BilinearFormIntegrators stored in dom, bdr, and skt are owned by another MixedBilinearForm.
Definition at line 518 of file bilinearform.hpp.
|
protected |
Owned.
Definition at line 511 of file bilinearform.hpp.
|
protected |
Trace face (skeleton) integrators.
Definition at line 525 of file bilinearform.hpp.
|
protected |
Not owned.
Definition at line 513 of file bilinearform.hpp.
|
protected |
Not owned.
Definition at line 513 of file bilinearform.hpp.