Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
block_preconditioner.cpp
1 #include "smith/numerics/block_preconditioner.hpp"
2 
3 #include <memory>
4 #include <utility>
5 #include <vector>
6 #include <stdexcept>
7 
8 #include "mfem.hpp"
9 #include "axom/slic/core/SimpleLogger.hpp"
10 #include "axom/fmt.hpp"
11 
12 namespace smith {
13 
14 namespace {
15 
16 void applyOverrides(int num_blocks, std::vector<std::unique_ptr<const mfem::Operator>>& block_op_overrides,
17  std::vector<BlockOverride> overrides)
18 {
19  for (auto& ov : overrides) {
20  const int i = ov.first;
21  auto& op = ov.second;
22 
23  if (i < 0 || i >= num_blocks) {
24  throw std::out_of_range("Override block index out of range");
25  }
26  if (!op) {
27  throw std::invalid_argument("Override operator must be non-null");
28  }
29  if (block_op_overrides[static_cast<size_t>(i)]) {
30  throw std::invalid_argument("Duplicate override for same block index");
31  }
32 
33  block_op_overrides[static_cast<size_t>(i)] = std::move(op);
34  }
35 }
36 
37 } // namespace
38 
39 BlockDiagonalPreconditioner::BlockDiagonalPreconditioner(std::vector<std::unique_ptr<mfem::Solver>> solvers,
40  std::vector<BlockOverride> overrides)
41  : block_offsets_(),
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_))
47 {
48  applyOverrides(num_blocks_, block_op_overrides_, std::move(overrides));
49 }
50 
51 void BlockDiagonalPreconditioner::Mult(const mfem::Vector& in, mfem::Vector& out) const { solver_diag_->Mult(in, out); }
52 
53 void BlockDiagonalPreconditioner::SetOperator(const mfem::Operator& jacobian)
54 {
55  height = jacobian.Height();
56  width = jacobian.Width();
57  // Cast the supplied jacobian to a block operator object
58  block_jacobian_ = dynamic_cast<const mfem::BlockOperator*>(&jacobian);
59  MFEM_VERIFY(block_jacobian_, "Jacobian must be a BlockOperator");
60 
61  SLIC_ERROR_ROOT_IF(
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()));
65 
66  block_offsets_.MakeRef(const_cast<mfem::Array<int>&>(block_jacobian_->RowOffsets()));
67  solver_diag_ = std::make_unique<mfem::BlockOperator>(block_offsets_);
68 
69  // For each diagonal block A_ii, configure the corresponding solver
70  for (int i = 0; i < num_blocks_; i++) {
71  // Attach operator to solver
72  const mfem::Operator* op = nullptr;
73  const size_t si = static_cast<size_t>(i);
74 
75  if (block_op_overrides_[si]) {
76  op = block_op_overrides_[si].get(); // use override
77  } else {
78  op = &block_jacobian_->GetBlock(i, i); // use Jacobian diagonal block
79  }
80 
81  mfem_solvers_[si]->SetOperator(*op);
82 
83  // Place the solver into the diagonal block of solver_diag_
84  solver_diag_->SetBlock(i, i, mfem_solvers_[static_cast<size_t>(i)].get());
85  }
86 }
87 
88 BlockDiagonalPreconditioner::~BlockDiagonalPreconditioner() {}
89 
90 BlockTriangularPreconditioner::BlockTriangularPreconditioner(std::vector<std::unique_ptr<mfem::Solver>> solvers,
92  std::vector<BlockOverride> overrides)
93  : block_offsets_(),
94  num_blocks_(static_cast<int>(solvers.size())),
95  block_jacobian_(nullptr),
96  mfem_solvers_(std::move(solvers)),
97  type_(type),
98  block_op_overrides_(static_cast<size_t>(num_blocks_))
99 {
100  applyOverrides(num_blocks_, block_op_overrides_, std::move(overrides));
101 }
102 
103 void BlockTriangularPreconditioner::LowerSweep(const mfem::Vector& in, mfem::Vector& out) const
104 {
105  mfem::BlockVector b(const_cast<mfem::Vector&>(in), block_offsets_);
106  mfem::BlockVector x(out, block_offsets_);
107 
108  // Forward sweep: i = 0 .. num_blocks_ - 1
109  for (int i = 0; i < num_blocks_; i++) {
110  mfem::Vector& bi = b.GetBlock(i);
111  mfem::Vector& xi = x.GetBlock(i);
112 
113  // rhs_i = b_i
114  mfem::Vector rhs_i(bi.Size());
115  rhs_i = bi;
116 
117  // Subtract sum_{j < i} A_ij x_j
118  for (int j = 0; j < i; j++) {
119  if (block_jacobian_->IsZeroBlock(i, j)) {
120  continue; // no coupling
121  }
122  const mfem::Operator& A_ij = block_jacobian_->GetBlock(i, j);
123 
124  mfem::Vector tmp(rhs_i.Size());
125  const mfem::Vector& xj = x.GetBlock(j);
126 
127  A_ij.Mult(xj, tmp); // tmp = A_ij x_j
128  rhs_i.Add(-1.0, tmp); // rhs_i -= A_ij x_j
129  }
130 
131  // Solve A_ii x_i = rhs_i with the i-th block solver
132  mfem_solvers_[static_cast<size_t>(i)]->Mult(rhs_i, xi);
133  }
134 }
135 
136 void BlockTriangularPreconditioner::UpperSweep(const mfem::Vector& in, mfem::Vector& out) const
137 {
138  mfem::BlockVector b(const_cast<mfem::Vector&>(in), block_offsets_);
139  mfem::BlockVector x(out, block_offsets_);
140 
141  // Backward sweep: i = num_blocks_ - 1 .. 0
142  for (int i = num_blocks_ - 1; i >= 0; i--) {
143  mfem::Vector& bi = b.GetBlock(i);
144  mfem::Vector& xi = x.GetBlock(i);
145 
146  // rhs_i = b_i
147  mfem::Vector rhs_i(bi.Size());
148  rhs_i = bi;
149 
150  // Subtract sum_{j > i} A_ij x_j
151  for (int j = i + 1; j < num_blocks_; j++) {
152  if (block_jacobian_->IsZeroBlock(i, j)) {
153  continue; // no coupling
154  }
155  const mfem::Operator& A_ij = block_jacobian_->GetBlock(i, j);
156 
157  mfem::Vector tmp(rhs_i.Size());
158  const mfem::Vector& xj = x.GetBlock(j);
159 
160  A_ij.Mult(xj, tmp); // tmp = A_ij x_j
161  rhs_i.Add(-1.0, tmp); // rhs_i -= A_ij x_j
162  }
163 
164  // Solve A_ii x_i = rhs_i
165  mfem_solvers_[static_cast<size_t>(i)]->Mult(rhs_i, xi);
166  }
167 }
168 
169 void BlockTriangularPreconditioner::Mult(const mfem::Vector& in, mfem::Vector& out) const
170 {
171  switch (type_) {
173  // x = P_lower^{-1} b
174  LowerSweep(in, out);
175  break;
176 
178  // x = P_upper^{-1} b
179  UpperSweep(in, out);
180  break;
181 
183  // Symmetric: x = P_upper^{-1} D P_lower^{-1} b
184  // 1) tmp = P_lower^{-1} b
185  mfem::Vector tmp(out.Size());
186  LowerSweep(in, tmp);
187 
188  // 2) tmp = D * tmp where D = diag(A_ii)
189  {
190  mfem::BlockVector tmp_block(tmp, block_offsets_);
191 
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());
195 
196  const mfem::Operator& A_ii = block_jacobian_->GetBlock(i, i);
197  A_ii.Mult(tmp_i, tmp_i_scaled); // tmp_i_scaled = A_ii * tmp_i
198 
199  tmp_i = tmp_i_scaled; // write back into block vector
200  }
201  }
202 
203  // 3) out = P_upper^{-1} tmp
204  UpperSweep(tmp, out);
205  break;
206  }
207  }
208 }
209 
210 void BlockTriangularPreconditioner::SetOperator(const mfem::Operator& jacobian)
211 {
212  height = jacobian.Height();
213  width = jacobian.Width();
214  // Cast the supplied jacobian to a block operator object
215  block_jacobian_ = dynamic_cast<const mfem::BlockOperator*>(&jacobian);
216  MFEM_VERIFY(block_jacobian_, "Jacobian must be a BlockOperator");
217 
218  SLIC_ERROR_ROOT_IF(
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()));
222 
223  block_offsets_.MakeRef(const_cast<mfem::Array<int>&>(block_jacobian_->RowOffsets()));
224 
225  // Configure all diagonal solves
226  for (int i = 0; i < num_blocks_; i++) {
227  // Attach operator to solver
228  const mfem::Operator* op = nullptr;
229  const size_t si = static_cast<size_t>(i);
230 
231  if (block_op_overrides_[si]) {
232  op = block_op_overrides_[si].get(); // use override
233  } else {
234  op = &block_jacobian_->GetBlock(i, i); // use Jacobian diagonal block
235  }
236 
237  mfem_solvers_[si]->SetOperator(*op);
238  }
239 }
240 
241 BlockTriangularPreconditioner::~BlockTriangularPreconditioner() {}
242 
243 BlockSchurPreconditioner::BlockSchurPreconditioner(std::vector<std::unique_ptr<mfem::Solver>> solvers,
244  BlockSchurType type, SchurApproxType approxType,
245  std::vector<BlockOverride> overrides)
246  : block_offsets_(),
247  block_jacobian_(nullptr),
248  solver_diag_(nullptr),
249  mfem_solvers_(std::move(solvers)),
250  type_(type),
251  approxType_(approxType),
252  block_op_overrides_(static_cast<size_t>(2))
253 {
254  SLIC_ERROR_IF(mfem_solvers_.size() != 2, "This precondition is specifically for 2X2 block systems");
255 
256  applyOverrides(2, block_op_overrides_, std::move(overrides));
257 
258  if (approxType_ == SchurApproxType::Custom && !block_op_overrides_[1]) {
259  throw std::invalid_argument(
260  "SchurApproxType::Custom requires an override operator for block index 1 (custom Schur operator)");
261  }
262 }
263 
264 void BlockSchurPreconditioner::LowerBlock(const mfem::Vector& in, mfem::Vector& out) const
265 {
266  // Interpret in, out as block vectors: in = [b1; b2], out = [x1; x2]
267  mfem::BlockVector b(const_cast<mfem::Vector&>(in), block_offsets_);
268  mfem::BlockVector x(out, block_offsets_);
269 
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);
274 
275  // 1) Solve A11 x1 = b1
276  mfem_solvers_[0]->Mult(b1, x1);
277 
278  // 2) Build x2 = b2 - A21 x1
279  A_21_->Mult(x1, x2); // x2 = A21 x1
280  x2.Neg(); // x2 = -A21 x1
281  x2 += b2; // x2 = b2 - A21 x1
282 
283  // 3) Reassign x1.
284  x1 = b1;
285 }
286 
287 void BlockSchurPreconditioner::UpperBlock(const mfem::Vector& in, mfem::Vector& out) const
288 {
289  // Interpret in, out as block vectors: in = [b1; b2], out = [x1; x2]
290  mfem::BlockVector b(const_cast<mfem::Vector&>(in), block_offsets_);
291  mfem::BlockVector x(out, block_offsets_);
292 
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);
297 
298  // 1) Build x1 = A12 b2
299  mfem::Vector rhs1(b1.Size());
300  A_12_->Mult(b2, rhs1); // rhs1 = A12 b2
301 
302  // 2) Solve A11 x1 = rhs1
303  mfem_solvers_[0]->Mult(rhs1, x1);
304 
305  // 3) Build b1 - A11^-1 A12 b2
306  x1.Neg(); // x1 = -x1
307  x1 += b1; // = b1 - A12 x2
308 
309  // 4) Assign x2
310  x2 = b2;
311 }
312 
313 void BlockSchurPreconditioner::Mult(const mfem::Vector& in, mfem::Vector& out) const
314 {
315  switch (type_) {
317  // x = [A11^-1, 0; 0, S^-1] b
318  solver_diag_->Mult(in, out);
319  break;
320  }
321 
322  case BlockSchurType::Lower: {
323  // x = [A11^-1, 0; 0, S^-1][I, 0; -A21 A11^-1, I] b
324  mfem::Vector tmp(out.Size());
325  LowerBlock(in, tmp);
326  solver_diag_->Mult(tmp, out);
327  break;
328  }
329 
330  case BlockSchurType::Upper: {
331  // x = [I, -A11^-1 A12; 0, I][A11^-1, 0; 0, S^-1] b
332  mfem::Vector tmp(out.Size());
333  solver_diag_->Mult(in, tmp);
334  UpperBlock(tmp, out);
335  break;
336  }
337 
338  case BlockSchurType::Full: {
339  // x = [I, -A11^-1 A12; 0, I][A11^-1, 0; 0, S^-1][I, 0; -A21 A11^-1, I] b
340  mfem::Vector tmp(out.Size());
341  mfem::Vector tmp2(out.Size());
342  LowerBlock(in, tmp);
343  solver_diag_->Mult(tmp, tmp2);
344  UpperBlock(tmp2, out);
345  break;
346  }
347  }
348 }
349 
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
367 {
368  // Extract diagonal of A11
369  auto* Md = new mfem::HypreParVector(A11.GetComm(), A11.GetGlobalNumRows(), A11.GetRowStarts());
370  A11.GetDiag(*Md);
371 
372  // Scale rows of A12 by diag(A11)^{-1}
373  auto* A12_scaled = new mfem::HypreParMatrix(A12);
374  A12_scaled->InvScaleRows(*Md);
375 
376  delete Md;
377  Md = nullptr;
378 
379  // Compute A21 * (diag(A11)^{-1} * A12)
380  mfem::HypreParMatrix* A21DinvA12 = mfem::ParMult(&A21, A12_scaled);
381  delete A12_scaled;
382  A12_scaled = nullptr;
383 
384  // S_approx = A22 - A21 * diag(A11)^{-1} * A12
385  mfem::HypreParMatrix* S = mfem::Add(1.0, A22, -1.0, *A21DinvA12);
386  delete A21DinvA12;
387  A21DinvA12 = nullptr;
388 
389  return S; // caller owns
390 }
391 
393 double matrixNorm(const mfem::HypreParMatrix& K)
394 {
395  const mfem::HypreParMatrix* H = &K;
396  hypre_ParCSRMatrix* Hhypre = static_cast<hypre_ParCSRMatrix*>(*H);
397  double Hfronorm;
398  hypre_ParCSRMatrixNormFro(Hhypre, &Hfronorm);
399  return Hfronorm;
400 }
401 
402 void BlockSchurPreconditioner::SetOperator(const mfem::Operator& jacobian)
403 {
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");
409 
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()));
413  SLIC_ERROR_ROOT_IF(
414  mfem_solvers_.size() != 2,
415  axom::fmt::format("BlockSchurPreconditioner requires exactly 2 solvers, got {}", mfem_solvers_.size()));
416 
417  block_offsets_.MakeRef(const_cast<mfem::Array<int>&>(block_jacobian_->RowOffsets()));
418  solver_diag_ = std::make_unique<mfem::BlockOperator>(block_offsets_);
419 
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));
424 
425  MFEM_VERIFY(A11 && A12 && A21 && A22,
426  "All blocks must be HypreParMatrix for assembled Schur complement preconditioner.");
427 
428  if (type_ == BlockSchurType::Lower || type_ == BlockSchurType::Full) {
429  A_21_ = A21;
430  }
431  if (type_ == BlockSchurType::Upper || type_ == BlockSchurType::Full) {
432  A_12_ = A12;
433  }
434  // Diagonal preconditioner for block (0,0)
435  const mfem::Operator* op = nullptr;
436  if (block_op_overrides_[0]) {
437  op = block_op_overrides_[0].get(); // use override
438  } else {
439  op = A11; // use Jacobian diagonal block
440  }
441  mfem_solvers_[0]->SetOperator(*op);
442  // Build Schur complement approximation
443  if (approxType_ == SchurApproxType::DiagInv) {
444  S_approx_owned_.reset(BuildSchurDiagApprox_(*A11, *A12, *A21, *A22));
445  S_approx_view_ = S_approx_owned_.get();
446  } else if (approxType_ == SchurApproxType::A22Only) {
447  S_approx_owned_.reset(new mfem::HypreParMatrix(*A22));
448  S_approx_view_ = S_approx_owned_.get();
449  } else {
450  S_approx_owned_.reset();
451  S_approx_view_ = block_op_overrides_[1].get();
452  }
453 
454  MFEM_VERIFY(S_approx_view_, "Schur complement approximation operator must be set");
455 
456  // Set the Schur complement preconditioner for block (1,1)
457  mfem_solvers_[1]->SetOperator(*S_approx_view_);
458 
459  // Set up block diagonal operator
460  solver_diag_->SetBlock(0, 0, mfem_solvers_[0].get());
461  solver_diag_->SetBlock(1, 1, mfem_solvers_[1].get());
462 }
463 
464 BlockSchurPreconditioner::~BlockSchurPreconditioner() {}
465 } // namespace smith
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.
Definition: smith.cpp:36
constexpr T & get(variant< T0, T1 > &v)
Returns the variant member of specified type.
Definition: variant.hpp:338
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
Definition: tensor.hpp:1939