12 #ifndef MFEM_BILINEARFORM
13 #define MFEM_BILINEARFORM
15 #include "../config/config.hpp"
16 #include "../linalg/linalg.hpp"
224 virtual double &
Elem(
int i,
int j);
227 virtual const double &
Elem(
int i,
int j)
const;
242 const double a = 1.0)
const
258 virtual void Finalize(
int skip_zeros = 1);
263 MFEM_VERIFY(
mat,
"mat is NULL and can't be dereferenced");
268 MFEM_VERIFY(
mat,
"mat is NULL and can't be dereferenced");
276 MFEM_VERIFY(
mat_e,
"mat_e is NULL and can't be dereferenced");
281 MFEM_VERIFY(
mat_e,
"mat_e is NULL and can't be dereferenced");
315 if (
mat != NULL) { *
mat = a; }
357 Vector &B,
int copy_interior = 0);
368 template <
typename OpType>
371 int copy_interior = 0)
375 OpType *A_ptr = Ah.
Is<OpType>();
376 MFEM_VERIFY(A_ptr,
"invalid OpType used");
392 template <
typename OpType>
397 OpType *A_ptr = Ah.
Is<OpType>();
398 MFEM_VERIFY(A_ptr,
"invalid OpType used");
555 virtual double &
Elem(
int i,
int j);
557 virtual const double &
Elem(
int i,
int j)
const;
562 const double a = 1.0)
const;
565 const double a = 1.0)
const;
572 virtual void Finalize(
int skip_zeros = 1);
691 virtual void Assemble(
int skip_zeros = 1);
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Pointer to an Operator of a specified type.
Data type dense matrix using column-major storage.
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Abstract data type for matrix inverse.
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
OpType * Is() const
Return the Operator pointer dynamically cast to a specified OpType.
DiscreteLinearOperator(FiniteElementSpace *domain_fes, FiniteElementSpace *range_fes)
Construct a DiscreteLinearOperator on the given FiniteElementSpaces domain_fes and range_fes...
Abstract data type matrix.
Set the diagonal value to one.
Dynamic 2D array using row-major layout.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void AddTraceFaceInterpolator(DiscreteInterpolator *di)
Adds a trace face interpolator. Assumes ownership of di.
int height
Dimension of the output / number of rows in the matrix.
FiniteElementSpace * GetTraceFESpace()
Return a pointer to the reduced/trace FE space.
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
Rank 3 tensor (array of matrices)
Array< BilinearFormIntegrator * > * GetDI()
Access all interpolators added with AddDomainInterpolator().
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.