Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
equation_solver.cpp
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
8 #include "smith/numerics/block_preconditioner.hpp"
9 
10 #include <cstdlib>
11 #include <iomanip>
12 #include <iostream>
13 #include <algorithm>
14 #include <cmath>
15 #include <exception>
16 #include <limits>
17 #include <string>
18 #include <tuple>
19 
20 #include "smith/smith_config.hpp"
22 #include "smith/numerics/trust_region_solver.hpp"
24 
25 namespace smith {
26 
27 namespace {
28 
29 class SolverWithPreconditioner : public mfem::Solver {
30  public:
31  SolverWithPreconditioner(std::unique_ptr<mfem::Solver> linear_solver, std::unique_ptr<mfem::Solver> preconditioner)
32  : linear_solver_(std::move(linear_solver)), preconditioner_(std::move(preconditioner))
33  {
34  SLIC_ERROR_IF(!linear_solver_, "SolverWithPreconditioner requires a non-null linear solver");
35  }
36 
37  void SetOperator(const mfem::Operator& op) override
38  {
39  height = op.Height();
40  width = op.Width();
41  linear_solver_->SetOperator(op);
42  }
43 
44  void Mult(const mfem::Vector& x, mfem::Vector& y) const override { linear_solver_->Mult(x, y); }
45 
46  private:
47  std::unique_ptr<mfem::Solver> linear_solver_;
48  std::unique_ptr<mfem::Solver> preconditioner_;
49 };
50 
51 bool preconditionerSupportsBlockOperator(Preconditioner preconditioner)
52 {
53  switch (preconditioner) {
58  return true;
59  default:
60  return false;
61  }
62 }
63 
64 bool linearSolverSupportsBlockOperator(LinearSolver linear_solver)
65 {
66  switch (linear_solver) {
67  case LinearSolver::CG:
70 #ifdef MFEM_USE_STRUMPACK
72 #endif
73 #ifdef SMITH_USE_PETSC
76 #endif
77  return true;
78  default:
79  return false;
80  }
81 }
82 
83 bool monolithicizeOperatorIfNeeded(const LinearSolverOptions& linear_options, mfem::Operator& assembled_gradient,
84  mfem::Operator*& gradient_operator)
85 {
86  auto* block_gradient = dynamic_cast<const mfem::BlockOperator*>(&assembled_gradient);
87  if (!block_gradient) {
88  gradient_operator = &assembled_gradient;
89  return false;
90  }
91 
92  if (!requiresMonolithicOperator(linear_options)) {
93  gradient_operator = &assembled_gradient;
94  return false;
95  }
96 
97  gradient_operator = buildMonolithicMatrix(*block_gradient).release();
98  SLIC_DEBUG_ROOT(
99  axom::fmt::format("Automatically monolithicizing block Jacobian for linear solver {} with "
100  "preconditioner {}",
101  linearName(linear_options.linear_solver), preconditionerName(linear_options.preconditioner)));
102  return true;
103 }
104 
105 } // namespace
106 
109 class NewtonSolver : public mfem::NewtonSolver, public ConvergenceManagedNonlinearSolver {
110  protected:
112  mutable mfem::Vector x0;
113 
115  NonlinearSolverOptions nonlinear_options;
116 
118  LinearSolverOptions linear_options;
119 
121  mutable size_t print_level = 0;
122 
124  mutable bool grad_monolithic = false;
125 
126  std::shared_ptr<EquationSolverConvergenceManager> convergence_manager_ = nullptr;
127 
128  public:
130  NewtonSolver(const NonlinearSolverOptions& nonlinear_opts, const LinearSolverOptions& linear_opts)
131  : nonlinear_options(nonlinear_opts), linear_options(linear_opts)
132  {
133  }
134 
135 #ifdef MFEM_USE_MPI
137  NewtonSolver(MPI_Comm comm_, const NonlinearSolverOptions& nonlinear_opts, const LinearSolverOptions& linear_opts)
138  : mfem::NewtonSolver(comm_), nonlinear_options(nonlinear_opts), linear_options(linear_opts)
139  {
140  }
141 #endif
142 
144  virtual ~NewtonSolver()
145  {
146  if (grad_monolithic) delete grad;
147  }
148 
149  void setConvergenceManager(std::shared_ptr<EquationSolverConvergenceManager> convergence_manager) override
150  {
151  convergence_manager_ = std::move(convergence_manager);
152  }
153 
155  ConvergenceStatus evaluateConvergence(const mfem::Vector& x, mfem::Vector& rOut) const
156  {
158  ConvergenceStatus status;
159  status.global_norm = std::numeric_limits<double>::max();
160  status.global_goal = std::numeric_limits<double>::max();
161  try {
162  oper->Mult(x, rOut);
163  if (convergence_manager_) {
164  status = convergence_manager_->evaluate(1.0, rOut);
165  } else {
166  status.block_norms = {Norm(rOut)};
167  status.global_norm = status.block_norms.front();
168  status.global_goal = std::max(rel_tol * initial_norm, abs_tol);
169  status.global_converged = status.global_norm <= status.global_goal;
170  status.converged = status.global_converged;
171  }
172  } catch (const std::exception&) {
173  status.global_norm = std::numeric_limits<double>::max();
174  status.global_goal = std::numeric_limits<double>::max();
175  }
176  return status;
177  }
178 
180  void assembleJacobian(const mfem::Vector& x) const
181  {
183  if (grad_monolithic) {
184  delete grad;
185  grad = nullptr;
186  grad_monolithic = false;
187  }
188  mfem::Operator& assembled_gradient = oper->GetGradient(x);
189  grad_monolithic = monolithicizeOperatorIfNeeded(linear_options, assembled_gradient, grad);
190  }
191 
193  void setPreconditioner() const
194  {
196  prec->SetOperator(*grad);
197  }
198 
200  void solveLinearSystem(const mfem::Vector& r_, mfem::Vector& c_) const
201  {
203  prec->Mult(r_, c_); // c = [DF(x_i)]^{-1} [F(x_i)-b]
204  }
205 
207  void Mult(const mfem::Vector&, mfem::Vector& x) const override
208  {
209  MFEM_ASSERT(oper != NULL, "the Operator is not set (use SetOperator).");
210  MFEM_ASSERT(prec != NULL, "the Solver is not set (use SetSolver).");
211 
212  print_level = print_options.iterations ? 1 : print_level;
213  print_level = print_options.summary ? 2 : print_level;
214 
215  using real_t = mfem::real_t;
216 
217  ConvergenceStatus status = evaluateConvergence(x, r);
218  real_t norm = status.global_norm;
219  initial_norm = norm;
220  if (norm == 0.0) return;
221 
222  if (print_level == 1) {
223  mfem::out << "Newton iteration " << std::setw(3) << 0 << " : ||r|| = " << std::setw(13) << norm << "\n";
224  }
225 
226  prec->iterative_mode = false;
227 
228  int it = 0;
229  for (; true; it++) {
230  MFEM_ASSERT(mfem::IsFinite(norm), "norm = " << norm);
231  if (print_level >= 2) {
232  mfem::out << "Newton iteration " << std::setw(3) << it << " : ||r|| = " << std::setw(13) << norm;
233  if (it > 0) {
234  mfem::out << ", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ? norm / initial_norm : norm);
235  }
236  mfem::out << '\n';
237  }
238 
239  if ((print_level >= 1) && (norm != norm)) {
240  mfem::out << "Initial residual for Newton iteration is undefined/nan.\n";
241  mfem::out << "Newton: No convergence!\n";
242  return;
243  }
244 
245  if (status.converged && it >= nonlinear_options.min_iterations) {
246  converged = true;
247  break;
248  } else if (it >= max_iter) {
249  converged = false;
250  break;
251  }
252 
253  real_t norm_nm1 = norm;
254 
255  assembleJacobian(x);
256  setPreconditioner();
257  solveLinearSystem(r, c);
258 
259  // there must be a better way to do this?
260  x0.SetSize(x.Size());
261  x0 = 0.0;
262  x0.Add(1.0, x);
263 
264  real_t stepScale = 1.0;
265  add(x0, -stepScale, c, x);
266  status = evaluateConvergence(x, r);
267  norm = status.global_norm;
268 
269  const int max_ls_iters = nonlinear_options.max_line_search_iterations;
270  static constexpr real_t reduction = 0.5;
271 
272  const double sufficientDecreaseParam = 0.0; // 1e-15;
273  const double cMagnitudeInR = sufficientDecreaseParam != 0.0 ? std::abs(Dot(c, r)) / norm_nm1 : 0.0;
274 
275  auto is_improved = [=](real_t currentNorm, real_t c_scale) {
276  return currentNorm < norm_nm1 - sufficientDecreaseParam * c_scale * cMagnitudeInR;
277  };
278 
279  // back-track linesearch
280  int ls_iter = 0;
281  int ls_iter_sum = 0;
282  for (; !is_improved(norm, stepScale) && ls_iter < max_ls_iters; ++ls_iter, ++ls_iter_sum) {
283  stepScale *= reduction;
284  add(x0, -stepScale, c, x);
285  status = evaluateConvergence(x, r);
286  norm = status.global_norm;
287  }
288 
289  // try the opposite direction and linesearch back from there
290  if (max_ls_iters > 0 && ls_iter == max_ls_iters && !is_improved(norm, stepScale)) {
291  stepScale = 1.0;
292  add(x0, stepScale, c, x);
293  status = evaluateConvergence(x, r);
294  norm = status.global_norm;
295 
296  ls_iter = 0;
297  for (; !is_improved(norm, stepScale) && ls_iter < max_ls_iters; ++ls_iter, ++ls_iter_sum) {
298  stepScale *= reduction;
299  add(x0, stepScale, c, x);
300  status = evaluateConvergence(x, r);
301  norm = status.global_norm;
302  }
303 
304  // ok, the opposite direction was also terrible, lets go back, cut in half 1 last time and accept it hoping for
305  // the best
306  if (ls_iter == max_ls_iters && !is_improved(norm, stepScale)) {
307  ++ls_iter_sum;
308  stepScale *= reduction;
309  add(x0, -stepScale, c, x);
310  status = evaluateConvergence(x, r);
311  norm = status.global_norm;
312  }
313  }
314 
315  if (ls_iter_sum) {
316  if (print_level >= 2) {
317  mfem::out << "Number of line search steps taken = " << ls_iter_sum << std::endl;
318  }
319  if (print_level >= 2 && (ls_iter_sum == 2 * max_ls_iters + 1)) {
320  mfem::out << "The maximum number of line search cut back have occurred, the resulting residual may not have "
321  "decreased. "
322  << std::endl;
323  }
324  }
325  }
326 
327  final_iter = it;
328  final_norm = norm;
329 
330  if (print_level == 1) {
331  mfem::out << "Newton iteration " << std::setw(3) << final_iter << " : ||r|| = " << std::setw(13) << norm << '\n';
332  }
333  if (!converged && print_level >= 1) { // (print_options.summary || print_options.warnings)) {
334  mfem::out << "Newton: No convergence!\n";
335  }
336  }
337 };
338 
340 struct TrustRegionSettings {
342  double cg_tol = 1e-8;
344  size_t min_cg_iterations = 0; //
346  size_t max_cg_iterations = 10000; //
348  size_t max_cumulative_iteration = 1;
350  double min_tr_size = 1e-13;
352  double t1 = 0.25;
354  double t2 = 1.75;
356  double eta1 = 1e-9;
358  double eta2 = 0.1;
360  double eta3 = 0.6;
362  double eta4 = 4.2;
363 };
364 
366 struct TrustRegionResults {
368  TrustRegionResults(int size)
369  {
370  z.SetSize(size);
371  H_z.SetSize(size);
372  d_old.SetSize(size);
373  H_d_old.SetSize(size);
374  d.SetSize(size);
375  H_d.SetSize(size);
376  Pr.SetSize(size);
377  cauchy_point.SetSize(size);
378  H_cauchy_point.SetSize(size);
379  }
380 
382  void reset()
383  {
384  z = 0.0;
385  cauchy_point = 0.0;
386  }
387 
389  enum class Status
390  {
391  Interior,
392  NegativeCurvature,
393  OnBoundary,
394  NonDescentDirection
395  };
396 
398  mfem::Vector z;
400  mfem::Vector H_z;
402  mfem::Vector d_old;
404  mfem::Vector H_d_old;
406  mfem::Vector d;
408  mfem::Vector H_d;
410  mfem::Vector Pr;
412  mfem::Vector cauchy_point;
414  mfem::Vector H_cauchy_point;
416  Status interior_status = Status::Interior;
418  size_t cg_iterations_count = 0;
419 };
420 
422 void printTrustRegionInfo(double realObjective, double modelObjective, size_t cgIters, double trSize, bool willAccept)
423 {
424  mfem::out << "real energy = " << std::setw(13) << realObjective << ", model energy = " << std::setw(13)
425  << modelObjective << ", cg iter = " << std::setw(7) << cgIters << ", next tr size = " << std::setw(8)
426  << trSize << ", accepting = " << willAccept << std::endl;
427 }
428 
438 class TrustRegion : public mfem::NewtonSolver, public ConvergenceManagedNonlinearSolver {
439  protected:
441  mutable mfem::Vector x_pred;
443  mutable mfem::Vector r_pred;
445  mutable mfem::Vector scratch;
447  mutable std::vector<std::shared_ptr<mfem::Vector>> left_mosts;
449  mutable std::vector<std::shared_ptr<mfem::Vector>> H_left_mosts;
450 
452  NonlinearSolverOptions nonlinear_options;
454  LinearSolverOptions linear_options;
455 
458  Solver& tr_precond;
459 
461  mutable size_t print_level = 0;
462 
464  mutable bool grad_monolithic = false;
465 
466  std::shared_ptr<EquationSolverConvergenceManager> convergence_manager_ = nullptr;
467 
468  public:
470  mutable size_t num_hess_vecs = 0;
472  mutable size_t num_preconds = 0;
474  mutable size_t num_residuals = 0;
476  mutable size_t num_subspace_solves = 0;
478  mutable size_t num_jacobian_assembles = 0;
479 
480 #ifdef MFEM_USE_MPI
482  TrustRegion(MPI_Comm comm_, const NonlinearSolverOptions& nonlinear_opts, const LinearSolverOptions& linear_opts,
483  Solver& tPrec)
484  : mfem::NewtonSolver(comm_), nonlinear_options(nonlinear_opts), linear_options(linear_opts), tr_precond(tPrec)
485  {
486  }
487 #endif
488 
490  virtual ~TrustRegion()
491  {
492  if (grad_monolithic) delete grad;
493  }
494 
495  void setConvergenceManager(std::shared_ptr<EquationSolverConvergenceManager> convergence_manager) override
496  {
497  convergence_manager_ = std::move(convergence_manager);
498  }
499 
501  void projectToBoundaryWithCoefs(mfem::Vector& z, const mfem::Vector& d, double delta, double zz, double zd,
502  double dd) const
503  {
504  // find z + tau d
505  double deltadelta_m_zz = delta * delta - zz;
506  if (deltadelta_m_zz == 0) return; // already on boundary
507  double tau = (std::sqrt(deltadelta_m_zz * dd + zd * zd) - zd) / dd;
508  z.Add(tau, d);
509  }
510 
512  template <typename HessVecFunc>
513  void solveTheSubspaceProblem([[maybe_unused]] mfem::Vector& z, [[maybe_unused]] const HessVecFunc& hess_vec_func,
514  [[maybe_unused]] const std::vector<const mfem::Vector*> ds,
515  [[maybe_unused]] const std::vector<const mfem::Vector*> Hds,
516  [[maybe_unused]] const mfem::Vector& g, [[maybe_unused]] double delta,
517  [[maybe_unused]] int num_leftmost) const
518  {
519 #ifdef SMITH_USE_SLEPC
521  ++num_subspace_solves;
522 
523  std::vector<const mfem::Vector*> directions;
524  for (auto& d : ds) {
525  directions.emplace_back(d);
526  }
527  for (auto& left : left_mosts) {
528  directions.emplace_back(left.get());
529  }
530 
531  std::vector<const mfem::Vector*> H_directions;
532  for (auto& Hd : Hds) {
533  H_directions.emplace_back(Hd);
534  }
535  for (auto& H_left : H_left_mosts) {
536  H_directions.emplace_back(H_left.get());
537  }
538 
539  try {
540  std::tie(directions, H_directions) = removeDependentDirections(directions, H_directions);
541  } catch (const std::exception& e) {
542  if (print_level >= 2) {
543  mfem::out << "remove dependent directions failed with " << e.what() << std::endl;
544  }
545  return;
546  }
547 
548  mfem::Vector b(g);
549  b *= -1;
550 
551  mfem::Vector sol;
552  std::vector<std::shared_ptr<mfem::Vector>> leftvecs;
553  std::vector<double> leftvals;
554  double energy_change;
555 
556  try {
557  std::tie(sol, leftvecs, leftvals, energy_change) =
558  solveSubspaceProblem(directions, H_directions, b, delta, num_leftmost);
559  } catch (const std::exception& e) {
560  if (print_level == 1) {
561  mfem::out << "subspace solve failed with " << e.what() << std::endl;
562  }
563  return;
564  }
565 
566  left_mosts.clear();
567  for (auto& lv : leftvecs) {
568  left_mosts.emplace_back(std::move(lv));
569  }
570 
571  double base_energy = computeEnergy(g, hess_vec_func, z);
572  double subspace_energy = computeEnergy(g, hess_vec_func, sol);
573 
574  if (print_level >= 2) {
575  double leftval = leftvals.size() ? leftvals[0] : 1.0;
576  mfem::out << "Energy using subspace solver from: " << base_energy << ", to: " << subspace_energy << " / "
577  << energy_change << ". Min eig: " << leftval << std::endl;
578  }
579 
580  if (subspace_energy < base_energy) {
581  z = sol;
582  }
583 #endif
584  }
585 
587  void projectToBoundaryBetweenWithCoefs(mfem::Vector& z, const mfem::Vector& y, double trSize, double zz, double zy,
588  double yy) const
589  {
590  double dd = yy - 2 * zy + zz;
591  double zd = zy - zz;
592  double tau = (std::sqrt((trSize * trSize - zz) * dd + zd * zd) - zd) / dd;
593  z.Add(-tau, z);
594  z.Add(tau, y);
595  }
596 
598  void doglegStep(const mfem::Vector& cp, const mfem::Vector& newtonP, double trSize, mfem::Vector& s) const
599  {
601  // MRT, could optimize some of these eventually, compute on the outside and save
602  double cc = Dot(cp, cp);
603  double nn = Dot(newtonP, newtonP);
604  double tt = trSize * trSize;
605 
606  s = 0.0;
607  if (cc >= tt) {
608  add(s, std::sqrt(tt / cc), cp, s);
609  } else if (cc > nn) {
610  if (print_level >= 2) {
611  mfem::out << "cp outside newton, preconditioner likely inaccurate\n";
612  }
613  add(s, 1.0, cp, s);
614  } else if (nn > tt) { // on the dogleg (we have nn >= cc, and tt >= cc)
615  add(s, 1.0, cp, s);
616  double cn = Dot(cp, newtonP);
617  projectToBoundaryBetweenWithCoefs(s, newtonP, trSize, cc, cn, nn);
618  } else {
619  s = newtonP;
620  }
621  }
622 
624  template <typename HessVecFunc>
625  double computeEnergy(const mfem::Vector& r_local, const HessVecFunc& H, const mfem::Vector& z) const
626  {
628  double rz = Dot(r_local, z);
629  mfem::Vector tmp(r_local);
630  tmp = 0.0;
631  H(z, tmp);
632  return rz + 0.5 * Dot(z, tmp);
633  }
634 
636  template <typename HessVecFunc, typename PrecondFunc>
637  void solveTrustRegionModelProblem(const mfem::Vector& r0, mfem::Vector& rCurrent, HessVecFunc hess_vec_func,
638  PrecondFunc precond, const TrustRegionSettings& settings, double& trSize,
639  TrustRegionResults& results) const
640  {
642  // minimize r0@z + 0.5*z@J@z
643  results.interior_status = TrustRegionResults::Status::Interior;
644  results.cg_iterations_count = 0;
645 
646  auto& z = results.z;
647  auto& cgIter = results.cg_iterations_count;
648  auto& d = results.d;
649  auto& Pr = results.Pr;
650  auto& Hd = results.H_d;
651 
652  const double cg_tol_squared = settings.cg_tol * settings.cg_tol;
653 
654  if (Dot(r0, r0) <= cg_tol_squared && settings.min_cg_iterations == 0) {
655  if (print_level >= 2) {
656  mfem::out << "Trust region solution state within tolerance on first iteration."
657  << "\n";
658  }
659  return;
660  }
661 
662  rCurrent = r0;
663  precond(rCurrent, Pr);
664 
665  // d = -Pr
666  d = Pr;
667  d *= -1.0;
668 
669  z = 0.0;
670  double zz = 0.;
671  double rPr = Dot(rCurrent, Pr);
672  double zd = 0.0;
673  double dd = Dot(d, d);
674 
675  // std::cout << "initial energy = " << computeEnergy(r0, hess_vec_func, z) << std::endl;
676 
677  for (cgIter = 1; cgIter <= settings.max_cg_iterations; ++cgIter) {
678  // check if this is a descent direction
679  if (Dot(d, rCurrent) > 0) {
680  d *= -1;
681  results.interior_status = TrustRegionResults::Status::NonDescentDirection;
682  }
683 
684  hess_vec_func(d, Hd);
685  const double curvature = Dot(d, Hd);
686  const double alphaCg = curvature != 0.0 ? rPr / curvature : 0.0;
687 
688  auto& zPred = Pr; // re-use Pr memory.
689  // This predicted step will no longer be used by the time Pr is, so we can avoid an extra
690  // vector floating around
691  add(z, alphaCg, d, zPred);
692  double zzNp1 = Dot(zPred, zPred);
693 
694  const bool go_to_boundary = curvature <= 0 || zzNp1 >= trSize * trSize;
695  if (go_to_boundary) {
696  projectToBoundaryWithCoefs(z, d, trSize, zz, zd, dd);
697  if (curvature <= 0) {
698  results.interior_status = TrustRegionResults::Status::NegativeCurvature;
699  } else {
700  results.interior_status = TrustRegionResults::Status::OnBoundary;
701  }
702  return;
703  }
704 
705  z = zPred;
706 
707  if (results.interior_status == TrustRegionResults::Status::NonDescentDirection) {
708  if (print_level >= 2) {
709  mfem::out << "Found a non descent direction\n";
710  }
711  return;
712  }
713 
714  add(rCurrent, alphaCg, Hd, rCurrent);
715 
716  precond(rCurrent, Pr);
717  double rPrNp1 = Dot(rCurrent, Pr);
718 
719  if (Dot(rCurrent, rCurrent) <= cg_tol_squared && cgIter >= settings.min_cg_iterations) {
720  return;
721  }
722 
723  double beta = rPrNp1 / rPr;
724  rPr = rPrNp1;
725  add(-1.0, Pr, beta, d, d);
726 
727  zz = zzNp1;
728  zd = Dot(z, d);
729  dd = Dot(d, d);
730  }
731  cgIter--; // if all cg iterations are taken, correct for output
732  }
733 
735  void assembleJacobian(const mfem::Vector& x) const
736  {
738  ++num_jacobian_assembles;
739  if (grad_monolithic) {
740  delete grad;
741  grad = nullptr;
742  grad_monolithic = false;
743  }
744  mfem::Operator& assembled_gradient = oper->GetGradient(x);
745  grad_monolithic = monolithicizeOperatorIfNeeded(linear_options, assembled_gradient, grad);
746  }
747 
749  mfem::real_t computeResidual(const mfem::Vector& x_, mfem::Vector& r_) const
750  {
752  ++num_residuals;
753  oper->Mult(x_, r_);
754  return Norm(r_);
755  }
756 
757  ConvergenceStatus evaluateConvergence(const mfem::Vector& x_, mfem::Vector& r_) const
758  {
759  ConvergenceStatus status;
760  status.global_norm = std::numeric_limits<double>::max();
761  status.global_goal = std::numeric_limits<double>::max();
762  try {
763  status.global_norm = computeResidual(x_, r_);
764  if (convergence_manager_) {
765  status = convergence_manager_->evaluate(1.0, r_);
766  } else {
767  status.block_norms = {status.global_norm};
768  status.global_goal = std::max(rel_tol * initial_norm, abs_tol);
769  status.global_converged = status.global_norm <= status.global_goal;
770  status.converged = status.global_converged;
771  }
772  } catch (const std::exception&) {
773  status.global_norm = std::numeric_limits<double>::max();
774  status.global_goal = std::numeric_limits<double>::max();
775  }
776  return status;
777  }
778 
780  void hessVec(const mfem::Vector& x_, mfem::Vector& v_) const
781  {
783  ++num_hess_vecs;
784  grad->Mult(x_, v_);
785  }
786 
788  void precond(const mfem::Vector& x_, mfem::Vector& v_) const
789  {
791  ++num_preconds;
792  tr_precond.Mult(x_, v_);
793  };
794 
796  void Mult(const mfem::Vector&, mfem::Vector& X) const override
797  {
798  MFEM_ASSERT(oper != NULL, "the Operator is not set (use SetOperator).");
799  MFEM_ASSERT(prec != NULL, "the Solver is not set (use SetSolver).");
800 
801  print_level = print_options.iterations ? 1 : print_level;
802  print_level = print_options.summary ? 2 : print_level;
803 
804  using real_t = mfem::real_t;
805 
806  num_hess_vecs = 0;
807  num_preconds = 0;
808  num_residuals = 0;
809  num_subspace_solves = 0;
810  num_jacobian_assembles = 0;
811 
812  ConvergenceStatus status = evaluateConvergence(X, r);
813  real_t norm = status.global_norm;
814  real_t norm_goal = status.global_goal;
815  initial_norm = norm;
816  if (norm == 0.0) return;
817 
818  if (print_level == 1) {
819  mfem::out << "TrustRegion iteration " << std::setw(3) << 0 << " : ||r|| = " << std::setw(13) << norm << "\n";
820  }
821 
822  prec->iterative_mode = false;
823  tr_precond.iterative_mode = false;
824 
825  // local arrays
826  x_pred.SetSize(X.Size());
827  x_pred = 0.0;
828  r_pred.SetSize(X.Size());
829  r_pred = 0.0;
830  scratch.SetSize(X.Size());
831  scratch = 0.0;
832 
833  TrustRegionResults trResults(X.Size());
834  TrustRegionSettings settings;
835  settings.min_cg_iterations = static_cast<size_t>(nonlinear_options.min_iterations);
836  settings.max_cg_iterations = static_cast<size_t>(linear_options.max_iterations);
837  settings.cg_tol = 0.5 * norm_goal;
838 
839  int subspace_option = nonlinear_options.subspace_option;
840  int num_leftmost = nonlinear_options.num_leftmost;
841 
842  scratch = 1.0;
843  double tr_size = nonlinear_options.trust_region_scaling * std::sqrt(Dot(scratch, scratch));
844  size_t cumulative_cg_iters_from_last_precond_update = 0;
845 
846  int it = 0;
847  for (; true; it++) {
848  MFEM_ASSERT(mfem::IsFinite(norm), "norm = " << norm);
849  if (print_level >= 2) {
850  mfem::out << "TrustRegion iteration " << std::setw(3) << it << " : ||r|| = " << std::setw(13) << norm;
851  if (it > 0) {
852  mfem::out << ", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ? norm / initial_norm : norm);
853  mfem::out << ", x_incr = " << std::setw(13) << trResults.d.Norml2();
854  } else {
855  mfem::out << ", norm goal = " << std::setw(13) << norm_goal;
856  }
857  mfem::out << '\n';
858  }
859 
860  if (print_level >= 1 && (norm != norm)) {
861  mfem::out << "Initial residual for trust-region iteration is undefined/nan." << std::endl;
862  mfem::out << "TrustRegion: No convergence!\n";
863  return;
864  }
865 
866  if (status.converged && it >= nonlinear_options.min_iterations) {
867  converged = true;
868  break;
869  } else if (it >= max_iter) {
870  converged = false;
871  break;
872  }
873 
874  assembleJacobian(X);
875 
876  if (it == 0 || (trResults.cg_iterations_count >= settings.max_cg_iterations ||
877  cumulative_cg_iters_from_last_precond_update >= settings.max_cumulative_iteration)) {
878  tr_precond.SetOperator(*grad);
879  cumulative_cg_iters_from_last_precond_update = 0;
880  }
881 
882  auto hess_vec_func = [&](const mfem::Vector& x_, mfem::Vector& v_) { hessVec(x_, v_); };
883  auto precond_func = [&](const mfem::Vector& x_, mfem::Vector& v_) { precond(x_, v_); };
884 
885  double cauchyPointNormSquared = tr_size * tr_size;
886  trResults.reset();
887 
888  hess_vec_func(r, trResults.H_d);
889  const double gKg = Dot(r, trResults.H_d);
890  if (gKg > 0) {
891  const double alphaCp = -Dot(r, r) / gKg;
892  add(trResults.cauchy_point, alphaCp, r, trResults.cauchy_point);
893  cauchyPointNormSquared = Dot(trResults.cauchy_point, trResults.cauchy_point);
894  } else {
895  const double alphaTr = -tr_size / std::sqrt(Dot(r, r));
896  add(trResults.cauchy_point, alphaTr, r, trResults.cauchy_point);
897  if (print_level >= 2) {
898  mfem::out << "Negative curvature un-preconditioned cauchy point direction found."
899  << "\n";
900  }
901  }
902 
903  if (cauchyPointNormSquared >= tr_size * tr_size) {
904  if (print_level >= 2) {
905  mfem::out << "Un-preconditioned gradient cauchy point outside trust region, step size = "
906  << std::sqrt(cauchyPointNormSquared) << "\n";
907  }
908  trResults.cauchy_point *= (tr_size / std::sqrt(cauchyPointNormSquared));
909  trResults.z = trResults.cauchy_point;
910 
911  trResults.cg_iterations_count = 1;
912  trResults.interior_status = TrustRegionResults::Status::OnBoundary;
913  } else {
914  settings.cg_tol = std::max(0.5 * norm_goal, 5e-5 * norm);
915  solveTrustRegionModelProblem(r, scratch, hess_vec_func, precond_func, settings, tr_size, trResults);
916  }
917  cumulative_cg_iters_from_last_precond_update += trResults.cg_iterations_count;
918 
919  bool have_computed_Hvs = false;
920 
921  int lineSearchIter = 0;
922  while (lineSearchIter <= nonlinear_options.max_line_search_iterations) {
923  ++lineSearchIter;
924 
925  doglegStep(trResults.cauchy_point, trResults.z, tr_size, trResults.d);
926 
927  bool use_with_option1 =
928  (subspace_option >= 1) && (trResults.interior_status == TrustRegionResults::Status::NonDescentDirection ||
929  trResults.interior_status == TrustRegionResults::Status::NegativeCurvature ||
930  ((Norm(trResults.d) > (1.0 - 1.0e-6) * tr_size) && lineSearchIter > 1));
931  bool use_with_option2 = (subspace_option >= 2) && (Norm(trResults.d) > (1.0 - 1.0e-6) * tr_size);
932  bool use_with_option3 = (subspace_option >= 3);
933 
934  if (use_with_option1 || use_with_option2 || use_with_option3) {
935  if (!have_computed_Hvs) {
936  have_computed_Hvs = true;
937  hess_vec_func(trResults.z, trResults.H_z);
938  hess_vec_func(trResults.d_old, trResults.H_d_old);
939  hess_vec_func(trResults.cauchy_point, trResults.H_cauchy_point);
940  }
941 
942  H_left_mosts.clear();
943  for (auto& left : left_mosts) {
944  H_left_mosts.emplace_back(std::make_shared<mfem::Vector>(*left));
945  hess_vec_func(*left, *H_left_mosts.back());
946  }
947 
948  std::vector<const mfem::Vector*> ds{&trResults.z, &trResults.d_old, &trResults.cauchy_point};
949  std::vector<const mfem::Vector*> H_ds{&trResults.H_z, &trResults.H_d_old, &trResults.H_cauchy_point};
950  solveTheSubspaceProblem(trResults.d, hess_vec_func, ds, H_ds, r, tr_size, num_leftmost);
951  }
952 
953  static constexpr double roundOffTol = 0.0; // 1e-14;
954 
955  hess_vec_func(trResults.d, trResults.H_d);
956  double dHd = Dot(trResults.d, trResults.H_d);
957  double modelObjective = Dot(r, trResults.d) + 0.5 * dHd - roundOffTol;
958 
959  add(X, trResults.d, x_pred);
960 
961  double realObjective = std::numeric_limits<double>::max();
962  double normPred = std::numeric_limits<double>::max();
963  try {
964  auto predicted_status = evaluateConvergence(x_pred, r_pred);
965  normPred = predicted_status.global_norm;
966  double obj1 = 0.5 * (Dot(r, trResults.d) + Dot(r_pred, trResults.d)) - roundOffTol;
967  realObjective = obj1;
968  if (predicted_status.converged) {
969  trResults.d_old = trResults.d;
970  X = x_pred;
971  r = r_pred;
972  status = predicted_status;
973  norm = status.global_norm;
974  if (print_level >= 2) {
975  printTrustRegionInfo(realObjective, modelObjective, trResults.cg_iterations_count, tr_size, true);
976  trResults.cg_iterations_count = 0;
977  }
978  break;
979  }
980  } catch (const std::exception&) {
981  realObjective = std::numeric_limits<double>::max();
983  }
984 
985  double modelImprove = -modelObjective;
986  double realImprove = -realObjective;
987 
988  double rho = realImprove / modelImprove;
989  if (modelObjective > 0) {
990  if (print_level >= 2) {
991  mfem::out << "Found a positive model objective increase. Debug if you see this.\n";
992  }
993  rho = realImprove / -modelImprove;
994  }
995 
996  // std::cout << "rho , stuff = " << rho << " " << settings.eta3 << std::endl;
997  // std::cout << "stat = "<< trResults.interior_status << std::endl;
998 
999  if (!(rho >= settings.eta2) ||
1000  rho > settings.eta4) { // not enough progress, decrease trust region. write it this way to handle NaNs.
1001  tr_size *= settings.t1;
1002  } else if ((rho > settings.eta3 && rho <= settings.eta4 &&
1003  trResults.interior_status == TrustRegionResults::Status::OnBoundary) ||
1004  (rho > 0.95 && rho < 1.05 &&
1005  trResults.interior_status ==
1006  TrustRegionResults::Status::NegativeCurvature)) { // good progress, on boundary, increase trust
1007  // region
1008  tr_size *= settings.t2;
1009  }
1010 
1011  // eventually extend to handle this case to handle occasional roundoff issues
1012  // modelRes = g + Jd
1013  // modelResNorm = np.linalg.norm(modelRes)
1014  // realResNorm = np.linalg.norm(gy)
1015  bool willAccept = rho >= settings.eta1 && rho <= settings.eta4; // or (rho >= -0 and realResNorm <= gNorm)
1016 
1017  if (print_level >= 2) {
1018  printTrustRegionInfo(realObjective, modelObjective, trResults.cg_iterations_count, tr_size, willAccept);
1019  trResults.cg_iterations_count =
1020  0; // zero this output so it doesn't look like the linesearch is doing cg iterations
1021  }
1022 
1023  if (willAccept) {
1024  trResults.d_old = trResults.d;
1025  X = x_pred;
1026  r = r_pred;
1027  status = convergence_manager_ ? convergence_manager_->evaluate(1.0, r_pred) : status;
1028  norm = normPred;
1029  break;
1030  }
1031  }
1032  }
1033 
1034  final_iter = it;
1035  final_norm = norm;
1036 
1037  if (print_level == 1) {
1038  mfem::out << "TrustRegion iteration " << std::setw(3) << final_iter << " : ||r|| = " << std::setw(13) << norm
1039  << '\n';
1040  }
1041  if (!converged && print_level >= 1) { // (print_options.summary || print_options.warnings)) {
1042  mfem::out << "TrustRegion: No convergence!\n";
1043  }
1044 
1045  if (false && print_level >= 2) {
1046  mfem::out << "num hess vecs = " << num_hess_vecs << "\n";
1047  mfem::out << "num preconds = " << num_preconds << "\n";
1048  mfem::out << "num residuals = " << num_residuals << "\n";
1049  mfem::out << "num subspace solves = " << num_subspace_solves << "\n";
1050  mfem::out << "num jacobian_assembles = " << num_jacobian_assembles << "\n";
1051  }
1052  }
1053 };
1055 
1057 {
1058  auto [lin_solver, preconditioner] = buildLinearSolverAndPreconditioner(lin_opts, comm);
1059 
1060  lin_solver_ = std::move(lin_solver);
1061  preconditioner_ = std::move(preconditioner);
1062  nonlin_solver_ = buildNonlinearSolver(nonlinear_opts, lin_opts, *preconditioner_, comm);
1063  convergence_manager_ = std::make_shared<EquationSolverConvergenceManager>(comm, nonlinear_opts.absolute_tol,
1064  nonlinear_opts.relative_tol);
1065  attachConvergenceManager();
1066 }
1067 
1068 EquationSolver::EquationSolver(std::unique_ptr<mfem::NewtonSolver> nonlinear_solver,
1069  std::unique_ptr<mfem::Solver> linear_solver,
1070  std::unique_ptr<mfem::Solver> preconditioner)
1071 {
1072  SLIC_ERROR_ROOT_IF(!nonlinear_solver, "Nonlinear solvers must be given to construct an EquationSolver");
1073  SLIC_ERROR_ROOT_IF(!linear_solver, "Linear solvers must be given to construct an EquationSolver");
1074 
1075  nonlin_solver_ = std::move(nonlinear_solver);
1076  lin_solver_ = std::move(linear_solver);
1077  preconditioner_ = std::move(preconditioner);
1078 }
1079 
1080 void EquationSolver::attachConvergenceManager() const
1081 {
1082  if (!convergence_manager_ || !nonlin_solver_) {
1083  return;
1084  }
1085 
1086  if (auto* managed_solver = dynamic_cast<ConvergenceManagedNonlinearSolver*>(nonlin_solver_.get())) {
1087  managed_solver->setConvergenceManager(convergence_manager_);
1088  }
1089 }
1090 
1091 void EquationSolver::initializeConvergenceManager(double abs_tol, double rel_tol, MPI_Comm comm) const
1092 {
1093  if (!convergence_manager_) {
1094  convergence_manager_ = std::make_shared<EquationSolverConvergenceManager>(comm, abs_tol, rel_tol);
1095  } else {
1096  convergence_manager_->setTolerances(abs_tol, rel_tol);
1097  }
1098  attachConvergenceManager();
1099 }
1100 
1101 void EquationSolver::setOperator(const mfem::Operator& op)
1102 {
1103  attachConvergenceManager();
1104  nonlin_solver_->SetOperator(op);
1105 
1106  // Now that the nonlinear solver knows about the operator, we can set its linear solver
1107  if (!nonlin_solver_set_solver_called_) {
1108  nonlin_solver_->SetSolver(linearSolver());
1109  nonlin_solver_set_solver_called_ = true;
1110  }
1111 }
1112 
1113 void EquationSolver::setConvergenceTolerances(double abs_tol, double rel_tol, MPI_Comm comm) const
1114 {
1115  initializeConvergenceManager(abs_tol, rel_tol, comm);
1116 }
1117 
1119 {
1120  if (convergence_manager_) {
1121  convergence_manager_->reset();
1122  }
1123 }
1124 
1125 void EquationSolver::solve(mfem::Vector& x) const
1126 {
1128  mfem::Vector zero(x);
1129  zero = 0.0;
1130  // KINSOL does not handle non-zero RHS, so we enforce that the RHS
1131  // of the nonlinear system is zero
1132  nonlin_solver_->Mult(zero, x);
1133 }
1134 
1135 void SuperLUSolver::Mult(const mfem::Vector& input, mfem::Vector& output) const
1136 {
1137  SLIC_ERROR_ROOT_IF(!superlu_mat_, "Operator must be set prior to solving with SuperLU");
1138 
1139  // Use the underlying MFEM-based solver and SuperLU matrix type to solve the system
1140  superlu_solver_.Mult(input, output);
1141 }
1142 
1158 std::unique_ptr<mfem::HypreParMatrix> buildMonolithicMatrix(const mfem::BlockOperator& block_operator)
1159 {
1160  int row_blocks = block_operator.NumRowBlocks();
1161  int col_blocks = block_operator.NumColBlocks();
1162 
1163  SLIC_ERROR_ROOT_IF(row_blocks != col_blocks, "Attempted to use a direct solver on a non-square block system.");
1164 
1165  mfem::Array2D<const mfem::HypreParMatrix*> hypre_blocks(row_blocks, col_blocks);
1166 
1167  for (int i = 0; i < row_blocks; ++i) {
1168  for (int j = 0; j < col_blocks; ++j) {
1169  // checks for presence of empty (null) blocks, which happen fairly common in multirank contact
1170  if (!block_operator.IsZeroBlock(i, j)) {
1171  auto* hypre_block = dynamic_cast<const mfem::HypreParMatrix*>(&block_operator.GetBlock(i, j));
1172  SLIC_ERROR_ROOT_IF(!hypre_block,
1173  "Trying to use SuperLU on a block operator that does not contain HypreParMatrix blocks.");
1174 
1175  hypre_blocks(i, j) = hypre_block;
1176  } else {
1177  hypre_blocks(i, j) = nullptr;
1178  }
1179  }
1180  }
1181 
1182  // Note that MFEM passes ownership of this matrix to the caller.
1183  // MFEM creates a new monolithic matrix (not a view), so this is a COPY operation.
1184  return std::unique_ptr<mfem::HypreParMatrix>(mfem::HypreParMatrixFromBlocks(hypre_blocks));
1185 }
1186 
1187 void SuperLUSolver::SetOperator(const mfem::Operator& op)
1188 {
1189  // Check if this is a block operator
1190  auto* block_operator = dynamic_cast<const mfem::BlockOperator*>(&op);
1191 
1192  // If it is, make a monolithic system from the underlying blocks
1193  if (block_operator) {
1194  monolithic_mat_ = buildMonolithicMatrix(*block_operator);
1195 
1196  superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*monolithic_mat_);
1197  } else {
1198  // If this is not a block system, check that the input operator is a HypreParMatrix as expected
1199  auto* matrix = dynamic_cast<const mfem::HypreParMatrix*>(&op);
1200 
1201  SLIC_ERROR_ROOT_IF(!matrix, "Matrix must be an assembled HypreParMatrix for use with SuperLU");
1202 
1203  superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*matrix);
1204  }
1205  height = op.Height();
1206  width = op.Width();
1207  superlu_solver_.SetOperator(*superlu_mat_);
1208 }
1209 
1210 #ifdef MFEM_USE_STRUMPACK
1211 
1212 void StrumpackSolver::Mult(const mfem::Vector& input, mfem::Vector& output) const
1213 {
1214  SLIC_ERROR_ROOT_IF(!strumpack_mat_, "Operator must be set prior to solving with Strumpack");
1215 
1216  // Use the underlying MFEM-based solver and Strumpack matrix type to solve the system
1217  strumpack_solver_.Mult(input, output);
1218 }
1219 
1220 void StrumpackSolver::SetOperator(const mfem::Operator& op)
1221 {
1222  // Check if this is a block operator
1223  auto* block_operator = dynamic_cast<const mfem::BlockOperator*>(&op);
1224 
1225  // If it is, make a monolithic system from the underlying blocks
1226  if (block_operator) {
1227  monolithic_mat_ = buildMonolithicMatrix(*block_operator);
1228 
1229  strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*monolithic_mat_);
1230  } else {
1231  // If this is not a block system, check that the input operator is a HypreParMatrix as expected
1232  auto* matrix = dynamic_cast<const mfem::HypreParMatrix*>(&op);
1233 
1234  SLIC_ERROR_ROOT_IF(!matrix, "Matrix must be an assembled HypreParMatrix for use with Strumpack");
1235 
1236  strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*matrix);
1237  }
1238  height = op.Height();
1239  width = op.Width();
1240  strumpack_solver_.SetOperator(*strumpack_mat_);
1241 }
1242 
1243 #endif
1244 
1245 std::unique_ptr<mfem::NewtonSolver> buildNonlinearSolver(NonlinearSolverOptions nonlinear_opts,
1246  const LinearSolverOptions& linear_opts, mfem::Solver& prec,
1247  MPI_Comm comm)
1248 {
1249  std::unique_ptr<mfem::NewtonSolver> nonlinear_solver;
1250 
1251  if (nonlinear_opts.nonlin_solver == NonlinearSolver::Newton) {
1252  nonlinear_opts.max_line_search_iterations = 0;
1253  SLIC_ERROR_ROOT_IF(nonlinear_opts.min_iterations != 0, "Newton's method does not support nonzero min_iterations");
1254  nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts, linear_opts);
1255  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::LBFGS) {
1256  nonlinear_opts.max_line_search_iterations = 0;
1257  SLIC_ERROR_ROOT_IF(nonlinear_opts.min_iterations != 0, "LBFGS does not support nonzero min_iterations");
1258  nonlinear_solver = std::make_unique<mfem::LBFGSSolver>(comm);
1259  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::NewtonLineSearch) {
1260  nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts, linear_opts);
1261  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::TrustRegion) {
1262  nonlinear_solver = std::make_unique<TrustRegion>(comm, nonlinear_opts, linear_opts, prec);
1263 #ifdef SMITH_USE_PETSC
1264  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscNewton) {
1265  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1266  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscNewtonBacktracking) {
1267  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1268  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscNewtonCriticalPoint) {
1269  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1270  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscTrustRegion) {
1271  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1272 #endif
1273  }
1274  // KINSOL
1275  else {
1276 #ifdef SMITH_USE_SUNDIALS
1277  nonlinear_opts.max_line_search_iterations = 0;
1278  SLIC_ERROR_ROOT_IF(nonlinear_opts.min_iterations != 0, "kinsol solvers do not support min_iterations");
1279 
1280  int kinsol_strat = KIN_NONE;
1281 
1282  switch (nonlinear_opts.nonlin_solver) {
1284  kinsol_strat = KIN_NONE;
1285  break;
1287  kinsol_strat = KIN_LINESEARCH;
1288  break;
1290  kinsol_strat = KIN_PICARD;
1291  break;
1292  default:
1293  kinsol_strat = KIN_NONE;
1294  SLIC_ERROR_ROOT("Unknown KINSOL nonlinear solver type given.");
1295  }
1296  auto kinsol_solver = std::make_unique<mfem::KINSolver>(comm, kinsol_strat, true);
1297  nonlinear_solver = std::move(kinsol_solver);
1298 #else
1299  SLIC_ERROR_ROOT("KINSOL was not enabled when MFEM was built");
1300 #endif
1301  }
1302 
1303  nonlinear_solver->SetRelTol(nonlinear_opts.relative_tol);
1304  nonlinear_solver->SetAbsTol(nonlinear_opts.absolute_tol);
1305  nonlinear_solver->SetMaxIter(nonlinear_opts.max_iterations);
1306  nonlinear_solver->SetPrintLevel(nonlinear_opts.print_level);
1307 
1308  // Iterative mode indicates we do not zero out the initial guess during the
1309  // nonlinear solver call. This is required as we apply the essential boundary
1310  // conditions before the nonlinear solver is applied.
1311  nonlinear_solver->iterative_mode = true;
1312 
1313  return nonlinear_solver;
1314 }
1315 
1316 std::pair<std::unique_ptr<mfem::Solver>, std::unique_ptr<mfem::Solver>> buildLinearSolverAndPreconditioner(
1317  LinearSolverOptions linear_opts, MPI_Comm comm)
1318 {
1319  auto preconditioner = buildPreconditioner(linear_opts, comm);
1320 
1321  if (linear_opts.linear_solver == LinearSolver::SuperLU) {
1322  auto lin_solver = std::make_unique<SuperLUSolver>(linear_opts.print_level, comm);
1323  return {std::move(lin_solver), std::move(preconditioner)};
1324  }
1325 
1326 #ifdef MFEM_USE_STRUMPACK
1327 
1328  if (linear_opts.linear_solver == LinearSolver::Strumpack) {
1329  auto lin_solver = std::make_unique<StrumpackSolver>(linear_opts.print_level, comm);
1330  return {std::move(lin_solver), std::move(preconditioner)};
1331  }
1332 
1333 #endif
1334 
1335  std::unique_ptr<mfem::IterativeSolver> iter_lin_solver;
1336 
1337  switch (linear_opts.linear_solver) {
1338  case LinearSolver::CG:
1339  iter_lin_solver = std::make_unique<mfem::CGSolver>(comm);
1340  break;
1341  case LinearSolver::GMRES:
1342  iter_lin_solver = std::make_unique<mfem::GMRESSolver>(comm);
1343  break;
1344 #ifdef SMITH_USE_PETSC
1345  case LinearSolver::PetscCG:
1346  iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPCG, std::string());
1347  break;
1349  iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPGMRES, std::string());
1350  break;
1351 #else
1352  case LinearSolver::PetscCG:
1354  SLIC_ERROR_ROOT("PETSc linear solver requested for non-PETSc build.");
1355  exit(1);
1356  break;
1357 #endif
1358  default:
1359  SLIC_ERROR_ROOT("Linear solver type not recognized.");
1360  exit(1);
1361  }
1362 
1363  iter_lin_solver->SetRelTol(linear_opts.relative_tol);
1364  iter_lin_solver->SetAbsTol(linear_opts.absolute_tol);
1365  iter_lin_solver->SetMaxIter(linear_opts.max_iterations);
1366  iter_lin_solver->SetPrintLevel(linear_opts.print_level);
1367 
1368  if (preconditioner) {
1369  iter_lin_solver->SetPreconditioner(*preconditioner);
1370  }
1371 
1372  return {std::move(iter_lin_solver), std::move(preconditioner)};
1373 }
1374 
1376 {
1377  return !linearSolverSupportsBlockOperator(linear_opts.linear_solver) ||
1378  !preconditionerSupportsBlockOperator(linear_opts.preconditioner);
1379 }
1380 
1381 #ifdef MFEM_USE_AMGX
1382 std::unique_ptr<mfem::AmgXSolver> buildAMGX(const AMGXOptions& options, const MPI_Comm comm)
1383 {
1384  auto amgx = std::make_unique<mfem::AmgXSolver>();
1385  conduit::Node options_node;
1386  options_node["config_version"] = 2;
1387  auto& solver_options = options_node["solver"];
1388  solver_options["solver"] = "AMG";
1389  solver_options["presweeps"] = 1;
1390  solver_options["postsweeps"] = 2;
1391  solver_options["interpolator"] = "D2";
1392  solver_options["max_iters"] = 2;
1393  solver_options["convergence"] = "ABSOLUTE";
1394  solver_options["cycle"] = "V";
1395 
1396  if (options.verbose) {
1397  options_node["solver/obtain_timings"] = 1;
1398  options_node["solver/monitor_residual"] = 1;
1399  options_node["solver/print_solve_stats"] = 1;
1400  }
1401 
1402  // TODO: Use magic_enum here when we can switch to GCC 9+
1403  // This is an immediately-invoked lambda so that the map
1404  // can be const without needed to initialize all the values
1405  // in the constructor
1406  static const auto solver_names = []() {
1407  std::unordered_map<AMGXSolver, std::string> names;
1408  names[AMGXSolver::AMG] = "AMG";
1409  names[AMGXSolver::PCGF] = "PCGF";
1410  names[AMGXSolver::CG] = "CG";
1411  names[AMGXSolver::PCG] = "PCG";
1412  names[AMGXSolver::PBICGSTAB] = "PBICGSTAB";
1413  names[AMGXSolver::BICGSTAB] = "BICGSTAB";
1414  names[AMGXSolver::FGMRES] = "FGMRES";
1415  names[AMGXSolver::JACOBI_L1] = "JACOBI_L1";
1416  names[AMGXSolver::GS] = "GS";
1417  names[AMGXSolver::POLYNOMIAL] = "POLYNOMIAL";
1418  names[AMGXSolver::KPZ_POLYNOMIAL] = "KPZ_POLYNOMIAL";
1419  names[AMGXSolver::BLOCK_JACOBI] = "BLOCK_JACOBI";
1420  names[AMGXSolver::MULTICOLOR_GS] = "MULTICOLOR_GS";
1421  names[AMGXSolver::MULTICOLOR_DILU] = "MULTICOLOR_DILU";
1422  return names;
1423  }();
1424 
1425  options_node["solver/solver"] = solver_names.at(options.solver);
1426  options_node["solver/smoother"] = solver_names.at(options.smoother);
1427 
1428  // Treat the string as the config (not a filename)
1429  amgx->ReadParameters(options_node.to_json(), mfem::AmgXSolver::INTERNAL);
1430  amgx->InitExclusiveGPU(comm);
1431 
1432  return amgx;
1433 }
1434 #endif
1435 
1436 std::unique_ptr<mfem::Solver> buildPreconditioner(LinearSolverOptions linear_opts, [[maybe_unused]] MPI_Comm comm)
1437 {
1438  std::unique_ptr<mfem::Solver> preconditioner_solver;
1439  auto preconditioner = linear_opts.preconditioner;
1440  auto preconditioner_print_level = linear_opts.preconditioner_print_level;
1441 
1442  // Handle the preconditioner - currently just BoomerAMG and HypreSmoother are supported
1443  if (preconditioner == Preconditioner::HypreAMG) {
1444  auto amg_preconditioner = std::make_unique<mfem::HypreBoomerAMG>();
1445  amg_preconditioner->SetPrintLevel(preconditioner_print_level);
1446  preconditioner_solver = std::move(amg_preconditioner);
1447  } else if (preconditioner == Preconditioner::HypreJacobi) {
1448  auto jac_preconditioner = std::make_unique<mfem::HypreSmoother>();
1449  jac_preconditioner->SetType(mfem::HypreSmoother::Type::Jacobi);
1450  preconditioner_solver = std::move(jac_preconditioner);
1451  } else if (preconditioner == Preconditioner::HypreL1Jacobi) {
1452  auto jacl1_preconditioner = std::make_unique<mfem::HypreSmoother>();
1453  jacl1_preconditioner->SetType(mfem::HypreSmoother::Type::l1Jacobi);
1454  preconditioner_solver = std::move(jacl1_preconditioner);
1455  } else if (preconditioner == Preconditioner::HypreGaussSeidel) {
1456  auto gs_preconditioner = std::make_unique<mfem::HypreSmoother>();
1457  gs_preconditioner->SetType(mfem::HypreSmoother::Type::GS);
1458  preconditioner_solver = std::move(gs_preconditioner);
1459  } else if (preconditioner == Preconditioner::HypreILU) {
1460  auto ilu_preconditioner = std::make_unique<mfem::HypreILU>();
1461  ilu_preconditioner->SetLevelOfFill(1);
1462  ilu_preconditioner->SetPrintLevel(preconditioner_print_level);
1463  preconditioner_solver = std::move(ilu_preconditioner);
1464  } else if (preconditioner == Preconditioner::AMGX) {
1465 #ifdef MFEM_USE_AMGX
1466  preconditioner_solver = buildAMGX(linear_opts.amgx_options, comm);
1467 #else
1468  SLIC_ERROR_ROOT("AMGX requested in non-GPU build");
1469 #endif
1470  } else if (preconditioner == Preconditioner::Petsc) {
1471 #ifdef SMITH_USE_PETSC
1472  preconditioner_solver = mfem_ext::buildPetscPreconditioner(linear_opts.petsc_preconditioner, comm);
1473 #else
1474  SLIC_ERROR_ROOT("PETSc preconditioner requested in non-PETSc build");
1475 #endif
1476  } else if (preconditioner == Preconditioner::AMGFContact) {
1477  auto amgfcontact_preconditioner = std::make_unique<mfem::AMGFSolver>();
1478  auto amgfcontact_opts = linear_opts.amgfcontact_options;
1479  amgfcontact_preconditioner->GetAMG().SetPrintLevel(preconditioner_print_level);
1480  amgfcontact_preconditioner->GetAMG().SetSystemsOptions(amgfcontact_opts.dim_systems_options);
1481  amgfcontact_preconditioner->GetAMG().SetRelaxType(amgfcontact_opts.relax_type);
1482  preconditioner_solver = std::move(amgfcontact_preconditioner);
1483  } else if (preconditioner == Preconditioner::BlockDiagonal || preconditioner == Preconditioner::BlockTriangular ||
1484  preconditioner == Preconditioner::BlockSchur) {
1485  std::vector<std::unique_ptr<mfem::Solver>> inner_solvers;
1486  for (const auto& opt : linear_opts.sub_block_linear_solver_options) {
1487  auto [lin, prec] = buildLinearSolverAndPreconditioner(opt, comm);
1488  inner_solvers.push_back(std::make_unique<SolverWithPreconditioner>(std::move(lin), std::move(prec)));
1489  }
1490 
1491  if (preconditioner == Preconditioner::BlockDiagonal) {
1492  preconditioner_solver = std::make_unique<BlockDiagonalPreconditioner>(std::move(inner_solvers));
1493  } else if (preconditioner == Preconditioner::BlockTriangular) {
1494  preconditioner_solver =
1495  std::make_unique<BlockTriangularPreconditioner>(std::move(inner_solvers), linear_opts.block_triangular_type);
1496  } else if (preconditioner == Preconditioner::BlockSchur) {
1497  preconditioner_solver = std::make_unique<BlockSchurPreconditioner>(
1498  std::move(inner_solvers), linear_opts.block_schur_type, linear_opts.schur_approx_type);
1499  }
1500  } else {
1501  SLIC_ERROR_ROOT_IF(preconditioner != Preconditioner::None, "Unknown preconditioner type requested");
1502  }
1503 
1504  return preconditioner_solver;
1505 }
1506 
1507 void EquationSolver::defineInputFileSchema(axom::inlet::Container& container)
1508 {
1509  auto& linear_container = container.addStruct("linear", "Linear Equation Solver Parameters");
1510  linear_container.required().registerVerifier([](const axom::inlet::Container& container_to_verify) {
1511  // Make sure that the provided options match the desired linear solver type
1512  const bool is_iterative = (container_to_verify["type"].get<std::string>() == "iterative") &&
1513  container_to_verify.contains("iterative_options");
1514  const bool is_direct =
1515  (container_to_verify["type"].get<std::string>() == "direct") && container_to_verify.contains("direct_options");
1516  return is_iterative || is_direct;
1517  });
1518 
1519  // Enforce the solver type - must be iterative or direct
1520  linear_container.addString("type", "The type of solver parameters to use (iterative|direct)")
1521  .required()
1522  .validValues({"iterative", "direct"});
1523 
1524  auto& iterative_container = linear_container.addStruct("iterative_options", "Iterative solver parameters");
1525  iterative_container.addDouble("rel_tol", "Relative tolerance for the linear solve.").defaultValue(1.0e-6);
1526  iterative_container.addDouble("abs_tol", "Absolute tolerance for the linear solve.").defaultValue(1.0e-8);
1527  iterative_container.addInt("max_iter", "Maximum iterations for the linear solve.").defaultValue(5000);
1528  iterative_container.addInt("print_level", "Linear print level.").defaultValue(0);
1529  iterative_container.addString("solver_type", "Solver type (gmres|minres|cg).").defaultValue("gmres");
1530  iterative_container.addString("prec_type", "Preconditioner type (JacobiSmoother|L1JacobiSmoother|AMG|ILU|Petsc).")
1531  .defaultValue("JacobiSmoother");
1532  iterative_container.addString("petsc_prec_type", "Type of PETSc preconditioner to use.").defaultValue("jacobi");
1533 
1534  auto& direct_container = linear_container.addStruct("direct_options", "Direct solver parameters");
1535  direct_container.addInt("print_level", "Linear print level.").defaultValue(0);
1536 
1537  // Only needed for nonlinear problems
1538  auto& nonlinear_container = container.addStruct("nonlinear", "Newton Equation Solver Parameters").required(false);
1539  nonlinear_container.addDouble("rel_tol", "Relative tolerance for the Newton solve.").defaultValue(1.0e-2);
1540  nonlinear_container.addDouble("abs_tol", "Absolute tolerance for the Newton solve.").defaultValue(1.0e-4);
1541  nonlinear_container.addInt("max_iter", "Maximum iterations for the Newton solve.").defaultValue(500);
1542  nonlinear_container.addInt("print_level", "Nonlinear print level.").defaultValue(0);
1543  nonlinear_container.addString("solver_type", "Solver type (Newton|KINFullStep|KINLineSearch)").defaultValue("Newton");
1544 }
1545 
1546 } // namespace smith
1547 
1548 using smith::EquationSolver;
1551 
1553 {
1554  LinearSolverOptions options;
1555  std::string type = base["type"];
1556 
1557  if (type == "direct") {
1559  options.print_level = base["direct_options/print_level"];
1560  return options;
1561  }
1562 
1563  auto config = base["iterative_options"];
1564  options.relative_tol = config["rel_tol"];
1565  options.absolute_tol = config["abs_tol"];
1566  options.max_iterations = config["max_iter"];
1567  options.print_level = config["print_level"];
1568  std::string solver_type = config["solver_type"];
1569  if (solver_type == "gmres") {
1571  } else if (solver_type == "cg") {
1573  } else {
1574  std::string msg = std::format("Unknown Linear solver type given: '{0}'", solver_type);
1575  SLIC_ERROR_ROOT(msg);
1576  }
1577  const std::string prec_type = config["prec_type"];
1578  if (prec_type == "JacobiSmoother") {
1580  } else if (prec_type == "L1JacobiSmoother") {
1582  } else if (prec_type == "HypreAMG") {
1584  } else if (prec_type == "ILU") {
1586 #ifdef MFEM_USE_AMGX
1587  } else if (prec_type == "AMGX") {
1589 #endif
1590  } else if (prec_type == "GaussSeidel") {
1592 #ifdef SMITH_USE_PETSC
1593  } else if (prec_type == "Petsc") {
1594  const std::string petsc_prec = config["petsc_prec_type"];
1596  options.petsc_preconditioner = smith::mfem_ext::stringToPetscPCType(petsc_prec);
1597 #endif
1598  } else if (prec_type == "AMGFContact") {
1600  } else {
1601  std::string msg = std::format("Unknown preconditioner type given: '{0}'", prec_type);
1602  SLIC_ERROR_ROOT(msg);
1603  }
1604 
1605  return options;
1606 }
1607 
1609 {
1610  NonlinearSolverOptions options;
1611  options.relative_tol = base["rel_tol"];
1612  options.absolute_tol = base["abs_tol"];
1613  options.max_iterations = base["max_iter"];
1614  options.print_level = base["print_level"];
1615  const std::string solver_type = base["solver_type"];
1616  if (solver_type == "Newton") {
1618  } else if (solver_type == "KINFullStep") {
1620  } else if (solver_type == "KINLineSearch") {
1622  } else if (solver_type == "KINPicard") {
1624  } else {
1625  SLIC_ERROR_ROOT(std::format("Unknown nonlinear solver type given: '{0}'", solver_type));
1626  }
1627  return options;
1628 }
1629 
1631 {
1632  auto lin = base["linear"].get<LinearSolverOptions>();
1633  auto nonlin = base["nonlinear"].get<NonlinearSolverOptions>();
1634 
1635  auto [linear_solver, preconditioner] = smith::buildLinearSolverAndPreconditioner(lin, MPI_COMM_WORLD);
1636 
1637  smith::EquationSolver eq_solver(smith::buildNonlinearSolver(nonlin, lin, *preconditioner, MPI_COMM_WORLD),
1638  std::move(linear_solver), std::move(preconditioner));
1639 
1640  return eq_solver;
1641 }
This class manages the objects typically required to solve a nonlinear set of equations arising from ...
void resetConvergenceState() const
Reset nonlinear convergence state stored across iterations of a single solve.
mfem::Solver & preconditioner()
void setConvergenceTolerances(double abs_tol, double rel_tol, MPI_Comm comm) const
Update scalar convergence tolerances for the managed nonlinear convergence path.
EquationSolver(std::unique_ptr< mfem::NewtonSolver > nonlinear_solver, std::unique_ptr< mfem::Solver > linear_solver, std::unique_ptr< mfem::Solver > preconditioner=nullptr)
void solve(mfem::Vector &x) const
static void defineInputFileSchema(axom::inlet::Container &container)
mfem::Solver & linearSolver()
void setOperator(const mfem::Operator &op)
void SetOperator(const mfem::Operator &op)
Set the underlying matrix operator to use in the solution algorithm.
void Mult(const mfem::Vector &input, mfem::Vector &output) const
Factor and solve the linear system y = Op^{-1} x using DSuperLU.
This file contains the declaration of an equation solver wrapper.
This file contains the all the necessary functions and macros required for logging as well as a helpe...
Accelerator functionality.
Definition: smith.cpp:36
std::unique_ptr< mfem::NewtonSolver > buildNonlinearSolver(NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions &linear_opts, mfem::Solver &prec, MPI_Comm comm)
Build a nonlinear solver using the nonlinear option struct.
bool requiresMonolithicOperator(const LinearSolverOptions &linear_opts)
Return true if the configured linear solve stack requires block operators to be merged.
std::unique_ptr< mfem::HypreParMatrix > buildMonolithicMatrix(const mfem::BlockOperator &block_operator)
Build a monolithic HypreParMatrix from a BlockOperator.
std::pair< std::unique_ptr< mfem::Solver >, std::unique_ptr< mfem::Solver > > buildLinearSolverAndPreconditioner(LinearSolverOptions linear_opts, MPI_Comm comm)
Build the linear solver and its associated preconditioner given a linear options struct.
std::string linearName(const LinearSolver &s)
Convert linear solver enums to their string names.
LinearSolver
Linear solution method indicator.
Preconditioner
The type of preconditioner to be used.
SMITH_HOST_DEVICE auto sqrt(dual< gradient_type > x)
implementation of square root for dual numbers
Definition: dual.hpp:308
constexpr SMITH_HOST_DEVICE auto norm(const isotropic_tensor< T, m, m > &I)
compute the Frobenius norm (sqrt(tr(dot(transpose(I), I)))) of an isotropic tensor
SMITH_HOST_DEVICE auto max(dual< gradient_type > a, double b)
Implementation of max for dual numbers.
Definition: dual.hpp:229
constexpr SMITH_HOST_DEVICE int size(const tensor< T, n... > &)
returns the total number of stored values in a tensor
Definition: tensor.hpp:1939
std::unique_ptr< mfem::Solver > buildPreconditioner(LinearSolverOptions linear_opts, [[maybe_unused]] MPI_Comm comm)
Build a preconditioner from the available options.
SMITH_HOST_DEVICE auto abs(dual< gradient_type > x)
Implementation of absolute value function for dual numbers.
Definition: dual.hpp:219
std::string preconditionerName(Preconditioner p)
Convert preconditioner enums to their string names.
Various helper functions and macros for profiling using Caliper.
#define SMITH_MARK_FUNCTION
Definition: profiling.hpp:90
smith::EquationSolver operator()(const axom::inlet::Container &base)
Returns created object from Inlet container.
smith::LinearSolverOptions operator()(const axom::inlet::Container &base)
Returns created object from Inlet container.
smith::NonlinearSolverOptions operator()(const axom::inlet::Container &base)
Returns created object from Inlet container.
Parameters for an iterative linear solution scheme.
int preconditioner_print_level
Debugging print level for the preconditioner.
int max_iterations
Maximum number of iterations.
double relative_tol
Relative tolerance.
int print_level
Debugging print level for the linear solver.
AMGFContactOptions amgfcontact_options
AMGFContact Options, used for Preconditioner::AMGFContact.
Preconditioner preconditioner
PreconditionerOptions selection.
SchurApproxType schur_approx_type
Schur approximation type.
PetscPCType petsc_preconditioner
PETSc preconditioner type.
std::vector< LinearSolverOptions > sub_block_linear_solver_options
Subblock linear solver options for block preconditioners.
BlockSchurType block_schur_type
Block Schur preconditioner factorization type.
double absolute_tol
Absolute tolerance.
AMGXOptions amgx_options
AMGX Options, used for Preconditioner::AMGX.
BlockTriangularType block_triangular_type
Block Triangular Preconditioner factorization type.
LinearSolver linear_solver
Linear solver selection.
Nonlinear solution scheme parameters.
int min_iterations
Minimum number of iterations.
int max_iterations
Maximum number of iterations.
double trust_region_scaling
Scaling for the initial trust region size.
int max_line_search_iterations
Maximum line search cutbacks.
double relative_tol
Relative tolerance.
int print_level
Debug print level.
int num_leftmost
Number of extra leftmost eigenvector to be stored between solves.
double absolute_tol
Absolute tolerance.
NonlinearSolver nonlin_solver
Nonlinear solver selection.
SubSpaceOptions subspace_option
Option for how when the subspace solver should be utilized within trust-region solver.
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:122