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