12 #ifndef HDIV_LINEAR_SOLVER_HPP 13 #define HDIV_LINEAR_SOLVER_HPP 59 const bool convert_map_type;
65 std::unique_ptr<HypreParMatrix> D, Dt, D_e;
66 std::shared_ptr<DGMassInverse> L_inv;
67 std::shared_ptr<Operator> A_11;
70 Vector L_diag, R_diag, L_diag_unweighted;
75 std::unique_ptr<OperatorJacobiSmoother> R_inv;
76 std::unique_ptr<HypreParMatrix> S;
81 std::unique_ptr<BlockOperator> A_block;
83 std::unique_ptr<BlockDiagonalPreconditioner> D_prec;
88 bool zero_l2_block =
false;
96 mutable Vector b_prime, x_prime, x_bc, w, z;
void EliminateBC(Vector &) const
Eliminates the BCs (called internally, not public interface).
Darcy/mixed Poisson problem.
A coefficient that is constant across space and time.
Pointer to an Operator of a specified type.
Abstract parallel finite element space.
void Setup()
Build the linear operator and solver. Must be called when the coefficients change.
The BoomerAMG solver in hypre.
void SetOperator(const Operator &op) override
No-op.
void Mult(const Vector &b, Vector &x) const override
Solve the linear system for L2 (scalar) and RT (flux) unknowns.
Quadrature function coefficient which requires that the quadrature rules used for this coefficient be...
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Change of basis operator between L2 spaces.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Mode
Which type of saddle-point problem is being solved?
Solve the H(div) saddle-point system using MINRES with matrix-free block-diagonal preconditioning...
Integrated GLL indicator functions.
void SetBC(const Vector &x_rt)
Sets the Dirichlet boundary conditions at the RT essential DOFs.
HdivSaddlePointSolver(ParMesh &mesh_, ParFiniteElementSpace &fes_rt_, ParFiniteElementSpace &fes_l2_, Coefficient &L_coeff_, Coefficient &R_coeff_, const Array< int > &ess_rt_dofs_, Mode mode_)
Creates a solver for the H(div) saddle-point system.
Class representing the storage layout of a QuadratureFunction.
MINRESSolver & GetMINRES()
Returns the internal MINRES solver.
const Array< int > & GetOffsets() const
Return the offsets of the block system.
Class for parallel meshes.
Represents values or vectors of values at quadrature points on a mesh.
int GetNumIterations() const
Get the number of MINRES iterations.
Arbitrary order "L2-conforming" discontinuous finite elements.