12#ifndef MFEM_CONSTRAINED
13#define MFEM_CONSTRAINED
23class FiniteElementSpace;
25class ParFiniteElementSpace;
94 Vector& x_and_lambda)
const;
236 { MFEM_ABORT(
"Operator cannot be reset!"); }
265 int dimension_=0,
bool reorder_=
false) :
267 dimension(dimension_), reorder(reorder_)
293 int dimension_=0,
bool reorder_=
false) :
295 dimension(dimension_), reorder(reorder_)
336 { MFEM_ABORT(
"Operator cannot be reset!"); }
552 bool parallel=
false);
A class to handle Block diagonal preconditioners in a matrix-free implementation.
A class to handle Block systems in a matrix-free implementation.
Conjugate gradient method.
An abstract class to solve the constrained system subject to the constraint .
virtual void LagrangeSystemMult(const Vector &f_and_r, Vector &x_and_lambda) const
Solve for (x, lambda) given (f, r)
virtual ~ConstrainedSolver()
void GetMultiplierSolution(Vector &lambda) const
Return the Lagrange multiplier solution in lambda.
virtual void SetConstraintRHS(const Vector &r)
Set the right-hand side r for the constraint B x = r.
ConstrainedSolver(MPI_Comm comm, Operator &A_, Operator &B_)
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
virtual void Mult(const Vector &f, Vector &x) const override
Solve for given .
Data type dense matrix using column-major storage.
virtual Solver * BuildPreconditioner() const override
Build preconditioner for eliminated system.
EliminationCGSolver(HypreParMatrix &A, SparseMatrix &B, Array< int > &constraint_rowstarts, int dimension_=0, bool reorder_=false)
virtual IterativeSolver * BuildKrylov() const override
Select krylov solver for eliminated system.
EliminationGMRESSolver(HypreParMatrix &A, SparseMatrix &B, Array< int > &constraint_rowstarts, int dimension_=0, bool reorder_=false)
virtual IterativeSolver * BuildKrylov() const override
Select krylov solver for eliminated system.
virtual Solver * BuildPreconditioner() const override
Build preconditioner for eliminated system.
SparseMatrix * AssembleExact() const
Assemble this projector as a (processor-local) SparseMatrix.
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void RecoverMultiplier(const Vector &primalrhs, const Vector &primalvars, Vector &lm) const
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
EliminationProjection(const Operator &A, Array< Eliminator * > &eliminators)
void BuildGTilde(const Vector &g, Vector >ilde) const
Solve constrained system by eliminating the constraint; see ConstrainedSolver.
void Mult(const Vector &x, Vector &y) const override
Solve for given .
void SetPreconditioner(Solver &precond) override
This should be called before SetOperator.
Array< Eliminator * > eliminators
void BuildExplicitOperator()
Internal utility routine; assembles eliminated matrix explicitly.
EliminationSolver(HypreParMatrix &A, SparseMatrix &B, Array< int > &primary_dofs, Array< int > &secondary_dofs)
Constructor, with explicit splitting into primary/secondary dofs.
EliminationProjection * projector
virtual Solver * BuildPreconditioner() const =0
Build preconditioner for eliminated system.
virtual IterativeSolver * BuildKrylov() const =0
Select krylov solver for eliminated system.
void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
HypreParMatrix * h_explicit_operator
Perform elimination of a single constraint.
void Eliminate(const Vector &in, Vector &out) const
void ExplicitAssembly(DenseMatrix &mat) const
Return explicitly assembled in mat.
void LagrangeSecondaryTranspose(const Vector &in, Vector &out) const
Transpose of LagrangeSecondary()
void LagrangeSecondary(const Vector &in, Vector &out) const
Maps Lagrange multipliers to secondary dofs, applies .
Eliminator(const SparseMatrix &B, const Array< int > &lagrange_dofs, const Array< int > &primary_tdofs, const Array< int > &secondary_tdofs)
const Array< int > & PrimaryDofs() const
const Array< int > & SecondaryDofs() const
void EliminateTranspose(const Vector &in, Vector &out) const
Transpose of Eliminate(), applies .
const Array< int > & LagrangeDofs() const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
The BoomerAMG solver in hypre.
void SetSystemsOptions(int dim, bool order_bynodes=false)
void SetPrintLevel(int print_level)
Wrapper for hypre's ParCSR matrix class.
Abstract base class for iterative solver.
MPI_Comm GetComm() const
Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
Abstract parallel finite element space.
Solve constrained system with penalty method; see ConstrainedSolver.
HypreParMatrix * penalized_mat
virtual IterativeSolver * BuildKrylov() const =0
Select krylov solver for penalized system.
void Mult(const Vector &x, Vector &y) const override
Solve for given .
void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
virtual Solver * BuildPreconditioner() const =0
Build preconditioner for penalized system.
~PenaltyConstrainedSolver()
void SetPreconditioner(Solver &precond) override
This should be called before SetOperator.
PenaltyConstrainedSolver(HypreParMatrix &A, SparseMatrix &B, real_t penalty_)
PenaltyGMRESSolver(HypreParMatrix &A, HypreParMatrix &B, real_t penalty_, int dimension=0, bool reorder=false)
virtual IterativeSolver * BuildKrylov() const override
Select krylov solver for penalized system.
virtual Solver * BuildPreconditioner() const override
Build preconditioner for penalized system.
PenaltyGMRESSolver(HypreParMatrix &A, HypreParMatrix &B, Vector &penalty_, int dimension=0, bool reorder=false)
PenaltyGMRESSolver(HypreParMatrix &A, SparseMatrix &B, real_t penalty_, int dimension=0, bool reorder=false)
PenaltyPCGSolver(HypreParMatrix &A, HypreParMatrix &B, Vector &penalty_, int dimension=0, bool reorder=false)
PenaltyPCGSolver(HypreParMatrix &A, HypreParMatrix &B, real_t penalty_, int dimension=0, bool reorder=false)
PenaltyPCGSolver(HypreParMatrix &A, SparseMatrix &B, real_t penalty_, int dimension=0, bool reorder=false)
virtual IterativeSolver * BuildKrylov() const override
Select krylov solver for penalized system.
virtual Solver * BuildPreconditioner() const override
Build preconditioner for penalized system.
Basic saddle-point solver with assembled blocks (ie, the operators are assembled HypreParMatrix objec...
SchurConstrainedHypreSolver(MPI_Comm comm, HypreParMatrix &hA_, HypreParMatrix &hB_, Solver *prec=nullptr, int dimension=0, bool reorder=false)
virtual ~SchurConstrainedHypreSolver()
Solve constrained system by solving original mixed system; see ConstrainedSolver.
BlockDiagonalPreconditioner * block_pc
SchurConstrainedSolver(MPI_Comm comm, Operator &A_, Operator &B_, Solver &primal_pc_)
virtual ~SchurConstrainedSolver()
virtual void LagrangeSystemMult(const Vector &x, Vector &y) const override
Solve for (x, lambda) given (f, r)
The transpose of a given operator. Switches the roles of the methods Mult() and MultTranspose().
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
SparseMatrix * BuildNormalConstraints(FiniteElementSpace &fespace, Array< int > &constrained_att, Array< int > &constraint_rowstarts, bool parallel)
Build a matrix constraining normal components to zero.
SparseMatrix * ParBuildNormalConstraints(ParFiniteElementSpace &fespace, Array< int > &constrained_att, Array< int > &constraint_rowstarts)
Parallel wrapper for BuildNormalConstraints.
std::function< real_t(const Vector &)> f(real_t mass_coeff)