MFEM
v4.0
Finite element discretization library
|
#include <bilinearform.hpp>
Public Member Functions | |
DiscreteLinearOperator (FiniteElementSpace *domain_fes, FiniteElementSpace *range_fes) | |
Construct a DiscreteLinearOperator on the given FiniteElementSpaces domain_fes and range_fes. More... | |
void | AddDomainInterpolator (DiscreteInterpolator *di) |
Adds a domain interpolator. Assumes ownership of di. More... | |
void | AddTraceFaceInterpolator (DiscreteInterpolator *di) |
Adds a trace face interpolator. Assumes ownership of di. More... | |
Array< BilinearFormIntegrator * > * | GetDI () |
Access all interpolators added with AddDomainInterpolator(). More... | |
virtual void | Assemble (int skip_zeros=1) |
Construct the internal matrix representation of the discrete linear operator. More... | |
Public Member Functions inherited from mfem::MixedBilinearForm | |
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... | |
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... | |
Protected Attributes inherited from mfem::MixedBilinearForm | |
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... | |
Class for constructing the matrix representation of a linear operator, v = L u
, from one FiniteElementSpace (domain) to another FiniteElementSpace (range). The constructed matrix A
is such that
V = A U
where U
and V
are the vectors of degrees of freedom representing the functions u
and v
, respectively. The dimensions of A
are
number of rows of A = dimension of the range space and number of cols of A = dimension of the domain space.
This class is very similar to MixedBilinearForm. One difference is that the linear operator L
is defined using a special kind of BilinearFormIntegrator (we reuse its functionality instead of defining a new class). The other difference with the MixedBilinearForm class is that the "assembly" process overwrites the global matrix entries using the local element matrices instead of adding them.
Note that if we define the bilinear form b(u,v) := (Lu,v)
using an inner product in the range space, then its matrix representation, B
, is
B = M A, (since V^t B U = b(u,v) = (Lu,v) = V^t M A U)
where M
denotes the mass matrix for the inner product in the range space: V1^t M V2 = (v1,v2)
. Similarly, if c(u,w) := (Lu,Lw)
then
C = A^t M A.
Definition at line 660 of file bilinearform.hpp.
|
inline |
Construct a DiscreteLinearOperator on the given FiniteElementSpaces domain_fes and range_fes.
The pointers domain_fes and range_fes are not owned by the newly constructed object.
Definition at line 674 of file bilinearform.hpp.
|
inline |
Adds a domain interpolator. Assumes ownership of di.
Definition at line 679 of file bilinearform.hpp.
|
inline |
Adds a trace face interpolator. Assumes ownership of di.
Definition at line 683 of file bilinearform.hpp.
|
virtual |
Construct the internal matrix representation of the discrete linear operator.
Definition at line 1269 of file bilinearform.cpp.
|
inline |
Access all interpolators added with AddDomainInterpolator().
Definition at line 687 of file bilinearform.hpp.