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