1 #include "smith/numerics/block_preconditioner.hpp"
9 #include "axom/slic/core/SimpleLogger.hpp"
10 #include "axom/fmt.hpp"
16 void applyOverrides(
int num_blocks, std::vector<std::unique_ptr<const mfem::Operator>>& block_op_overrides,
17 std::vector<BlockOverride> overrides)
19 for (
auto& ov : overrides) {
20 const int i = ov.first;
23 if (i < 0 || i >= num_blocks) {
24 throw std::out_of_range(
"Override block index out of range");
27 throw std::invalid_argument(
"Override operator must be non-null");
29 if (block_op_overrides[
static_cast<size_t>(i)]) {
30 throw std::invalid_argument(
"Duplicate override for same block index");
33 block_op_overrides[
static_cast<size_t>(i)] = std::move(op);
40 std::vector<BlockOverride> overrides)
42 num_blocks_(static_cast<int>(solvers.
size())),
43 block_jacobian_(nullptr),
44 solver_diag_(nullptr),
45 mfem_solvers_(std::move(solvers)),
46 block_op_overrides_(static_cast<size_t>(num_blocks_))
48 applyOverrides(num_blocks_, block_op_overrides_, std::move(overrides));
55 height = jacobian.Height();
56 width = jacobian.Width();
58 block_jacobian_ =
dynamic_cast<const mfem::BlockOperator*
>(&jacobian);
59 MFEM_VERIFY(block_jacobian_,
"Jacobian must be a BlockOperator");
62 block_jacobian_->NumRowBlocks() != num_blocks_ || block_jacobian_->NumColBlocks() != num_blocks_,
63 axom::fmt::format(
"BlockDiagonalPreconditioner solver count ({}) must match block operator size ({}x{})",
64 num_blocks_, block_jacobian_->NumRowBlocks(), block_jacobian_->NumColBlocks()));
66 block_offsets_.MakeRef(
const_cast<mfem::Array<int>&
>(block_jacobian_->RowOffsets()));
67 solver_diag_ = std::make_unique<mfem::BlockOperator>(block_offsets_);
70 for (
int i = 0; i < num_blocks_; i++) {
72 const mfem::Operator* op =
nullptr;
73 const size_t si =
static_cast<size_t>(i);
75 if (block_op_overrides_[si]) {
76 op = block_op_overrides_[si].get();
78 op = &block_jacobian_->GetBlock(i, i);
81 mfem_solvers_[si]->SetOperator(*op);
84 solver_diag_->SetBlock(i, i, mfem_solvers_[
static_cast<size_t>(i)].get());
88 BlockDiagonalPreconditioner::~BlockDiagonalPreconditioner() {}
92 std::vector<BlockOverride> overrides)
94 num_blocks_(static_cast<int>(solvers.
size())),
95 block_jacobian_(nullptr),
96 mfem_solvers_(std::move(solvers)),
98 block_op_overrides_(static_cast<size_t>(num_blocks_))
100 applyOverrides(num_blocks_, block_op_overrides_, std::move(overrides));
103 void BlockTriangularPreconditioner::LowerSweep(
const mfem::Vector& in, mfem::Vector& out)
const
105 mfem::BlockVector b(
const_cast<mfem::Vector&
>(in), block_offsets_);
106 mfem::BlockVector x(out, block_offsets_);
109 for (
int i = 0; i < num_blocks_; i++) {
110 mfem::Vector& bi = b.GetBlock(i);
111 mfem::Vector& xi = x.GetBlock(i);
114 mfem::Vector rhs_i(bi.Size());
118 for (
int j = 0; j < i; j++) {
119 if (block_jacobian_->IsZeroBlock(i, j)) {
122 const mfem::Operator& A_ij = block_jacobian_->GetBlock(i, j);
124 mfem::Vector tmp(rhs_i.Size());
125 const mfem::Vector& xj = x.GetBlock(j);
128 rhs_i.Add(-1.0, tmp);
132 mfem_solvers_[
static_cast<size_t>(i)]->
Mult(rhs_i, xi);
136 void BlockTriangularPreconditioner::UpperSweep(
const mfem::Vector& in, mfem::Vector& out)
const
138 mfem::BlockVector b(
const_cast<mfem::Vector&
>(in), block_offsets_);
139 mfem::BlockVector x(out, block_offsets_);
142 for (
int i = num_blocks_ - 1; i >= 0; i--) {
143 mfem::Vector& bi = b.GetBlock(i);
144 mfem::Vector& xi = x.GetBlock(i);
147 mfem::Vector rhs_i(bi.Size());
151 for (
int j = i + 1; j < num_blocks_; j++) {
152 if (block_jacobian_->IsZeroBlock(i, j)) {
155 const mfem::Operator& A_ij = block_jacobian_->GetBlock(i, j);
157 mfem::Vector tmp(rhs_i.Size());
158 const mfem::Vector& xj = x.GetBlock(j);
161 rhs_i.Add(-1.0, tmp);
165 mfem_solvers_[
static_cast<size_t>(i)]->
Mult(rhs_i, xi);
185 mfem::Vector tmp(out.Size());
190 mfem::BlockVector tmp_block(tmp, block_offsets_);
192 for (
int i = 0; i < num_blocks_; i++) {
193 mfem::Vector& tmp_i = tmp_block.GetBlock(i);
194 mfem::Vector tmp_i_scaled(tmp_i.Size());
196 const mfem::Operator& A_ii = block_jacobian_->GetBlock(i, i);
197 A_ii.Mult(tmp_i, tmp_i_scaled);
199 tmp_i = tmp_i_scaled;
204 UpperSweep(tmp, out);
212 height = jacobian.Height();
213 width = jacobian.Width();
215 block_jacobian_ =
dynamic_cast<const mfem::BlockOperator*
>(&jacobian);
216 MFEM_VERIFY(block_jacobian_,
"Jacobian must be a BlockOperator");
219 block_jacobian_->NumRowBlocks() != num_blocks_ || block_jacobian_->NumColBlocks() != num_blocks_,
220 axom::fmt::format(
"BlockTriangularPreconditioner solver count ({}) must match block operator size ({}x{})",
221 num_blocks_, block_jacobian_->NumRowBlocks(), block_jacobian_->NumColBlocks()));
223 block_offsets_.MakeRef(
const_cast<mfem::Array<int>&
>(block_jacobian_->RowOffsets()));
226 for (
int i = 0; i < num_blocks_; i++) {
228 const mfem::Operator* op =
nullptr;
229 const size_t si =
static_cast<size_t>(i);
231 if (block_op_overrides_[si]) {
232 op = block_op_overrides_[si].get();
234 op = &block_jacobian_->GetBlock(i, i);
237 mfem_solvers_[si]->SetOperator(*op);
241 BlockTriangularPreconditioner::~BlockTriangularPreconditioner() {}
245 std::vector<BlockOverride> overrides)
247 block_jacobian_(nullptr),
248 solver_diag_(nullptr),
249 mfem_solvers_(std::move(solvers)),
251 approxType_(approxType),
252 block_op_overrides_(static_cast<size_t>(2))
254 SLIC_ERROR_IF(mfem_solvers_.size() != 2,
"This precondition is specifically for 2X2 block systems");
256 applyOverrides(2, block_op_overrides_, std::move(overrides));
259 throw std::invalid_argument(
260 "SchurApproxType::Custom requires an override operator for block index 1 (custom Schur operator)");
264 void BlockSchurPreconditioner::LowerBlock(
const mfem::Vector& in, mfem::Vector& out)
const
267 mfem::BlockVector b(
const_cast<mfem::Vector&
>(in), block_offsets_);
268 mfem::BlockVector x(out, block_offsets_);
270 mfem::Vector& b1 = b.GetBlock(0);
271 mfem::Vector& b2 = b.GetBlock(1);
272 mfem::Vector& x1 = x.GetBlock(0);
273 mfem::Vector& x2 = x.GetBlock(1);
276 mfem_solvers_[0]->Mult(b1, x1);
287 void BlockSchurPreconditioner::UpperBlock(
const mfem::Vector& in, mfem::Vector& out)
const
290 mfem::BlockVector b(
const_cast<mfem::Vector&
>(in), block_offsets_);
291 mfem::BlockVector x(out, block_offsets_);
293 mfem::Vector& b1 = b.GetBlock(0);
294 mfem::Vector& b2 = b.GetBlock(1);
295 mfem::Vector& x1 = x.GetBlock(0);
296 mfem::Vector& x2 = x.GetBlock(1);
299 mfem::Vector rhs1(b1.Size());
300 A_12_->Mult(b2, rhs1);
303 mfem_solvers_[0]->Mult(rhs1, x1);
318 solver_diag_->Mult(in, out);
324 mfem::Vector tmp(out.Size());
326 solver_diag_->Mult(tmp, out);
332 mfem::Vector tmp(out.Size());
333 solver_diag_->Mult(in, tmp);
334 UpperBlock(tmp, out);
340 mfem::Vector tmp(out.Size());
341 mfem::Vector tmp2(out.Size());
343 solver_diag_->Mult(tmp, tmp2);
344 UpperBlock(tmp2, out);
363 mfem::HypreParMatrix* BlockSchurPreconditioner::BuildSchurDiagApprox_(
const mfem::HypreParMatrix& A11,
364 const mfem::HypreParMatrix& A12,
365 const mfem::HypreParMatrix& A21,
366 const mfem::HypreParMatrix& A22)
const
369 auto* Md =
new mfem::HypreParVector(A11.GetComm(), A11.GetGlobalNumRows(), A11.GetRowStarts());
373 auto* A12_scaled =
new mfem::HypreParMatrix(A12);
374 A12_scaled->InvScaleRows(*Md);
380 mfem::HypreParMatrix* A21DinvA12 = mfem::ParMult(&A21, A12_scaled);
382 A12_scaled =
nullptr;
385 mfem::HypreParMatrix* S = mfem::Add(1.0, A22, -1.0, *A21DinvA12);
387 A21DinvA12 =
nullptr;
395 const mfem::HypreParMatrix* H = &K;
396 hypre_ParCSRMatrix* Hhypre =
static_cast<hypre_ParCSRMatrix*
>(*H);
398 hypre_ParCSRMatrixNormFro(Hhypre, &Hfronorm);
404 S_approx_view_ =
nullptr;
405 height = jacobian.Height();
406 width = jacobian.Width();
407 block_jacobian_ =
dynamic_cast<const mfem::BlockOperator*
>(&jacobian);
408 MFEM_VERIFY(block_jacobian_,
"Jacobian must be a BlockOperator");
410 SLIC_ERROR_ROOT_IF(block_jacobian_->NumRowBlocks() != 2 || block_jacobian_->NumColBlocks() != 2,
411 axom::fmt::format(
"BlockSchurPreconditioner requires a 2x2 block operator, got {}x{}",
412 block_jacobian_->NumRowBlocks(), block_jacobian_->NumColBlocks()));
414 mfem_solvers_.size() != 2,
415 axom::fmt::format(
"BlockSchurPreconditioner requires exactly 2 solvers, got {}", mfem_solvers_.size()));
417 block_offsets_.MakeRef(
const_cast<mfem::Array<int>&
>(block_jacobian_->RowOffsets()));
418 solver_diag_ = std::make_unique<mfem::BlockOperator>(block_offsets_);
420 auto* A11 =
dynamic_cast<const mfem::HypreParMatrix*
>(&block_jacobian_->GetBlock(0, 0));
421 auto* A12 =
dynamic_cast<const mfem::HypreParMatrix*
>(&block_jacobian_->GetBlock(0, 1));
422 auto* A21 =
dynamic_cast<const mfem::HypreParMatrix*
>(&block_jacobian_->GetBlock(1, 0));
423 auto* A22 =
dynamic_cast<const mfem::HypreParMatrix*
>(&block_jacobian_->GetBlock(1, 1));
425 MFEM_VERIFY(A11 && A12 && A21 && A22,
426 "All blocks must be HypreParMatrix for assembled Schur complement preconditioner.");
435 const mfem::Operator* op =
nullptr;
436 if (block_op_overrides_[0]) {
437 op = block_op_overrides_[0].get();
441 mfem_solvers_[0]->SetOperator(*op);
444 S_approx_owned_.reset(BuildSchurDiagApprox_(*A11, *A12, *A21, *A22));
445 S_approx_view_ = S_approx_owned_.get();
447 S_approx_owned_.reset(
new mfem::HypreParMatrix(*A22));
448 S_approx_view_ = S_approx_owned_.get();
450 S_approx_owned_.reset();
451 S_approx_view_ = block_op_overrides_[1].get();
454 MFEM_VERIFY(S_approx_view_,
"Schur complement approximation operator must be set");
457 mfem_solvers_[1]->SetOperator(*S_approx_view_);
460 solver_diag_->SetBlock(0, 0, mfem_solvers_[0].
get());
461 solver_diag_->SetBlock(1, 1, mfem_solvers_[1].
get());
464 BlockSchurPreconditioner::~BlockSchurPreconditioner() {}
virtual void Mult(const mfem::Vector &in, mfem::Vector &out) const
The action of the precondition on the block vector (b_1, ..., b_n)
virtual void SetOperator(const mfem::Operator &jacobian)
Set the preconditioner to use the supplied linearized block Jacobian.
BlockDiagonalPreconditioner(std::vector< std::unique_ptr< mfem::Solver >> solvers, std::vector< BlockOverride > overrides={})
Construct a new N by N block diagonal preconditioner.
virtual void Mult(const mfem::Vector &in, mfem::Vector &out) const
The action of the precondition on the block vector (b_1, b_2)
virtual void SetOperator(const mfem::Operator &jacobian)
Set the preconditioner to use the supplied linearized block Jacobian.
BlockSchurPreconditioner(std::vector< std::unique_ptr< mfem::Solver >> solvers, BlockSchurType type=BlockSchurType::Diagonal, SchurApproxType approxType=SchurApproxType::DiagInv, std::vector< BlockOverride > overrides={})
Construct a new 2x2 block Schur complement preconditioner.
virtual void Mult(const mfem::Vector &in, mfem::Vector &out) const
The action of the precondition on the block vector (b_1, ..., b_n)
BlockTriangularPreconditioner(std::vector< std::unique_ptr< mfem::Solver >> solvers, BlockTriangularType type=BlockTriangularType::Lower, std::vector< BlockOverride > overrides={})
Construct a new nxn block triangular preconditioner.
virtual void SetOperator(const mfem::Operator &jacobian)
Set the preconditioner to use the supplied linearized block Jacobian.
Accelerator functionality.
constexpr T & get(variant< T0, T1 > &v)
Returns the variant member of specified type.
SchurApproxType
Selects how the (1,1) Schur operator is approximated.
double matrixNorm(std::unique_ptr< mfem::HypreParMatrix > &K)
Utility to compute the matrix norm.
BlockTriangularType
Selects the block triangular sweep used by BlockTriangularPreconditioner.
BlockSchurType
Selects the block Schur preconditioner variant.
constexpr SMITH_HOST_DEVICE int size(const tensor< T, n... > &)
returns the total number of stored values in a tensor