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 
9 #include <cstdlib>
10 #include <iomanip>
11 #include <iostream>
12 #include <algorithm>
13 #include <cmath>
14 #include <exception>
15 #include <limits>
16 #include <string>
17 #include <tuple>
18 
19 #include "smith/smith_config.hpp"
21 #include "smith/numerics/trust_region_solver.hpp"
23 
24 namespace smith {
25 
27 class NewtonSolver : public mfem::NewtonSolver {
28  protected:
30  mutable mfem::Vector x0;
31 
34 
36  mutable size_t print_level = 0;
37 
38  public:
40  NewtonSolver(const NonlinearSolverOptions& nonlinear_opts) : nonlinear_options(nonlinear_opts) {}
41 
42 #ifdef MFEM_USE_MPI
44  NewtonSolver(MPI_Comm comm_, const NonlinearSolverOptions& nonlinear_opts)
45  : mfem::NewtonSolver(comm_), nonlinear_options(nonlinear_opts)
46  {
47  }
48 #endif
49 
51  double evaluateNorm(const mfem::Vector& x, mfem::Vector& rOut) const
52  {
54  double normEval = std::numeric_limits<double>::max();
55  try {
56  oper->Mult(x, rOut);
57  normEval = Norm(rOut);
58  } catch (const std::exception&) {
60  }
61  return normEval;
62  }
63 
65  void assembleJacobian(const mfem::Vector& x) const
66  {
68  grad = &oper->GetGradient(x);
70  auto* grad_blocked = dynamic_cast<mfem::BlockOperator*>(grad);
71  if (grad_blocked) grad = buildMonolithicMatrix(*grad_blocked).release();
72  }
73  }
74 
76  void setPreconditioner() const
77  {
79  prec->SetOperator(*grad);
80  }
81 
83  void solveLinearSystem(const mfem::Vector& r_, mfem::Vector& c_) const
84  {
86  prec->Mult(r_, c_); // c = [DF(x_i)]^{-1} [F(x_i)-b]
87  }
88 
90  void Mult(const mfem::Vector&, mfem::Vector& x) const
91  {
92  MFEM_ASSERT(oper != NULL, "the Operator is not set (use SetOperator).");
93  MFEM_ASSERT(prec != NULL, "the Solver is not set (use SetSolver).");
94 
95  print_level = print_options.iterations ? 1 : print_level;
96  print_level = print_options.summary ? 2 : print_level;
97 
98  using real_t = mfem::real_t;
99 
100  real_t norm, norm_goal = 0;
101  norm = initial_norm = evaluateNorm(x, r);
102  if (norm == 0.0) return;
103 
104  if (print_level == 1) {
105  mfem::out << "Newton iteration " << std::setw(3) << 0 << " : ||r|| = " << std::setw(13) << norm << "\n";
106  }
107 
108  norm_goal = std::max(rel_tol * initial_norm, abs_tol);
109  prec->iterative_mode = false;
110 
111  int it = 0;
112  for (; true; it++) {
113  MFEM_ASSERT(mfem::IsFinite(norm), "norm = " << norm);
114  if (print_level >= 2) {
115  mfem::out << "Newton iteration " << std::setw(3) << it << " : ||r|| = " << std::setw(13) << norm;
116  if (it > 0) {
117  mfem::out << ", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ? norm / initial_norm : norm);
118  }
119  mfem::out << '\n';
120  }
121 
122  if ((print_level >= 1) && (norm != norm)) {
123  mfem::out << "Initial residual for Newton iteration is undefined/nan.\n";
124  mfem::out << "Newton: No convergence!\n";
125  return;
126  }
127 
128  if (norm <= norm_goal && it >= nonlinear_options.min_iterations) {
129  converged = true;
130  break;
131  } else if (it >= max_iter) {
132  converged = false;
133  break;
134  }
135 
136  real_t norm_nm1 = norm;
137 
138  assembleJacobian(x);
140  solveLinearSystem(r, c);
141 
142  // there must be a better way to do this?
143  x0.SetSize(x.Size());
144  x0 = 0.0;
145  x0.Add(1.0, x);
146 
147  real_t stepScale = 1.0;
148  add(x0, -stepScale, c, x);
149  norm = evaluateNorm(x, r);
150 
151  const int max_ls_iters = nonlinear_options.max_line_search_iterations;
152  static constexpr real_t reduction = 0.5;
153 
154  const double sufficientDecreaseParam = 0.0; // 1e-15;
155  const double cMagnitudeInR = sufficientDecreaseParam != 0.0 ? std::abs(Dot(c, r)) / norm_nm1 : 0.0;
156 
157  auto is_improved = [=](real_t currentNorm, real_t c_scale) {
158  return currentNorm < norm_nm1 - sufficientDecreaseParam * c_scale * cMagnitudeInR;
159  };
160 
161  // back-track linesearch
162  int ls_iter = 0;
163  int ls_iter_sum = 0;
164  for (; !is_improved(norm, stepScale) && ls_iter < max_ls_iters; ++ls_iter, ++ls_iter_sum) {
165  stepScale *= reduction;
166  add(x0, -stepScale, c, x);
167  norm = evaluateNorm(x, r);
168  }
169 
170  // try the opposite direction and linesearch back from there
171  if (max_ls_iters > 0 && ls_iter == max_ls_iters && !is_improved(norm, stepScale)) {
172  stepScale = 1.0;
173  add(x0, stepScale, c, x);
174  norm = evaluateNorm(x, r);
175 
176  ls_iter = 0;
177  for (; !is_improved(norm, stepScale) && ls_iter < max_ls_iters; ++ls_iter, ++ls_iter_sum) {
178  stepScale *= reduction;
179  add(x0, stepScale, c, x);
180  norm = evaluateNorm(x, r);
181  }
182 
183  // ok, the opposite direction was also terrible, lets go back, cut in half 1 last time and accept it hoping for
184  // the best
185  if (ls_iter == max_ls_iters && !is_improved(norm, stepScale)) {
186  ++ls_iter_sum;
187  stepScale *= reduction;
188  add(x0, -stepScale, c, x);
189  norm = evaluateNorm(x, r);
190  }
191  }
192 
193  if (ls_iter_sum) {
194  if (print_level >= 2) {
195  mfem::out << "Number of line search steps taken = " << ls_iter_sum << std::endl;
196  }
197  if (print_level >= 2 && (ls_iter_sum == 2 * max_ls_iters + 1)) {
198  mfem::out << "The maximum number of line search cut back have occurred, the resulting residual may not have "
199  "decreased. "
200  << std::endl;
201  }
202  }
203  }
204 
205  final_iter = it;
206  final_norm = norm;
207 
208  if (print_level == 1) {
209  mfem::out << "Newton iteration " << std::setw(3) << final_iter << " : ||r|| = " << std::setw(13) << norm << '\n';
210  }
211  if (!converged && print_level >= 1) { // (print_options.summary || print_options.warnings)) {
212  mfem::out << "Newton: No convergence!\n";
213  }
214  }
215 };
216 
220  double cg_tol = 1e-8;
222  size_t min_cg_iterations = 0; //
224  size_t max_cg_iterations = 10000; //
228  double min_tr_size = 1e-13;
230  double t1 = 0.25;
232  double t2 = 1.75;
234  double eta1 = 1e-9;
236  double eta2 = 0.1;
238  double eta3 = 0.6;
240  double eta4 = 4.2;
241 };
242 
247  {
248  z.SetSize(size);
249  H_z.SetSize(size);
250  d_old.SetSize(size);
251  H_d_old.SetSize(size);
252  d.SetSize(size);
253  H_d.SetSize(size);
254  Pr.SetSize(size);
255  cauchy_point.SetSize(size);
256  H_cauchy_point.SetSize(size);
257  }
258 
260  void reset()
261  {
262  z = 0.0;
263  cauchy_point = 0.0;
264  }
265 
267  enum class Status
268  {
269  Interior,
270  NegativeCurvature,
271  OnBoundary,
272  NonDescentDirection
273  };
274 
276  mfem::Vector z;
278  mfem::Vector H_z;
280  mfem::Vector d_old;
282  mfem::Vector H_d_old;
284  mfem::Vector d;
286  mfem::Vector H_d;
288  mfem::Vector Pr;
290  mfem::Vector cauchy_point;
292  mfem::Vector H_cauchy_point;
294  Status interior_status = Status::Interior;
297 };
298 
300 void printTrustRegionInfo(double realObjective, double modelObjective, size_t cgIters, double trSize, bool willAccept)
301 {
302  mfem::out << "real energy = " << std::setw(13) << realObjective << ", model energy = " << std::setw(13)
303  << modelObjective << ", cg iter = " << std::setw(7) << cgIters << ", next tr size = " << std::setw(8)
304  << trSize << ", accepting = " << willAccept << std::endl;
305 }
306 
316 class TrustRegion : public mfem::NewtonSolver {
317  protected:
319  mutable mfem::Vector x_pred;
321  mutable mfem::Vector r_pred;
323  mutable mfem::Vector scratch;
325  mutable std::vector<std::shared_ptr<mfem::Vector>> left_mosts;
327  mutable std::vector<std::shared_ptr<mfem::Vector>> H_left_mosts;
328 
333 
336  Solver& tr_precond;
337 
339  mutable size_t print_level = 0;
340 
341  public:
343  mutable size_t num_hess_vecs = 0;
345  mutable size_t num_preconds = 0;
347  mutable size_t num_residuals = 0;
349  mutable size_t num_subspace_solves = 0;
351  mutable size_t num_jacobian_assembles = 0;
352 
353 #ifdef MFEM_USE_MPI
355  TrustRegion(MPI_Comm comm_, const NonlinearSolverOptions& nonlinear_opts, const LinearSolverOptions& linear_opts,
356  Solver& tPrec)
357  : mfem::NewtonSolver(comm_), nonlinear_options(nonlinear_opts), linear_options(linear_opts), tr_precond(tPrec)
358  {
359  }
360 #endif
361 
363  void projectToBoundaryWithCoefs(mfem::Vector& z, const mfem::Vector& d, double delta, double zz, double zd,
364  double dd) const
365  {
366  // find z + tau d
367  double deltadelta_m_zz = delta * delta - zz;
368  if (deltadelta_m_zz == 0) return; // already on boundary
369  double tau = (std::sqrt(deltadelta_m_zz * dd + zd * zd) - zd) / dd;
370  z.Add(tau, d);
371  }
372 
374  template <typename HessVecFunc>
375  void solveTheSubspaceProblem([[maybe_unused]] mfem::Vector& z, [[maybe_unused]] const HessVecFunc& hess_vec_func,
376  [[maybe_unused]] const std::vector<const mfem::Vector*> ds,
377  [[maybe_unused]] const std::vector<const mfem::Vector*> Hds,
378  [[maybe_unused]] const mfem::Vector& g, [[maybe_unused]] double delta,
379  [[maybe_unused]] int num_leftmost) const
380  {
381 #ifdef SMITH_USE_SLEPC
384 
385  std::vector<const mfem::Vector*> directions;
386  for (auto& d : ds) {
387  directions.emplace_back(d);
388  }
389  for (auto& left : left_mosts) {
390  directions.emplace_back(left.get());
391  }
392 
393  std::vector<const mfem::Vector*> H_directions;
394  for (auto& Hd : Hds) {
395  H_directions.emplace_back(Hd);
396  }
397  for (auto& H_left : H_left_mosts) {
398  H_directions.emplace_back(H_left.get());
399  }
400 
401  try {
402  std::tie(directions, H_directions) = removeDependentDirections(directions, H_directions);
403  } catch (const std::exception& e) {
404  if (print_level >= 2) {
405  mfem::out << "remove dependent directions failed with " << e.what() << std::endl;
406  }
407  return;
408  }
409 
410  mfem::Vector b(g);
411  b *= -1;
412 
413  mfem::Vector sol;
414  std::vector<std::shared_ptr<mfem::Vector>> leftvecs;
415  std::vector<double> leftvals;
416  double energy_change;
417 
418  try {
419  std::tie(sol, leftvecs, leftvals, energy_change) =
420  solveSubspaceProblem(directions, H_directions, b, delta, num_leftmost);
421  } catch (const std::exception& e) {
422  if (print_level == 1) {
423  mfem::out << "subspace solve failed with " << e.what() << std::endl;
424  }
425  return;
426  }
427 
428  left_mosts.clear();
429  for (auto& lv : leftvecs) {
430  left_mosts.emplace_back(std::move(lv));
431  }
432 
433  double base_energy = computeEnergy(g, hess_vec_func, z);
434  double subspace_energy = computeEnergy(g, hess_vec_func, sol);
435 
436  if (print_level >= 2) {
437  double leftval = leftvals.size() ? leftvals[0] : 1.0;
438  mfem::out << "Energy using subspace solver from: " << base_energy << ", to: " << subspace_energy << " / "
439  << energy_change << ". Min eig: " << leftval << std::endl;
440  }
441 
442  if (subspace_energy < base_energy) {
443  z = sol;
444  }
445 #endif
446  }
447 
449  void projectToBoundaryBetweenWithCoefs(mfem::Vector& z, const mfem::Vector& y, double trSize, double zz, double zy,
450  double yy) const
451  {
452  double dd = yy - 2 * zy + zz;
453  double zd = zy - zz;
454  double tau = (std::sqrt((trSize * trSize - zz) * dd + zd * zd) - zd) / dd;
455  z.Add(-tau, z);
456  z.Add(tau, y);
457  }
458 
460  void doglegStep(const mfem::Vector& cp, const mfem::Vector& newtonP, double trSize, mfem::Vector& s) const
461  {
463  // MRT, could optimize some of these eventually, compute on the outside and save
464  double cc = Dot(cp, cp);
465  double nn = Dot(newtonP, newtonP);
466  double tt = trSize * trSize;
467 
468  s = 0.0;
469  if (cc >= tt) {
470  add(s, std::sqrt(tt / cc), cp, s);
471  } else if (cc > nn) {
472  if (print_level >= 2) {
473  mfem::out << "cp outside newton, preconditioner likely inaccurate\n";
474  }
475  add(s, 1.0, cp, s);
476  } else if (nn > tt) { // on the dogleg (we have nn >= cc, and tt >= cc)
477  add(s, 1.0, cp, s);
478  double cn = Dot(cp, newtonP);
479  projectToBoundaryBetweenWithCoefs(s, newtonP, trSize, cc, cn, nn);
480  } else {
481  s = newtonP;
482  }
483  }
484 
486  template <typename HessVecFunc>
487  double computeEnergy(const mfem::Vector& r_local, const HessVecFunc& H, const mfem::Vector& z) const
488  {
490  double rz = Dot(r_local, z);
491  mfem::Vector tmp(r_local);
492  tmp = 0.0;
493  H(z, tmp);
494  return rz + 0.5 * Dot(z, tmp);
495  }
496 
498  template <typename HessVecFunc, typename PrecondFunc>
499  void solveTrustRegionModelProblem(const mfem::Vector& r0, mfem::Vector& rCurrent, HessVecFunc hess_vec_func,
500  PrecondFunc precond, const TrustRegionSettings& settings, double& trSize,
501  TrustRegionResults& results) const
502  {
504  // minimize r0@z + 0.5*z@J@z
505  results.interior_status = TrustRegionResults::Status::Interior;
506  results.cg_iterations_count = 0;
507 
508  auto& z = results.z;
509  auto& cgIter = results.cg_iterations_count;
510  auto& d = results.d;
511  auto& Pr = results.Pr;
512  auto& Hd = results.H_d;
513 
514  const double cg_tol_squared = settings.cg_tol * settings.cg_tol;
515 
516  if (Dot(r0, r0) <= cg_tol_squared && settings.min_cg_iterations == 0) {
517  if (print_level >= 2) {
518  mfem::out << "Trust region solution state within tolerance on first iteration."
519  << "\n";
520  }
521  return;
522  }
523 
524  rCurrent = r0;
525  precond(rCurrent, Pr);
526 
527  // d = -Pr
528  d = Pr;
529  d *= -1.0;
530 
531  z = 0.0;
532  double zz = 0.;
533  double rPr = Dot(rCurrent, Pr);
534  double zd = 0.0;
535  double dd = Dot(d, d);
536 
537  // std::cout << "initial energy = " << computeEnergy(r0, hess_vec_func, z) << std::endl;
538 
539  for (cgIter = 1; cgIter <= settings.max_cg_iterations; ++cgIter) {
540  // check if this is a descent direction
541  if (Dot(d, rCurrent) > 0) {
542  d *= -1;
543  results.interior_status = TrustRegionResults::Status::NonDescentDirection;
544  }
545 
546  hess_vec_func(d, Hd);
547  const double curvature = Dot(d, Hd);
548  const double alphaCg = curvature != 0.0 ? rPr / curvature : 0.0;
549 
550  auto& zPred = Pr; // re-use Pr memory.
551  // This predicted step will no longer be used by the time Pr is, so we can avoid an extra
552  // vector floating around
553  add(z, alphaCg, d, zPred);
554  double zzNp1 = Dot(zPred, zPred);
555 
556  const bool go_to_boundary = curvature <= 0 || zzNp1 >= trSize * trSize;
557  if (go_to_boundary) {
558  projectToBoundaryWithCoefs(z, d, trSize, zz, zd, dd);
559  if (curvature <= 0) {
560  results.interior_status = TrustRegionResults::Status::NegativeCurvature;
561  } else {
562  results.interior_status = TrustRegionResults::Status::OnBoundary;
563  }
564  return;
565  }
566 
567  z = zPred;
568 
569  if (results.interior_status == TrustRegionResults::Status::NonDescentDirection) {
570  if (print_level >= 2) {
571  mfem::out << "Found a non descent direction\n";
572  }
573  return;
574  }
575 
576  add(rCurrent, alphaCg, Hd, rCurrent);
577 
578  precond(rCurrent, Pr);
579  double rPrNp1 = Dot(rCurrent, Pr);
580 
581  if (Dot(rCurrent, rCurrent) <= cg_tol_squared && cgIter >= settings.min_cg_iterations) {
582  return;
583  }
584 
585  double beta = rPrNp1 / rPr;
586  rPr = rPrNp1;
587  add(-1.0, Pr, beta, d, d);
588 
589  zz = zzNp1;
590  zd = Dot(z, d);
591  dd = Dot(d, d);
592  }
593  cgIter--; // if all cg iterations are taken, correct for output
594  }
595 
597  void assembleJacobian(const mfem::Vector& x) const
598  {
601  grad = &oper->GetGradient(x);
603  auto* grad_blocked = dynamic_cast<mfem::BlockOperator*>(grad);
604  if (grad_blocked) grad = buildMonolithicMatrix(*grad_blocked).release();
605  }
606  }
607 
609  mfem::real_t computeResidual(const mfem::Vector& x_, mfem::Vector& r_) const
610  {
612  ++num_residuals;
613  oper->Mult(x_, r_);
614  return Norm(r_);
615  }
616 
618  void hessVec(const mfem::Vector& x_, mfem::Vector& v_) const
619  {
621  ++num_hess_vecs;
622  grad->Mult(x_, v_);
623  }
624 
626  void precond(const mfem::Vector& x_, mfem::Vector& v_) const
627  {
629  ++num_preconds;
630  tr_precond.Mult(x_, v_);
631  };
632 
634  void Mult(const mfem::Vector&, mfem::Vector& X) const
635  {
636  MFEM_ASSERT(oper != NULL, "the Operator is not set (use SetOperator).");
637  MFEM_ASSERT(prec != NULL, "the Solver is not set (use SetSolver).");
638 
639  print_level = print_options.iterations ? 1 : print_level;
640  print_level = print_options.summary ? 2 : print_level;
641 
642  using real_t = mfem::real_t;
643 
644  num_hess_vecs = 0;
645  num_preconds = 0;
646  num_residuals = 0;
649 
650  real_t norm, norm_goal = 0.0;
651  norm = initial_norm = computeResidual(X, r);
652  if (norm == 0.0) return;
653 
654  norm_goal = std::max(rel_tol * initial_norm, abs_tol);
655 
656  if (print_level == 1) {
657  mfem::out << "TrustRegion iteration " << std::setw(3) << 0 << " : ||r|| = " << std::setw(13) << norm << "\n";
658  }
659 
660  prec->iterative_mode = false;
661  tr_precond.iterative_mode = false;
662 
663  // local arrays
664  x_pred.SetSize(X.Size());
665  x_pred = 0.0;
666  r_pred.SetSize(X.Size());
667  r_pred = 0.0;
668  scratch.SetSize(X.Size());
669  scratch = 0.0;
670 
671  TrustRegionResults trResults(X.Size());
672  TrustRegionSettings settings;
673  settings.min_cg_iterations = static_cast<size_t>(nonlinear_options.min_iterations);
674  settings.max_cg_iterations = static_cast<size_t>(linear_options.max_iterations);
675  settings.cg_tol = 0.5 * norm_goal;
676 
677  int subspace_option = nonlinear_options.subspace_option;
678  int num_leftmost = nonlinear_options.num_leftmost;
679 
680  scratch = 1.0;
682  size_t cumulative_cg_iters_from_last_precond_update = 0;
683 
684  int it = 0;
685  for (; true; it++) {
686  MFEM_ASSERT(mfem::IsFinite(norm), "norm = " << norm);
687  if (print_level >= 2) {
688  mfem::out << "TrustRegion iteration " << std::setw(3) << it << " : ||r|| = " << std::setw(13) << norm;
689  if (it > 0) {
690  mfem::out << ", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ? norm / initial_norm : norm);
691  mfem::out << ", x_incr = " << std::setw(13) << trResults.d.Norml2();
692  } else {
693  mfem::out << ", norm goal = " << std::setw(13) << norm_goal;
694  }
695  mfem::out << '\n';
696  }
697 
698  if (print_level >= 1 && (norm != norm)) {
699  mfem::out << "Initial residual for trust-region iteration is undefined/nan." << std::endl;
700  mfem::out << "TrustRegion: No convergence!\n";
701  return;
702  }
703 
704  if (norm <= norm_goal && it >= nonlinear_options.min_iterations) {
705  converged = true;
706  break;
707  } else if (it >= max_iter) {
708  converged = false;
709  break;
710  }
711 
712  assembleJacobian(X);
713 
714  if (it == 0 || (trResults.cg_iterations_count >= settings.max_cg_iterations ||
715  cumulative_cg_iters_from_last_precond_update >= settings.max_cumulative_iteration)) {
716  tr_precond.SetOperator(*grad);
717  cumulative_cg_iters_from_last_precond_update = 0;
718  }
719 
720  auto hess_vec_func = [&](const mfem::Vector& x_, mfem::Vector& v_) { hessVec(x_, v_); };
721  auto precond_func = [&](const mfem::Vector& x_, mfem::Vector& v_) { precond(x_, v_); };
722 
723  double cauchyPointNormSquared = tr_size * tr_size;
724  trResults.reset();
725 
726  hess_vec_func(r, trResults.H_d);
727  const double gKg = Dot(r, trResults.H_d);
728  if (gKg > 0) {
729  const double alphaCp = -Dot(r, r) / gKg;
730  add(trResults.cauchy_point, alphaCp, r, trResults.cauchy_point);
731  cauchyPointNormSquared = Dot(trResults.cauchy_point, trResults.cauchy_point);
732  } else {
733  const double alphaTr = -tr_size / std::sqrt(Dot(r, r));
734  add(trResults.cauchy_point, alphaTr, r, trResults.cauchy_point);
735  if (print_level >= 2) {
736  mfem::out << "Negative curvature un-preconditioned cauchy point direction found."
737  << "\n";
738  }
739  }
740 
741  if (cauchyPointNormSquared >= tr_size * tr_size) {
742  if (print_level >= 2) {
743  mfem::out << "Un-preconditioned gradient cauchy point outside trust region, step size = "
744  << std::sqrt(cauchyPointNormSquared) << "\n";
745  }
746  trResults.cauchy_point *= (tr_size / std::sqrt(cauchyPointNormSquared));
747  trResults.z = trResults.cauchy_point;
748 
749  trResults.cg_iterations_count = 1;
750  trResults.interior_status = TrustRegionResults::Status::OnBoundary;
751  } else {
752  settings.cg_tol = std::max(0.5 * norm_goal, 5e-5 * norm);
753  solveTrustRegionModelProblem(r, scratch, hess_vec_func, precond_func, settings, tr_size, trResults);
754  }
755  cumulative_cg_iters_from_last_precond_update += trResults.cg_iterations_count;
756 
757  bool have_computed_Hvs = false;
758 
759  int lineSearchIter = 0;
760  while (lineSearchIter <= nonlinear_options.max_line_search_iterations) {
761  ++lineSearchIter;
762 
763  doglegStep(trResults.cauchy_point, trResults.z, tr_size, trResults.d);
764 
765  bool use_with_option1 =
766  (subspace_option >= 1) && (trResults.interior_status == TrustRegionResults::Status::NonDescentDirection ||
767  trResults.interior_status == TrustRegionResults::Status::NegativeCurvature ||
768  ((Norm(trResults.d) > (1.0 - 1.0e-6) * tr_size) && lineSearchIter > 1));
769  bool use_with_option2 = (subspace_option >= 2) && (Norm(trResults.d) > (1.0 - 1.0e-6) * tr_size);
770  bool use_with_option3 = (subspace_option >= 3);
771 
772  if (use_with_option1 || use_with_option2 || use_with_option3) {
773  if (!have_computed_Hvs) {
774  have_computed_Hvs = true;
775  hess_vec_func(trResults.z, trResults.H_z);
776  hess_vec_func(trResults.d_old, trResults.H_d_old);
777  hess_vec_func(trResults.cauchy_point, trResults.H_cauchy_point);
778  }
779 
780  H_left_mosts.clear();
781  for (auto& left : left_mosts) {
782  H_left_mosts.emplace_back(std::make_shared<mfem::Vector>(*left));
783  hess_vec_func(*left, *H_left_mosts.back());
784  }
785 
786  std::vector<const mfem::Vector*> ds{&trResults.z, &trResults.d_old, &trResults.cauchy_point};
787  std::vector<const mfem::Vector*> H_ds{&trResults.H_z, &trResults.H_d_old, &trResults.H_cauchy_point};
788  solveTheSubspaceProblem(trResults.d, hess_vec_func, ds, H_ds, r, tr_size, num_leftmost);
789  }
790 
791  static constexpr double roundOffTol = 0.0; // 1e-14;
792 
793  hess_vec_func(trResults.d, trResults.H_d);
794  double dHd = Dot(trResults.d, trResults.H_d);
795  double modelObjective = Dot(r, trResults.d) + 0.5 * dHd - roundOffTol;
796 
797  add(X, trResults.d, x_pred);
798 
799  double realObjective = std::numeric_limits<double>::max();
800  double normPred = std::numeric_limits<double>::max();
801  try {
802  normPred = computeResidual(x_pred, r_pred);
803  double obj1 = 0.5 * (Dot(r, trResults.d) + Dot(r_pred, trResults.d)) - roundOffTol;
804  realObjective = obj1;
805  } catch (const std::exception&) {
806  realObjective = std::numeric_limits<double>::max();
808  }
809 
810  if (normPred <= norm_goal) {
811  trResults.d_old = trResults.d;
812  X = x_pred;
813  r = r_pred;
814  norm = normPred;
815  if (print_level >= 2) {
816  printTrustRegionInfo(realObjective, modelObjective, trResults.cg_iterations_count, tr_size, true);
817  trResults.cg_iterations_count =
818  0; // zero this output so it doesn't look like the linesearch is doing cg iterations
819  }
820  break;
821  }
822 
823  double modelImprove = -modelObjective;
824  double realImprove = -realObjective;
825 
826  double rho = realImprove / modelImprove;
827  if (modelObjective > 0) {
828  if (print_level >= 2) {
829  mfem::out << "Found a positive model objective increase. Debug if you see this.\n";
830  }
831  rho = realImprove / -modelImprove;
832  }
833 
834  // std::cout << "rho , stuff = " << rho << " " << settings.eta3 << std::endl;
835  // std::cout << "stat = "<< trResults.interior_status << std::endl;
836 
837  if (!(rho >= settings.eta2) ||
838  rho > settings.eta4) { // not enough progress, decrease trust region. write it this way to handle NaNs.
839  tr_size *= settings.t1;
840  } else if ((rho > settings.eta3 && rho <= settings.eta4 &&
841  trResults.interior_status == TrustRegionResults::Status::OnBoundary) ||
842  (rho > 0.95 && rho < 1.05 &&
843  trResults.interior_status ==
844  TrustRegionResults::Status::NegativeCurvature)) { // good progress, on boundary, increase trust
845  // region
846  tr_size *= settings.t2;
847  }
848 
849  // eventually extend to handle this case to handle occasional roundoff issues
850  // modelRes = g + Jd
851  // modelResNorm = np.linalg.norm(modelRes)
852  // realResNorm = np.linalg.norm(gy)
853  bool willAccept = rho >= settings.eta1 && rho <= settings.eta4; // or (rho >= -0 and realResNorm <= gNorm)
854 
855  if (print_level >= 2) {
856  printTrustRegionInfo(realObjective, modelObjective, trResults.cg_iterations_count, tr_size, willAccept);
857  trResults.cg_iterations_count =
858  0; // zero this output so it doesn't look like the linesearch is doing cg iterations
859  }
860 
861  if (willAccept) {
862  trResults.d_old = trResults.d;
863  X = x_pred;
864  r = r_pred;
865  norm = normPred;
866  break;
867  }
868  }
869  }
870 
871  final_iter = it;
872  final_norm = norm;
873 
874  if (print_level == 1) {
875  mfem::out << "TrustRegion iteration " << std::setw(3) << final_iter << " : ||r|| = " << std::setw(13) << norm
876  << '\n';
877  }
878  if (!converged && print_level >= 1) { // (print_options.summary || print_options.warnings)) {
879  mfem::out << "TrustRegion: No convergence!\n";
880  }
881 
882  if (false && print_level >= 2) {
883  mfem::out << "num hess vecs = " << num_hess_vecs << "\n";
884  mfem::out << "num preconds = " << num_preconds << "\n";
885  mfem::out << "num residuals = " << num_residuals << "\n";
886  mfem::out << "num subspace solves = " << num_subspace_solves << "\n";
887  mfem::out << "num jacobian_assembles = " << num_jacobian_assembles << "\n";
888  }
889  }
890 };
891 
893 {
894  auto [lin_solver, preconditioner] = buildLinearSolverAndPreconditioner(lin_opts, comm);
895 
896  lin_solver_ = std::move(lin_solver);
897  preconditioner_ = std::move(preconditioner);
898  nonlin_solver_ = buildNonlinearSolver(nonlinear_opts, lin_opts, *preconditioner_, comm);
899 }
900 
901 EquationSolver::EquationSolver(std::unique_ptr<mfem::NewtonSolver> nonlinear_solver,
902  std::unique_ptr<mfem::Solver> linear_solver,
903  std::unique_ptr<mfem::Solver> preconditioner)
904 {
905  SLIC_ERROR_ROOT_IF(!nonlinear_solver, "Nonlinear solvers must be given to construct an EquationSolver");
906  SLIC_ERROR_ROOT_IF(!linear_solver, "Linear solvers must be given to construct an EquationSolver");
907 
908  nonlin_solver_ = std::move(nonlinear_solver);
909  lin_solver_ = std::move(linear_solver);
910  preconditioner_ = std::move(preconditioner);
911 }
912 
913 void EquationSolver::setOperator(const mfem::Operator& op)
914 {
915  nonlin_solver_->SetOperator(op);
916 
917  // Now that the nonlinear solver knows about the operator, we can set its linear solver
918  if (!nonlin_solver_set_solver_called_) {
919  nonlin_solver_->SetSolver(linearSolver());
920  nonlin_solver_set_solver_called_ = true;
921  }
922 }
923 
924 void EquationSolver::solve(mfem::Vector& x) const
925 {
926  mfem::Vector zero(x);
927  zero = 0.0;
928  // KINSOL does not handle non-zero RHS, so we enforce that the RHS
929  // of the nonlinear system is zero
930  nonlin_solver_->Mult(zero, x);
931 }
932 
933 void SuperLUSolver::Mult(const mfem::Vector& input, mfem::Vector& output) const
934 {
935  SLIC_ERROR_ROOT_IF(!superlu_mat_, "Operator must be set prior to solving with SuperLU");
936 
937  // Use the underlying MFEM-based solver and SuperLU matrix type to solve the system
938  superlu_solver_.Mult(input, output);
939 }
940 
941 std::unique_ptr<mfem::HypreParMatrix> buildMonolithicMatrix(const mfem::BlockOperator& block_operator)
942 {
943  int row_blocks = block_operator.NumRowBlocks();
944  int col_blocks = block_operator.NumColBlocks();
945 
946  SLIC_ERROR_ROOT_IF(row_blocks != col_blocks, "Attempted to use a direct solver on a non-square block system.");
947 
948  mfem::Array2D<const mfem::HypreParMatrix*> hypre_blocks(row_blocks, col_blocks);
949 
950  for (int i = 0; i < row_blocks; ++i) {
951  for (int j = 0; j < col_blocks; ++j) {
952  // checks for presence of empty (null) blocks, which happen fairly common in multirank contact
953  if (!block_operator.IsZeroBlock(i, j)) {
954  auto* hypre_block = dynamic_cast<const mfem::HypreParMatrix*>(&block_operator.GetBlock(i, j));
955  SLIC_ERROR_ROOT_IF(!hypre_block,
956  "Trying to use SuperLU on a block operator that does not contain HypreParMatrix blocks.");
957 
958  hypre_blocks(i, j) = hypre_block;
959  } else {
960  hypre_blocks(i, j) = nullptr;
961  }
962  }
963  }
964 
965  // Note that MFEM passes ownership of this matrix to the caller
966  return std::unique_ptr<mfem::HypreParMatrix>(mfem::HypreParMatrixFromBlocks(hypre_blocks));
967 }
968 
969 void SuperLUSolver::SetOperator(const mfem::Operator& op)
970 {
971  // Check if this is a block operator
972  auto* block_operator = dynamic_cast<const mfem::BlockOperator*>(&op);
973 
974  // If it is, make a monolithic system from the underlying blocks
975  if (block_operator) {
976  auto monolithic_mat = buildMonolithicMatrix(*block_operator);
977 
978  superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*monolithic_mat);
979  } else {
980  // If this is not a block system, check that the input operator is a HypreParMatrix as expected
981  auto* matrix = dynamic_cast<const mfem::HypreParMatrix*>(&op);
982 
983  SLIC_ERROR_ROOT_IF(!matrix, "Matrix must be an assembled HypreParMatrix for use with SuperLU");
984 
985  superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*matrix);
986  }
987 
988  superlu_solver_.SetOperator(*superlu_mat_);
989 }
990 
991 #ifdef MFEM_USE_STRUMPACK
992 
993 void StrumpackSolver::Mult(const mfem::Vector& input, mfem::Vector& output) const
994 {
995  SLIC_ERROR_ROOT_IF(!strumpack_mat_, "Operator must be set prior to solving with Strumpack");
996 
997  // Use the underlying MFEM-based solver and Strumpack matrix type to solve the system
998  strumpack_solver_.Mult(input, output);
999 }
1000 
1001 void StrumpackSolver::SetOperator(const mfem::Operator& op)
1002 {
1003  // Check if this is a block operator
1004  auto* block_operator = dynamic_cast<const mfem::BlockOperator*>(&op);
1005 
1006  // If it is, make a monolithic system from the underlying blocks
1007  if (block_operator) {
1008  auto monolithic_mat = buildMonolithicMatrix(*block_operator);
1009 
1010  strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*monolithic_mat);
1011  } else {
1012  // If this is not a block system, check that the input operator is a HypreParMatrix as expected
1013  auto* matrix = dynamic_cast<const mfem::HypreParMatrix*>(&op);
1014 
1015  SLIC_ERROR_ROOT_IF(!matrix, "Matrix must be an assembled HypreParMatrix for use with Strumpack");
1016 
1017  strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*matrix);
1018  }
1019  height = op.Height();
1020  width = op.Width();
1021  strumpack_solver_.SetOperator(*strumpack_mat_);
1022 }
1023 
1024 #endif
1025 
1026 std::unique_ptr<mfem::NewtonSolver> buildNonlinearSolver(NonlinearSolverOptions nonlinear_opts,
1027  const LinearSolverOptions& linear_opts, mfem::Solver& prec,
1028  MPI_Comm comm)
1029 {
1030  std::unique_ptr<mfem::NewtonSolver> nonlinear_solver;
1031 
1032  if (nonlinear_opts.nonlin_solver == NonlinearSolver::Newton) {
1033  nonlinear_opts.max_line_search_iterations = 0;
1034  SLIC_ERROR_ROOT_IF(nonlinear_opts.min_iterations != 0, "Newton's method does not support nonzero min_iterations");
1035  nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts);
1036  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::LBFGS) {
1037  nonlinear_opts.max_line_search_iterations = 0;
1038  SLIC_ERROR_ROOT_IF(nonlinear_opts.min_iterations != 0, "LBFGS does not support nonzero min_iterations");
1039  nonlinear_solver = std::make_unique<mfem::LBFGSSolver>(comm);
1040  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::NewtonLineSearch) {
1041  nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts);
1042  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::TrustRegion) {
1043  nonlinear_solver = std::make_unique<TrustRegion>(comm, nonlinear_opts, linear_opts, prec);
1044 #ifdef SMITH_USE_PETSC
1045  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscNewton) {
1046  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1047  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscNewtonBacktracking) {
1048  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1049  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscNewtonCriticalPoint) {
1050  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1051  } else if (nonlinear_opts.nonlin_solver == NonlinearSolver::PetscTrustRegion) {
1052  nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1053 #endif
1054  }
1055  // KINSOL
1056  else {
1057 #ifdef SMITH_USE_SUNDIALS
1058  nonlinear_opts.max_line_search_iterations = 0;
1059  SLIC_ERROR_ROOT_IF(nonlinear_opts.min_iterations != 0, "kinsol solvers do not support min_iterations");
1060 
1061  int kinsol_strat = KIN_NONE;
1062 
1063  switch (nonlinear_opts.nonlin_solver) {
1065  kinsol_strat = KIN_NONE;
1066  break;
1068  kinsol_strat = KIN_LINESEARCH;
1069  break;
1071  kinsol_strat = KIN_PICARD;
1072  break;
1073  default:
1074  kinsol_strat = KIN_NONE;
1075  SLIC_ERROR_ROOT("Unknown KINSOL nonlinear solver type given.");
1076  }
1077  auto kinsol_solver = std::make_unique<mfem::KINSolver>(comm, kinsol_strat, true);
1078  nonlinear_solver = std::move(kinsol_solver);
1079 #else
1080  SLIC_ERROR_ROOT("KINSOL was not enabled when MFEM was built");
1081 #endif
1082  }
1083 
1084  nonlinear_solver->SetRelTol(nonlinear_opts.relative_tol);
1085  nonlinear_solver->SetAbsTol(nonlinear_opts.absolute_tol);
1086  nonlinear_solver->SetMaxIter(nonlinear_opts.max_iterations);
1087  nonlinear_solver->SetPrintLevel(nonlinear_opts.print_level);
1088 
1089  // Iterative mode indicates we do not zero out the initial guess during the
1090  // nonlinear solver call. This is required as we apply the essential boundary
1091  // conditions before the nonlinear solver is applied.
1092  nonlinear_solver->iterative_mode = true;
1093 
1094  return nonlinear_solver;
1095 }
1096 
1097 std::pair<std::unique_ptr<mfem::Solver>, std::unique_ptr<mfem::Solver>> buildLinearSolverAndPreconditioner(
1098  LinearSolverOptions linear_opts, MPI_Comm comm)
1099 {
1100  auto preconditioner = buildPreconditioner(linear_opts, comm);
1101 
1102  if (linear_opts.linear_solver == LinearSolver::SuperLU) {
1103  auto lin_solver = std::make_unique<SuperLUSolver>(linear_opts.print_level, comm);
1104  return {std::move(lin_solver), std::move(preconditioner)};
1105  }
1106 
1107 #ifdef MFEM_USE_STRUMPACK
1108 
1109  if (linear_opts.linear_solver == LinearSolver::Strumpack) {
1110  auto lin_solver = std::make_unique<StrumpackSolver>(linear_opts.print_level, comm);
1111  return {std::move(lin_solver), std::move(preconditioner)};
1112  }
1113 
1114 #endif
1115 
1116  std::unique_ptr<mfem::IterativeSolver> iter_lin_solver;
1117 
1118  switch (linear_opts.linear_solver) {
1119  case LinearSolver::CG:
1120  iter_lin_solver = std::make_unique<mfem::CGSolver>(comm);
1121  break;
1122  case LinearSolver::GMRES:
1123  iter_lin_solver = std::make_unique<mfem::GMRESSolver>(comm);
1124  break;
1125 #ifdef SMITH_USE_PETSC
1126  case LinearSolver::PetscCG:
1127  iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPCG, std::string());
1128  break;
1130  iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPGMRES, std::string());
1131  break;
1132 #else
1133  case LinearSolver::PetscCG:
1135  SLIC_ERROR_ROOT("PETSc linear solver requested for non-PETSc build.");
1136  exit(1);
1137  break;
1138 #endif
1139  default:
1140  SLIC_ERROR_ROOT("Linear solver type not recognized.");
1141  exit(1);
1142  }
1143 
1144  iter_lin_solver->SetRelTol(linear_opts.relative_tol);
1145  iter_lin_solver->SetAbsTol(linear_opts.absolute_tol);
1146  iter_lin_solver->SetMaxIter(linear_opts.max_iterations);
1147  iter_lin_solver->SetPrintLevel(linear_opts.print_level);
1148 
1149  if (preconditioner) {
1150  iter_lin_solver->SetPreconditioner(*preconditioner);
1151  }
1152 
1153  return {std::move(iter_lin_solver), std::move(preconditioner)};
1154 }
1155 
1156 #ifdef MFEM_USE_AMGX
1157 std::unique_ptr<mfem::AmgXSolver> buildAMGX(const AMGXOptions& options, const MPI_Comm comm)
1158 {
1159  auto amgx = std::make_unique<mfem::AmgXSolver>();
1160  conduit::Node options_node;
1161  options_node["config_version"] = 2;
1162  auto& solver_options = options_node["solver"];
1163  solver_options["solver"] = "AMG";
1164  solver_options["presweeps"] = 1;
1165  solver_options["postsweeps"] = 2;
1166  solver_options["interpolator"] = "D2";
1167  solver_options["max_iters"] = 2;
1168  solver_options["convergence"] = "ABSOLUTE";
1169  solver_options["cycle"] = "V";
1170 
1171  if (options.verbose) {
1172  options_node["solver/obtain_timings"] = 1;
1173  options_node["solver/monitor_residual"] = 1;
1174  options_node["solver/print_solve_stats"] = 1;
1175  }
1176 
1177  // TODO: Use magic_enum here when we can switch to GCC 9+
1178  // This is an immediately-invoked lambda so that the map
1179  // can be const without needed to initialize all the values
1180  // in the constructor
1181  static const auto solver_names = []() {
1182  std::unordered_map<AMGXSolver, std::string> names;
1183  names[AMGXSolver::AMG] = "AMG";
1184  names[AMGXSolver::PCGF] = "PCGF";
1185  names[AMGXSolver::CG] = "CG";
1186  names[AMGXSolver::PCG] = "PCG";
1187  names[AMGXSolver::PBICGSTAB] = "PBICGSTAB";
1188  names[AMGXSolver::BICGSTAB] = "BICGSTAB";
1189  names[AMGXSolver::FGMRES] = "FGMRES";
1190  names[AMGXSolver::JACOBI_L1] = "JACOBI_L1";
1191  names[AMGXSolver::GS] = "GS";
1192  names[AMGXSolver::POLYNOMIAL] = "POLYNOMIAL";
1193  names[AMGXSolver::KPZ_POLYNOMIAL] = "KPZ_POLYNOMIAL";
1194  names[AMGXSolver::BLOCK_JACOBI] = "BLOCK_JACOBI";
1195  names[AMGXSolver::MULTICOLOR_GS] = "MULTICOLOR_GS";
1196  names[AMGXSolver::MULTICOLOR_DILU] = "MULTICOLOR_DILU";
1197  return names;
1198  }();
1199 
1200  options_node["solver/solver"] = solver_names.at(options.solver);
1201  options_node["solver/smoother"] = solver_names.at(options.smoother);
1202 
1203  // Treat the string as the config (not a filename)
1204  amgx->ReadParameters(options_node.to_json(), mfem::AmgXSolver::INTERNAL);
1205  amgx->InitExclusiveGPU(comm);
1206 
1207  return amgx;
1208 }
1209 #endif
1210 
1211 std::unique_ptr<mfem::Solver> buildPreconditioner(LinearSolverOptions linear_opts, [[maybe_unused]] MPI_Comm comm)
1212 {
1213  std::unique_ptr<mfem::Solver> preconditioner_solver;
1214  auto preconditioner = linear_opts.preconditioner;
1215  auto preconditioner_print_level = linear_opts.preconditioner_print_level;
1216 
1217  // Handle the preconditioner - currently just BoomerAMG and HypreSmoother are supported
1218  if (preconditioner == Preconditioner::HypreAMG) {
1219  auto amg_preconditioner = std::make_unique<mfem::HypreBoomerAMG>();
1220  amg_preconditioner->SetPrintLevel(preconditioner_print_level);
1221  preconditioner_solver = std::move(amg_preconditioner);
1222  } else if (preconditioner == Preconditioner::HypreJacobi) {
1223  auto jac_preconditioner = std::make_unique<mfem::HypreSmoother>();
1224  jac_preconditioner->SetType(mfem::HypreSmoother::Type::Jacobi);
1225  preconditioner_solver = std::move(jac_preconditioner);
1226  } else if (preconditioner == Preconditioner::HypreL1Jacobi) {
1227  auto jacl1_preconditioner = std::make_unique<mfem::HypreSmoother>();
1228  jacl1_preconditioner->SetType(mfem::HypreSmoother::Type::l1Jacobi);
1229  preconditioner_solver = std::move(jacl1_preconditioner);
1230  } else if (preconditioner == Preconditioner::HypreGaussSeidel) {
1231  auto gs_preconditioner = std::make_unique<mfem::HypreSmoother>();
1232  gs_preconditioner->SetType(mfem::HypreSmoother::Type::GS);
1233  preconditioner_solver = std::move(gs_preconditioner);
1234  } else if (preconditioner == Preconditioner::HypreILU) {
1235  auto ilu_preconditioner = std::make_unique<mfem::HypreILU>();
1236  ilu_preconditioner->SetLevelOfFill(1);
1237  ilu_preconditioner->SetPrintLevel(preconditioner_print_level);
1238  preconditioner_solver = std::move(ilu_preconditioner);
1239  } else if (preconditioner == Preconditioner::AMGX) {
1240 #ifdef MFEM_USE_AMGX
1241  preconditioner_solver = buildAMGX(linear_opts.amgx_options, comm);
1242 #else
1243  SLIC_ERROR_ROOT("AMGX requested in non-GPU build");
1244 #endif
1245  } else if (preconditioner == Preconditioner::Petsc) {
1246 #ifdef SMITH_USE_PETSC
1247  preconditioner_solver = mfem_ext::buildPetscPreconditioner(linear_opts.petsc_preconditioner, comm);
1248 #else
1249  SLIC_ERROR_ROOT("PETSc preconditioner requested in non-PETSc build");
1250 #endif
1251  } else if (preconditioner == Preconditioner::AMGFContact) {
1252  auto amgfcontact_preconditioner = std::make_unique<mfem::AMGFSolver>();
1253  auto amgfcontact_opts = linear_opts.amgfcontact_options;
1254  amgfcontact_preconditioner->GetAMG().SetPrintLevel(preconditioner_print_level);
1255  amgfcontact_preconditioner->GetAMG().SetSystemsOptions(amgfcontact_opts.dim_systems_options);
1256  amgfcontact_preconditioner->GetAMG().SetRelaxType(amgfcontact_opts.relax_type);
1257  preconditioner_solver = std::move(amgfcontact_preconditioner);
1258  } else {
1259  SLIC_ERROR_ROOT_IF(preconditioner != Preconditioner::None, "Unknown preconditioner type requested");
1260  }
1261 
1262  return preconditioner_solver;
1263 }
1264 
1265 void EquationSolver::defineInputFileSchema(axom::inlet::Container& container)
1266 {
1267  auto& linear_container = container.addStruct("linear", "Linear Equation Solver Parameters");
1268  linear_container.required().registerVerifier([](const axom::inlet::Container& container_to_verify) {
1269  // Make sure that the provided options match the desired linear solver type
1270  const bool is_iterative = (container_to_verify["type"].get<std::string>() == "iterative") &&
1271  container_to_verify.contains("iterative_options");
1272  const bool is_direct =
1273  (container_to_verify["type"].get<std::string>() == "direct") && container_to_verify.contains("direct_options");
1274  return is_iterative || is_direct;
1275  });
1276 
1277  // Enforce the solver type - must be iterative or direct
1278  linear_container.addString("type", "The type of solver parameters to use (iterative|direct)")
1279  .required()
1280  .validValues({"iterative", "direct"});
1281 
1282  auto& iterative_container = linear_container.addStruct("iterative_options", "Iterative solver parameters");
1283  iterative_container.addDouble("rel_tol", "Relative tolerance for the linear solve.").defaultValue(1.0e-6);
1284  iterative_container.addDouble("abs_tol", "Absolute tolerance for the linear solve.").defaultValue(1.0e-8);
1285  iterative_container.addInt("max_iter", "Maximum iterations for the linear solve.").defaultValue(5000);
1286  iterative_container.addInt("print_level", "Linear print level.").defaultValue(0);
1287  iterative_container.addString("solver_type", "Solver type (gmres|minres|cg).").defaultValue("gmres");
1288  iterative_container.addString("prec_type", "Preconditioner type (JacobiSmoother|L1JacobiSmoother|AMG|ILU|Petsc).")
1289  .defaultValue("JacobiSmoother");
1290  iterative_container.addString("petsc_prec_type", "Type of PETSc preconditioner to use.").defaultValue("jacobi");
1291 
1292  auto& direct_container = linear_container.addStruct("direct_options", "Direct solver parameters");
1293  direct_container.addInt("print_level", "Linear print level.").defaultValue(0);
1294 
1295  // Only needed for nonlinear problems
1296  auto& nonlinear_container = container.addStruct("nonlinear", "Newton Equation Solver Parameters").required(false);
1297  nonlinear_container.addDouble("rel_tol", "Relative tolerance for the Newton solve.").defaultValue(1.0e-2);
1298  nonlinear_container.addDouble("abs_tol", "Absolute tolerance for the Newton solve.").defaultValue(1.0e-4);
1299  nonlinear_container.addInt("max_iter", "Maximum iterations for the Newton solve.").defaultValue(500);
1300  nonlinear_container.addInt("print_level", "Nonlinear print level.").defaultValue(0);
1301  nonlinear_container.addString("solver_type", "Solver type (Newton|KINFullStep|KINLineSearch)").defaultValue("Newton");
1302 }
1303 
1304 } // namespace smith
1305 
1306 using smith::EquationSolver;
1309 
1311 {
1312  LinearSolverOptions options;
1313  std::string type = base["type"];
1314 
1315  if (type == "direct") {
1317  options.print_level = base["direct_options/print_level"];
1318  return options;
1319  }
1320 
1321  auto config = base["iterative_options"];
1322  options.relative_tol = config["rel_tol"];
1323  options.absolute_tol = config["abs_tol"];
1324  options.max_iterations = config["max_iter"];
1325  options.print_level = config["print_level"];
1326  std::string solver_type = config["solver_type"];
1327  if (solver_type == "gmres") {
1329  } else if (solver_type == "cg") {
1331  } else {
1332  std::string msg = std::format("Unknown Linear solver type given: '{0}'", solver_type);
1333  SLIC_ERROR_ROOT(msg);
1334  }
1335  const std::string prec_type = config["prec_type"];
1336  if (prec_type == "JacobiSmoother") {
1338  } else if (prec_type == "L1JacobiSmoother") {
1340  } else if (prec_type == "HypreAMG") {
1342  } else if (prec_type == "ILU") {
1344 #ifdef MFEM_USE_AMGX
1345  } else if (prec_type == "AMGX") {
1347 #endif
1348  } else if (prec_type == "GaussSeidel") {
1350 #ifdef SMITH_USE_PETSC
1351  } else if (prec_type == "Petsc") {
1352  const std::string petsc_prec = config["petsc_prec_type"];
1354  options.petsc_preconditioner = smith::mfem_ext::stringToPetscPCType(petsc_prec);
1355 #endif
1356  } else if (prec_type == "AMGFContact") {
1358  } else {
1359  std::string msg = std::format("Unknown preconditioner type given: '{0}'", prec_type);
1360  SLIC_ERROR_ROOT(msg);
1361  }
1362 
1363  return options;
1364 }
1365 
1367 {
1368  NonlinearSolverOptions options;
1369  options.relative_tol = base["rel_tol"];
1370  options.absolute_tol = base["abs_tol"];
1371  options.max_iterations = base["max_iter"];
1372  options.print_level = base["print_level"];
1373  const std::string solver_type = base["solver_type"];
1374  if (solver_type == "Newton") {
1376  } else if (solver_type == "KINFullStep") {
1378  } else if (solver_type == "KINLineSearch") {
1380  } else if (solver_type == "KINPicard") {
1382  } else {
1383  SLIC_ERROR_ROOT(std::format("Unknown nonlinear solver type given: '{0}'", solver_type));
1384  }
1385  return options;
1386 }
1387 
1389 {
1390  auto lin = base["linear"].get<LinearSolverOptions>();
1391  auto nonlin = base["nonlinear"].get<NonlinearSolverOptions>();
1392 
1393  auto [linear_solver, preconditioner] = smith::buildLinearSolverAndPreconditioner(lin, MPI_COMM_WORLD);
1394 
1395  smith::EquationSolver eq_solver(smith::buildNonlinearSolver(nonlin, lin, *preconditioner, MPI_COMM_WORLD),
1396  std::move(linear_solver), std::move(preconditioner));
1397 
1398  return eq_solver;
1399 }
This class manages the objects typically required to solve a nonlinear set of equations arising from ...
mfem::Solver & preconditioner()
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)
Newton solver with a 2-way line-search. Reverts to regular Newton if max_line_search_iterations is se...
size_t print_level
reconstructed smith print level
mfem::Vector x0
initial solution vector to do line-search off of
void setPreconditioner() const
set the preconditioner for the linear solver
void Mult(const mfem::Vector &, mfem::Vector &x) const
This is an overloaded member function, provided for convenience. It differs from the above function o...
void assembleJacobian(const mfem::Vector &x) const
assemble the jacobian
void solveLinearSystem(const mfem::Vector &r_, mfem::Vector &c_) const
solve the linear system
NonlinearSolverOptions nonlinear_options
nonlinear solver options
double evaluateNorm(const mfem::Vector &x, mfem::Vector &rOut) const
Evaluate the residual, put in rOut and return its norm.
NewtonSolver(const NonlinearSolverOptions &nonlinear_opts)
constructor
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.
Equation solver class based on a standard preconditioned trust-region algorithm.
std::vector< std::shared_ptr< mfem::Vector > > H_left_mosts
the action of the stiffness/hessian (H) on the left most eigenvectors
size_t num_residuals
internal counter for residuals
void assembleJacobian(const mfem::Vector &x) const
assemble the jacobian
void Mult(const mfem::Vector &, mfem::Vector &X) const
This is an overloaded member function, provided for convenience. It differs from the above function o...
void solveTrustRegionModelProblem(const mfem::Vector &r0, mfem::Vector &rCurrent, HessVecFunc hess_vec_func, PrecondFunc precond, const TrustRegionSettings &settings, double &trSize, TrustRegionResults &results) const
Minimize quadratic sub-problem given residual vector, the action of the stiffness and a preconditione...
void hessVec(const mfem::Vector &x_, mfem::Vector &v_) const
apply the action of the assembled Jacobian matrix to a vector
NonlinearSolverOptions nonlinear_options
nonlinear solution options
void solveTheSubspaceProblem([[maybe_unused]] mfem::Vector &z, [[maybe_unused]] const HessVecFunc &hess_vec_func, [[maybe_unused]] const std::vector< const mfem::Vector * > ds, [[maybe_unused]] const std::vector< const mfem::Vector * > Hds, [[maybe_unused]] const mfem::Vector &g, [[maybe_unused]] double delta, [[maybe_unused]] int num_leftmost) const
solve the exact trust-region subspace problem with directions ds, and the leftmosts
std::vector< std::shared_ptr< mfem::Vector > > left_mosts
left most eigenvectors
mfem::Vector x_pred
predicted solution
size_t num_subspace_solves
internal counter for subspace solves
void projectToBoundaryBetweenWithCoefs(mfem::Vector &z, const mfem::Vector &y, double trSize, double zz, double zy, double yy) const
finds tau s.t. (z + tau*(y-z))^2 = trSize^2
void precond(const mfem::Vector &x_, mfem::Vector &v_) const
apply trust region specific preconditioner
mfem::Vector scratch
scratch
size_t num_jacobian_assembles
internal counter for matrix assembles
void projectToBoundaryWithCoefs(mfem::Vector &z, const mfem::Vector &d, double delta, double zz, double zd, double dd) const
finds tau s.t. (z + tau*d)^2 = trSize^2
size_t num_preconds
internal counter for preconditions
mfem::Vector r_pred
predicted residual
mfem::real_t computeResidual(const mfem::Vector &x_, mfem::Vector &r_) const
evaluate the nonlinear residual
size_t num_hess_vecs
internal counter for hess-vecs
size_t print_level
reconstructed smith print level
LinearSolverOptions linear_options
linear solution options
void doglegStep(const mfem::Vector &cp, const mfem::Vector &newtonP, double trSize, mfem::Vector &s) const
take a dogleg step in direction s, solution norm must be within trSize
double computeEnergy(const mfem::Vector &r_local, const HessVecFunc &H, const mfem::Vector &z) const
compute the energy of the linearized system for a given solution vector z
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.
std::unique_ptr< mfem::HypreParMatrix > buildMonolithicMatrix(const mfem::BlockOperator &block_operator)
Function for building a monolithic parallel Hypre matrix from a block system of smaller Hypre matrice...
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.
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
constexpr SMITH_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
Definition: tuple.hpp:376
void printTrustRegionInfo(double realObjective, double modelObjective, size_t cgIters, double trSize, bool willAccept)
trust region printing utility function
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
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.
PetscPCType petsc_preconditioner
PETSc preconditioner type.
double absolute_tol
Absolute tolerance.
AMGXOptions amgx_options
AMGX Options, used for Preconditioner::AMGX.
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.
bool force_monolithic
Should the gradient be converted to a monolithic matrix.
Internal structure for storing trust region stateful data.
size_t cg_iterations_count
iteration counter
mfem::Vector cauchy_point
cauchy point
mfem::Vector H_z
action of hessian on current step z
mfem::Vector z
step direction
mfem::Vector d
incrementalCG direction
mfem::Vector H_cauchy_point
action of hessian on direction of cauchy point
mfem::Vector H_d_old
action of hessian on previous step z_old
void reset()
resets trust region results for a new outer iteration
mfem::Vector H_d
action of hessian on direction d
TrustRegionResults(int size)
Constructor takes the size of the solution vector.
mfem::Vector Pr
preconditioned residual
Status
enumerates the possible final status of the trust region steps
Status interior_status
specifies if step is interior, exterior, negative curvature, etc.
mfem::Vector d_old
old step direction
Internal structure for storing trust region settings.
double eta2
non-ideal energy drop ratio. trust region decreases if energy drop is worse than this.
size_t max_cg_iterations
max cg iters should be around # of system dofs
double eta1
worse case energy drop ratio. trust region accepted if energy drop is better than this.
double eta3
ideal energy drop ratio. trust region increases if energy drop is better than this.
size_t min_cg_iterations
min cg iters
size_t max_cumulative_iteration
max cumulative iterations
double min_tr_size
minimum trust region size
double t1
trust region decrease factor
double eta4
parameter limiting how fast the energy can drop relative to the prediction (in case the energy surrog...
double t2
trust region increase factor
A sentinel struct for eliding no-op tensor operations.
Definition: tensor.hpp:122