19 #include "serac/serac_config.hpp"
21 #include "serac/numerics/trust_region_solver.hpp"
30 mutable mfem::Vector
x0;
47 double evaluateNorm(
const mfem::Vector& x, mfem::Vector& rOut)
const
53 normEval = Norm(rOut);
54 }
catch (
const std::exception&) {
64 grad = &oper->GetGradient(x);
66 auto* grad_blocked =
dynamic_cast<mfem::BlockOperator*
>(grad);
75 prec->SetOperator(*grad);
86 void Mult(
const mfem::Vector&, mfem::Vector& x)
const
88 MFEM_ASSERT(oper != NULL,
"the Operator is not set (use SetOperator).");
89 MFEM_ASSERT(prec != NULL,
"the Solver is not set (use SetSolver).");
91 using real_t = mfem::real_t;
93 real_t
norm, norm_goal = 0;
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";
100 norm_goal =
std::max(rel_tol * initial_norm, abs_tol);
101 prec->iterative_mode =
false;
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;
109 mfem::out <<
", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ?
norm / initial_norm :
norm);
115 mfem::out <<
"Initial residual for Newton iteration is undefined/nan.\n";
116 mfem::out <<
"Newton: No convergence!\n";
123 }
else if (it >= max_iter) {
128 real_t norm_nm1 =
norm;
135 x0.SetSize(x.Size());
139 real_t stepScale = 1.0;
140 add(
x0, -stepScale, c, x);
144 static constexpr real_t reduction = 0.5;
146 const double sufficientDecreaseParam = 0.0;
147 const double cMagnitudeInR = sufficientDecreaseParam != 0.0 ?
std::abs(Dot(c, r)) / norm_nm1 : 0.0;
149 auto is_improved = [=](real_t currentNorm, real_t c_scale) {
150 return currentNorm < norm_nm1 - sufficientDecreaseParam * c_scale * cMagnitudeInR;
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);
163 if (max_ls_iters > 0 && ls_iter == max_ls_iters && !is_improved(
norm, stepScale)) {
165 add(
x0, stepScale, c, x);
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);
177 if (ls_iter == max_ls_iters && !is_improved(
norm, stepScale)) {
179 stepScale *= reduction;
180 add(
x0, -stepScale, c, x);
186 if (print_options.iterations) {
187 mfem::out <<
"Number of line search steps taken = " << ls_iter_sum << std::endl;
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 "
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';
203 if (!converged && (print_options.summary || print_options.warnings)) {
204 mfem::out <<
"Newton: No convergence!\n";
292 void printTrustRegionInfo(
double realObjective,
double modelObjective,
size_t cgIters,
double trSize,
bool willAccept)
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;
317 mutable std::vector<std::shared_ptr<mfem::Vector>>
left_mosts;
355 double deltadelta_m_zz = delta * delta - zz;
356 if (deltadelta_m_zz == 0)
return;
357 double tau = (
std::sqrt(deltadelta_m_zz * dd + zd * zd) - zd) / dd;
362 template <
typename HessVecFunc>
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
369 #ifdef SERAC_USE_SLEPC
373 std::vector<const mfem::Vector*> directions;
375 directions.emplace_back(d);
378 directions.emplace_back(left.get());
381 std::vector<const mfem::Vector*> H_directions;
382 for (
auto& Hd : Hds) {
383 H_directions.emplace_back(Hd);
386 H_directions.emplace_back(H_left.get());
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;
402 std::vector<std::shared_ptr<mfem::Vector>> leftvecs;
403 std::vector<double> leftvals;
404 double energy_change;
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;
417 for (
auto& lv : leftvecs) {
422 double subspace_energy =
computeEnergy(g, hess_vec_func, sol);
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;
430 if (subspace_energy < base_energy) {
440 double dd = yy - 2 * zy + zz;
442 double tau = (
std::sqrt((trSize * trSize - zz) * dd + zd * zd) - zd) / dd;
448 void doglegStep(
const mfem::Vector& cp,
const mfem::Vector& newtonP,
double trSize, mfem::Vector& s)
const
452 double cc = Dot(cp, cp);
453 double nn = Dot(newtonP, newtonP);
454 double tt = trSize * trSize;
459 }
else if (cc > nn) {
460 if (print_options.warnings) {
461 mfem::out <<
"cp outside newton, preconditioner likely inaccurate\n";
464 }
else if (nn > tt) {
466 double cn = Dot(cp, newtonP);
474 template <
typename HessVecFunc>
475 double computeEnergy(
const mfem::Vector& r_local,
const HessVecFunc& H,
const mfem::Vector& z)
const
478 double rz = Dot(r_local, z);
479 mfem::Vector tmp(r_local);
482 return rz + 0.5 * Dot(z, tmp);
486 template <
typename HessVecFunc,
typename PrecondFunc>
499 auto& Pr = results.
Pr;
500 auto& Hd = results.
H_d;
502 const double cg_tol_squared = settings.
cg_tol * settings.
cg_tol;
505 mfem::out <<
"Trust region solution state within tolerance on first iteration."
519 double rPr = Dot(rCurrent, Pr);
521 double dd = Dot(d, d);
527 if (Dot(d, rCurrent) > 0) {
529 results.
interior_status = TrustRegionResults::Status::NonDescentDirection;
532 hess_vec_func(d, Hd);
533 const double curvature = Dot(d, Hd);
534 const double alphaCg = curvature != 0.0 ? rPr / curvature : 0.0;
539 add(z, alphaCg, d, zPred);
540 double zzNp1 = Dot(zPred, zPred);
542 const bool go_to_boundary = curvature <= 0 || zzNp1 >= trSize * trSize;
543 if (go_to_boundary) {
545 if (curvature <= 0) {
546 results.
interior_status = TrustRegionResults::Status::NegativeCurvature;
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";
562 add(rCurrent, alphaCg, Hd, rCurrent);
565 double rPrNp1 = Dot(rCurrent, Pr);
567 if (Dot(rCurrent, rCurrent) <= cg_tol_squared && cgIter >= settings.
min_cg_iterations) {
571 double beta = rPrNp1 / rPr;
573 add(-1.0, Pr, beta, d, d);
587 grad = &oper->GetGradient(x);
589 auto* grad_blocked =
dynamic_cast<mfem::BlockOperator*
>(grad);
604 void hessVec(
const mfem::Vector& x_, mfem::Vector& v_)
const
612 void precond(
const mfem::Vector& x_, mfem::Vector& v_)
const
620 void Mult(
const mfem::Vector&, mfem::Vector& X)
const
622 MFEM_ASSERT(oper != NULL,
"the Operator is not set (use SetOperator).");
623 MFEM_ASSERT(prec != NULL,
"the Solver is not set (use SetSolver).");
625 using real_t = mfem::real_t;
633 real_t
norm, norm_goal = 0.0;
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";
639 prec->iterative_mode =
false;
654 settings.
cg_tol = 0.5 * norm_goal;
661 size_t cumulative_cg_iters_from_last_precond_update = 0;
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;
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();
672 mfem::out <<
", norm goal = " << std::setw(13) << norm_goal;
678 mfem::out <<
"Initial residual for trust-region iteration is undefined/nan." << std::endl;
679 mfem::out <<
"Newton: No convergence!\n";
686 }
else if (it >= max_iter) {
696 cumulative_cg_iters_from_last_precond_update = 0;
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_); };
702 double cauchyPointNormSquared = tr_size * tr_size;
705 hess_vec_func(r, trResults.
H_d);
706 const double gKg = Dot(r, trResults.
H_d);
708 const double alphaCp = -Dot(r, r) / gKg;
712 const double alphaTr = -tr_size /
std::sqrt(Dot(r, r));
714 if (print_options.iterations) {
715 mfem::out <<
"Negative curvature un-preconditioned cauchy point direction found."
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";
736 bool have_computed_Hvs =
false;
738 int lineSearchIter = 0;
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);
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);
761 H_left_mosts.emplace_back(std::make_shared<mfem::Vector>(*left));
765 std::vector<const mfem::Vector*> ds{&trResults.
z, &trResults.
d_old, &trResults.
cauchy_point};
770 static constexpr
double roundOffTol = 0.0;
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;
782 double obj1 = 0.5 * (Dot(r, trResults.
d) + Dot(
r_pred, trResults.
d)) - roundOffTol;
783 realObjective = obj1;
784 }
catch (
const std::exception&) {
789 if (normPred <= norm_goal) {
790 trResults.
d_old = trResults.
d;
794 if (print_options.iterations) {
802 double modelImprove = -modelObjective;
803 double realImprove = -realObjective;
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";
810 rho = realImprove / -modelImprove;
816 if (!(rho >= settings.
eta2) ||
817 rho > settings.
eta4) {
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 &&
823 TrustRegionResults::Status::NegativeCurvature)) {
825 tr_size *= settings.
t2;
832 bool willAccept = rho >= settings.
eta1 && rho <= settings.
eta4;
834 if (print_options.iterations) {
841 trResults.
d_old = trResults.
d;
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';
856 if (!converged && (print_options.summary || print_options.warnings)) {
857 mfem::out <<
"Newton: No convergence!\n";
860 if (
false && (print_options.summary || print_options.warnings)) {
874 lin_solver_ = std::move(lin_solver);
880 std::unique_ptr<mfem::Solver> linear_solver,
881 std::unique_ptr<mfem::Solver> preconditioner)
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");
886 nonlin_solver_ = std::move(nonlinear_solver);
887 lin_solver_ = std::move(linear_solver);
893 nonlin_solver_->SetOperator(op);
896 if (!nonlin_solver_set_solver_called_) {
898 nonlin_solver_set_solver_called_ =
true;
904 mfem::Vector
zero(x);
908 nonlin_solver_->Mult(
zero, x);
913 SLIC_ERROR_ROOT_IF(!superlu_mat_,
"Operator must be set prior to solving with SuperLU");
916 superlu_solver_.Mult(input, output);
921 int row_blocks = block_operator.NumRowBlocks();
922 int col_blocks = block_operator.NumColBlocks();
924 SLIC_ERROR_ROOT_IF(row_blocks != col_blocks,
"Attempted to use a direct solver on a non-square block system.");
926 mfem::Array2D<const mfem::HypreParMatrix*> hypre_blocks(row_blocks, col_blocks);
928 for (
int i = 0; i < row_blocks; ++i) {
929 for (
int j = 0; j < col_blocks; ++j) {
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.");
936 hypre_blocks(i, j) = hypre_block;
938 hypre_blocks(i, j) =
nullptr;
944 return std::unique_ptr<mfem::HypreParMatrix>(mfem::HypreParMatrixFromBlocks(hypre_blocks));
950 auto* block_operator =
dynamic_cast<const mfem::BlockOperator*
>(&op);
953 if (block_operator) {
956 superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*monolithic_mat);
959 auto* matrix =
dynamic_cast<const mfem::HypreParMatrix*
>(&op);
961 SLIC_ERROR_ROOT_IF(!matrix,
"Matrix must be an assembled HypreParMatrix for use with SuperLU");
963 superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*matrix);
966 superlu_solver_.SetOperator(*superlu_mat_);
969 #ifdef MFEM_USE_STRUMPACK
971 void StrumpackSolver::Mult(
const mfem::Vector& input, mfem::Vector& output)
const
973 SLIC_ERROR_ROOT_IF(!strumpack_mat_,
"Operator must be set prior to solving with Strumpack");
976 strumpack_solver_.Mult(input, output);
979 void StrumpackSolver::SetOperator(
const mfem::Operator& op)
982 auto* block_operator =
dynamic_cast<const mfem::BlockOperator*
>(&op);
985 if (block_operator) {
988 strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*monolithic_mat);
991 auto* matrix =
dynamic_cast<const mfem::HypreParMatrix*
>(&op);
993 SLIC_ERROR_ROOT_IF(!matrix,
"Matrix must be an assembled HypreParMatrix for use with Strumpack");
995 strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*matrix);
998 strumpack_solver_.SetOperator(*strumpack_mat_);
1007 std::unique_ptr<mfem::NewtonSolver> nonlinear_solver;
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);
1016 "LBFGS does not support nonzero min_iterations or max_line_search_iterations");
1017 nonlinear_solver = std::make_unique<mfem::LBFGSSolver>(comm);
1019 nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts);
1021 nonlinear_solver = std::make_unique<TrustRegion>(comm, nonlinear_opts, linear_opts, prec);
1022 #ifdef SERAC_USE_PETSC
1024 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1026 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1028 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1030 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1035 #ifdef SERAC_USE_SUNDIALS
1038 "kinsol solvers do not support min_iterations or max_line_search_iterations");
1040 int kinsol_strat = KIN_NONE;
1044 kinsol_strat = KIN_NONE;
1047 kinsol_strat = KIN_LINESEARCH;
1050 kinsol_strat = KIN_PICARD;
1053 kinsol_strat = KIN_NONE;
1054 SLIC_ERROR_ROOT(
"Unknown KINSOL nonlinear solver type given.");
1056 auto kinsol_solver = std::make_unique<mfem::KINSolver>(comm, kinsol_strat,
true);
1057 nonlinear_solver = std::move(kinsol_solver);
1059 SLIC_ERROR_ROOT(
"KINSOL was not enabled when MFEM was built");
1063 nonlinear_solver->SetRelTol(nonlinear_opts.
relative_tol);
1064 nonlinear_solver->SetAbsTol(nonlinear_opts.
absolute_tol);
1066 nonlinear_solver->SetPrintLevel(nonlinear_opts.
print_level);
1071 nonlinear_solver->iterative_mode =
true;
1073 return nonlinear_solver;
1082 auto lin_solver = std::make_unique<SuperLUSolver>(linear_opts.
print_level, comm);
1083 return {std::move(lin_solver), std::move(preconditioner)};
1086 #ifdef MFEM_USE_STRUMPACK
1089 auto lin_solver = std::make_unique<StrumpackSolver>(linear_opts.
print_level, comm);
1090 return {std::move(lin_solver), std::move(preconditioner)};
1095 std::unique_ptr<mfem::IterativeSolver> iter_lin_solver;
1099 iter_lin_solver = std::make_unique<mfem::CGSolver>(comm);
1102 iter_lin_solver = std::make_unique<mfem::GMRESSolver>(comm);
1104 #ifdef SERAC_USE_PETSC
1106 iter_lin_solver = std::make_unique<serac::mfem_ext::PetscKSPSolver>(comm, KSPCG, std::string());
1109 iter_lin_solver = std::make_unique<serac::mfem_ext::PetscKSPSolver>(comm, KSPGMRES, std::string());
1114 SLIC_ERROR_ROOT(
"PETSc linear solver requested for non-PETSc build.");
1119 SLIC_ERROR_ROOT(
"Linear solver type not recognized.");
1126 iter_lin_solver->SetPrintLevel(linear_opts.
print_level);
1128 if (preconditioner) {
1129 iter_lin_solver->SetPreconditioner(*preconditioner);
1132 return {std::move(iter_lin_solver), std::move(preconditioner)};
1135 #ifdef MFEM_USE_AMGX
1136 std::unique_ptr<mfem::AmgXSolver> buildAMGX(
const AMGXOptions& options,
const MPI_Comm comm)
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";
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;
1160 static const auto solver_names = []() {
1161 std::unordered_map<AMGXSolver, std::string> names;
1179 options_node[
"solver/solver"] = solver_names.at(options.solver);
1180 options_node[
"solver/smoother"] = solver_names.at(options.smoother);
1183 amgx->ReadParameters(options_node.to_json(), mfem::AmgXSolver::INTERNAL);
1184 amgx->InitExclusiveGPU(comm);
1192 std::unique_ptr<mfem::Solver> preconditioner_solver;
1198 auto amg_preconditioner = std::make_unique<mfem::HypreBoomerAMG>();
1199 amg_preconditioner->SetPrintLevel(print_level);
1200 preconditioner_solver = std::move(amg_preconditioner);
1202 auto jac_preconditioner = std::make_unique<mfem::HypreSmoother>();
1203 jac_preconditioner->SetType(mfem::HypreSmoother::Type::Jacobi);
1204 preconditioner_solver = std::move(jac_preconditioner);
1206 auto jacl1_preconditioner = std::make_unique<mfem::HypreSmoother>();
1207 jacl1_preconditioner->SetType(mfem::HypreSmoother::Type::l1Jacobi);
1208 preconditioner_solver = std::move(jacl1_preconditioner);
1210 auto gs_preconditioner = std::make_unique<mfem::HypreSmoother>();
1211 gs_preconditioner->SetType(mfem::HypreSmoother::Type::GS);
1212 preconditioner_solver = std::move(gs_preconditioner);
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);
1219 #ifdef MFEM_USE_AMGX
1220 preconditioner_solver = buildAMGX(linear_opts.
amgx_options, comm);
1222 SLIC_ERROR_ROOT(
"AMGX requested in non-GPU build");
1225 #ifdef SERAC_USE_PETSC
1226 preconditioner_solver = mfem_ext::buildPetscPreconditioner(linear_opts.
petsc_preconditioner, comm);
1228 SLIC_ERROR_ROOT(
"PETSc preconditioner requested in non-PETSc build");
1231 SLIC_ERROR_ROOT_IF(preconditioner !=
Preconditioner::None,
"Unknown preconditioner type requested");
1234 return preconditioner_solver;
1239 auto& linear_container = container.addStruct(
"linear",
"Linear Equation Solver Parameters");
1240 linear_container.required().registerVerifier([](
const axom::inlet::Container& container_to_verify) {
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;
1250 linear_container.addString(
"type",
"The type of solver parameters to use (iterative|direct)")
1252 .validValues({
"iterative",
"direct"});
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");
1264 auto& direct_container = linear_container.addStruct(
"direct_options",
"Direct solver parameters");
1265 direct_container.addInt(
"print_level",
"Linear print level.").defaultValue(0);
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");
1285 std::string
type = base[
"type"];
1287 if (
type ==
"direct") {
1289 options.
print_level = base[
"direct_options/print_level"];
1293 auto config = base[
"iterative_options"];
1298 std::string solver_type = config[
"solver_type"];
1299 if (solver_type ==
"gmres") {
1301 }
else if (solver_type ==
"cg") {
1304 std::string msg = axom::fmt::format(
"Unknown Linear solver type given: '{0}'", solver_type);
1305 SLIC_ERROR_ROOT(msg);
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") {
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"];
1329 std::string msg = axom::fmt::format(
"Unknown preconditioner type given: '{0}'", prec_type);
1330 SLIC_ERROR_ROOT(msg);
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") {
1353 SLIC_ERROR_ROOT(axom::fmt::format(
"Unknown nonlinear solver type given: '{0}'", solver_type));
1366 std::move(linear_solver), std::move(preconditioner));
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.
void printTrustRegionInfo(double realObjective, double modelObjective, size_t cgIters, double trSize, bool willAccept)
trust region printing utility function
@ PetscNewtonBacktracking
@ PetscNewtonCriticalPoint
@ KINBacktrackingLineSearch
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.
SERAC_HOST_DEVICE auto sqrt(dual< gradient_type > x)
implementation of square root for dual numbers
constexpr SERAC_HOST_DEVICE int size(const tensor< T, n... > &)
returns the total number of stored values in a tensor
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.
constexpr SERAC_HOST_DEVICE auto type(const tuple< T... > &values)
a function intended to be used for extracting the ith type from a tuple.
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
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.