12#ifndef MFEM_DGMASSINV_HPP
13#define MFEM_DGMASSINV_HPP
37 std::unique_ptr<class BilinearForm>
M;
111 template<
int DIM,
int D1D = 0,
int Q1D = 0>
@ GaussLegendre
Open type.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Solver for the discontinuous Galerkin mass matrix.
Vector diag_inv
Jacobi preconditioner.
void MultTranspose(const Vector &b, Vector &u) const override
Same as Mult() since the mass matrix is symmetric.
class MassIntegrator * m
Mass integrator, owned by the form M.
void(DGMassInverse::*)(const Vector &b_, Vector &u) const CGKernelType
real_t rel_tol
Relative CG tolerance.
int max_iter
Maximum number of CG iterations;.
Array< real_t > B_
Inverse of change of basis.
std::unique_ptr< class BilinearForm > M
Mass bilinear form.
real_t abs_tol
Absolute CG tolerance.
DGMassInverse(const FiniteElementSpace &fes_, Coefficient *coeff, const IntegrationRule *ir, int btype)
Protected constructor, used internally.
void SetMaxIter(const int max_iter_)
Set the maximum number of iterations.
Array< real_t > Bt_
Inverse of change of basis, transposed.
DG_FECollection fec
FE collection in requested basis.
void SetOperator(const Operator &op) override
Not implemented. Aborts.
void Update()
Recompute operator and preconditioner (when coefficient or mesh changes).
void DGMassCGIteration(const Vector &b_, Vector &u_) const
Solve the system M b = u. Not part of the public interface.
MFEM_REGISTER_KERNELS(CGKernels, CGKernelType,(int, int, int))
void Mult(const Vector &b, Vector &u) const override
Solve the system M b = u.
void SetRelTol(const real_t rel_tol_)
Set the relative tolerance.
FiniteElementSpace fes
FE space in requested basis.
const DofToQuad * d2q
Change of basis. Not owned.
void SetAbsTol(const real_t abs_tol_)
Set the absolute tolerance.
Structure representing the matrices/tensors needed to evaluate (in reference space) the values,...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Class for an integration rule - an Array of IntegrationPoint.
Arbitrary order "L2-conforming" discontinuous finite elements.
real_t u(const Vector &xvec)