19 #include "smith/smith_config.hpp"
21 #include "smith/numerics/trust_region_solver.hpp"
30 mutable mfem::Vector
x0;
51 double evaluateNorm(
const mfem::Vector& x, mfem::Vector& rOut)
const
57 normEval = Norm(rOut);
58 }
catch (
const std::exception&) {
68 grad = &oper->GetGradient(x);
70 auto* grad_blocked =
dynamic_cast<mfem::BlockOperator*
>(grad);
79 prec->SetOperator(*grad);
90 void Mult(
const mfem::Vector&, mfem::Vector& x)
const
92 MFEM_ASSERT(oper != NULL,
"the Operator is not set (use SetOperator).");
93 MFEM_ASSERT(prec != NULL,
"the Solver is not set (use SetSolver).");
98 using real_t = mfem::real_t;
100 real_t
norm, norm_goal = 0;
104 mfem::out <<
"Newton iteration " << std::setw(3) << 0 <<
" : ||r|| = " << std::setw(13) <<
norm <<
"\n";
107 norm_goal =
std::max(rel_tol * initial_norm, abs_tol);
108 prec->iterative_mode =
false;
112 MFEM_ASSERT(mfem::IsFinite(
norm),
"norm = " <<
norm);
114 mfem::out <<
"Newton iteration " << std::setw(3) << it <<
" : ||r|| = " << std::setw(13) <<
norm;
116 mfem::out <<
", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ?
norm / initial_norm :
norm);
122 mfem::out <<
"Initial residual for Newton iteration is undefined/nan.\n";
123 mfem::out <<
"Newton: No convergence!\n";
130 }
else if (it >= max_iter) {
135 real_t norm_nm1 =
norm;
142 x0.SetSize(x.Size());
146 real_t stepScale = 1.0;
147 add(
x0, -stepScale, c, x);
151 static constexpr real_t reduction = 0.5;
153 const double sufficientDecreaseParam = 0.0;
154 const double cMagnitudeInR = sufficientDecreaseParam != 0.0 ?
std::abs(Dot(c, r)) / norm_nm1 : 0.0;
156 auto is_improved = [=](real_t currentNorm, real_t c_scale) {
157 return currentNorm < norm_nm1 - sufficientDecreaseParam * c_scale * cMagnitudeInR;
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);
170 if (max_ls_iters > 0 && ls_iter == max_ls_iters && !is_improved(
norm, stepScale)) {
172 add(
x0, stepScale, c, x);
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);
184 if (ls_iter == max_ls_iters && !is_improved(
norm, stepScale)) {
186 stepScale *= reduction;
187 add(
x0, -stepScale, c, x);
194 mfem::out <<
"Number of line search steps taken = " << ls_iter_sum << std::endl;
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 "
208 mfem::out <<
"Newton iteration " << std::setw(3) << final_iter <<
" : ||r|| = " << std::setw(13) <<
norm <<
'\n';
211 mfem::out <<
"Newton: No convergence!\n";
299 void printTrustRegionInfo(
double realObjective,
double modelObjective,
size_t cgIters,
double trSize,
bool willAccept)
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;
324 mutable std::vector<std::shared_ptr<mfem::Vector>>
left_mosts;
366 double deltadelta_m_zz = delta * delta - zz;
367 if (deltadelta_m_zz == 0)
return;
368 double tau = (
std::sqrt(deltadelta_m_zz * dd + zd * zd) - zd) / dd;
373 template <
typename HessVecFunc>
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
380 #ifdef SMITH_USE_SLEPC
384 std::vector<const mfem::Vector*> directions;
386 directions.emplace_back(d);
389 directions.emplace_back(left.get());
392 std::vector<const mfem::Vector*> H_directions;
393 for (
auto& Hd : Hds) {
394 H_directions.emplace_back(Hd);
397 H_directions.emplace_back(H_left.get());
401 std::tie(directions, H_directions) = removeDependentDirections(directions, H_directions);
402 }
catch (
const std::exception& e) {
404 mfem::out <<
"remove dependent directions failed with " << e.what() << std::endl;
413 std::vector<std::shared_ptr<mfem::Vector>> leftvecs;
414 std::vector<double> leftvals;
415 double energy_change;
418 std::tie(sol, leftvecs, leftvals, energy_change) =
419 solveSubspaceProblem(directions, H_directions, b, delta, num_leftmost);
420 }
catch (
const std::exception& e) {
422 mfem::out <<
"subspace solve failed with " << e.what() << std::endl;
428 for (
auto& lv : leftvecs) {
433 double subspace_energy =
computeEnergy(g, hess_vec_func, sol);
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;
441 if (subspace_energy < base_energy) {
451 double dd = yy - 2 * zy + zz;
453 double tau = (
std::sqrt((trSize * trSize - zz) * dd + zd * zd) - zd) / dd;
459 void doglegStep(
const mfem::Vector& cp,
const mfem::Vector& newtonP,
double trSize, mfem::Vector& s)
const
463 double cc = Dot(cp, cp);
464 double nn = Dot(newtonP, newtonP);
465 double tt = trSize * trSize;
470 }
else if (cc > nn) {
472 mfem::out <<
"cp outside newton, preconditioner likely inaccurate\n";
475 }
else if (nn > tt) {
477 double cn = Dot(cp, newtonP);
485 template <
typename HessVecFunc>
486 double computeEnergy(
const mfem::Vector& r_local,
const HessVecFunc& H,
const mfem::Vector& z)
const
489 double rz = Dot(r_local, z);
490 mfem::Vector tmp(r_local);
493 return rz + 0.5 * Dot(z, tmp);
497 template <
typename HessVecFunc,
typename PrecondFunc>
510 auto& Pr = results.
Pr;
511 auto& Hd = results.
H_d;
513 const double cg_tol_squared = settings.
cg_tol * settings.
cg_tol;
517 mfem::out <<
"Trust region solution state within tolerance on first iteration."
532 double rPr = Dot(rCurrent, Pr);
534 double dd = Dot(d, d);
540 if (Dot(d, rCurrent) > 0) {
542 results.
interior_status = TrustRegionResults::Status::NonDescentDirection;
545 hess_vec_func(d, Hd);
546 const double curvature = Dot(d, Hd);
547 const double alphaCg = curvature != 0.0 ? rPr / curvature : 0.0;
552 add(z, alphaCg, d, zPred);
553 double zzNp1 = Dot(zPred, zPred);
555 const bool go_to_boundary = curvature <= 0 || zzNp1 >= trSize * trSize;
556 if (go_to_boundary) {
558 if (curvature <= 0) {
559 results.
interior_status = TrustRegionResults::Status::NegativeCurvature;
568 if (results.
interior_status == TrustRegionResults::Status::NonDescentDirection) {
570 mfem::out <<
"Found a non descent direction\n";
575 add(rCurrent, alphaCg, Hd, rCurrent);
578 double rPrNp1 = Dot(rCurrent, Pr);
580 if (Dot(rCurrent, rCurrent) <= cg_tol_squared && cgIter >= settings.
min_cg_iterations) {
584 double beta = rPrNp1 / rPr;
586 add(-1.0, Pr, beta, d, d);
600 grad = &oper->GetGradient(x);
602 auto* grad_blocked =
dynamic_cast<mfem::BlockOperator*
>(grad);
617 void hessVec(
const mfem::Vector& x_, mfem::Vector& v_)
const
625 void precond(
const mfem::Vector& x_, mfem::Vector& v_)
const
633 void Mult(
const mfem::Vector&, mfem::Vector& X)
const
635 MFEM_ASSERT(oper != NULL,
"the Operator is not set (use SetOperator).");
636 MFEM_ASSERT(prec != NULL,
"the Solver is not set (use SetSolver).");
638 using real_t = mfem::real_t;
646 real_t
norm, norm_goal = 0.0;
648 norm_goal =
std::max(rel_tol * initial_norm, abs_tol);
651 mfem::out <<
"TrustRegion iteration " << std::setw(3) << 0 <<
" : ||r|| = " << std::setw(13) <<
norm <<
"\n";
654 prec->iterative_mode =
false;
669 settings.
cg_tol = 0.5 * norm_goal;
676 size_t cumulative_cg_iters_from_last_precond_update = 0;
680 MFEM_ASSERT(mfem::IsFinite(
norm),
"norm = " <<
norm);
682 mfem::out <<
"TrustRegion iteration " << std::setw(3) << it <<
" : ||r|| = " << std::setw(13) <<
norm;
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();
687 mfem::out <<
", norm goal = " << std::setw(13) << norm_goal;
693 mfem::out <<
"Initial residual for trust-region iteration is undefined/nan." << std::endl;
694 mfem::out <<
"TrustRegion: No convergence!\n";
701 }
else if (it >= max_iter) {
711 cumulative_cg_iters_from_last_precond_update = 0;
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_); };
717 double cauchyPointNormSquared = tr_size * tr_size;
720 hess_vec_func(r, trResults.
H_d);
721 const double gKg = Dot(r, trResults.
H_d);
723 const double alphaCp = -Dot(r, r) / gKg;
727 const double alphaTr = -tr_size /
std::sqrt(Dot(r, r));
730 mfem::out <<
"Negative curvature un-preconditioned cauchy point direction found."
735 if (cauchyPointNormSquared >= tr_size * tr_size) {
737 mfem::out <<
"Un-preconditioned gradient cauchy point outside trust region, step size = "
738 <<
std::sqrt(cauchyPointNormSquared) <<
"\n";
751 bool have_computed_Hvs =
false;
753 int lineSearchIter = 0;
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);
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);
776 H_left_mosts.emplace_back(std::make_shared<mfem::Vector>(*left));
780 std::vector<const mfem::Vector*> ds{&trResults.
z, &trResults.
d_old, &trResults.
cauchy_point};
785 static constexpr
double roundOffTol = 0.0;
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;
797 double obj1 = 0.5 * (Dot(r, trResults.
d) + Dot(
r_pred, trResults.
d)) - roundOffTol;
798 realObjective = obj1;
799 }
catch (
const std::exception&) {
804 if (normPred <= norm_goal) {
805 trResults.
d_old = trResults.
d;
809 if (print_options.iterations) {
817 double modelImprove = -modelObjective;
818 double realImprove = -realObjective;
820 double rho = realImprove / modelImprove;
821 if (modelObjective > 0) {
823 mfem::out <<
"Found a positive model objective increase. Debug if you see this.\n";
825 rho = realImprove / -modelImprove;
831 if (!(rho >= settings.
eta2) ||
832 rho > settings.
eta4) {
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 &&
838 TrustRegionResults::Status::NegativeCurvature)) {
840 tr_size *= settings.
t2;
847 bool willAccept = rho >= settings.
eta1 && rho <= settings.
eta4;
849 if (print_options.iterations) {
856 trResults.
d_old = trResults.
d;
869 mfem::out <<
"TrustRegion iteration " << std::setw(3) << final_iter <<
" : ||r|| = " << std::setw(13) <<
norm
873 mfem::out <<
"TrustRegion: No convergence!\n";
890 lin_solver_ = std::move(lin_solver);
896 std::unique_ptr<mfem::Solver> linear_solver,
897 std::unique_ptr<mfem::Solver> preconditioner)
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");
902 nonlin_solver_ = std::move(nonlinear_solver);
903 lin_solver_ = std::move(linear_solver);
909 nonlin_solver_->SetOperator(op);
912 if (!nonlin_solver_set_solver_called_) {
914 nonlin_solver_set_solver_called_ =
true;
920 mfem::Vector
zero(x);
924 nonlin_solver_->Mult(
zero, x);
929 SLIC_ERROR_ROOT_IF(!superlu_mat_,
"Operator must be set prior to solving with SuperLU");
932 superlu_solver_.Mult(input, output);
937 int row_blocks = block_operator.NumRowBlocks();
938 int col_blocks = block_operator.NumColBlocks();
940 SLIC_ERROR_ROOT_IF(row_blocks != col_blocks,
"Attempted to use a direct solver on a non-square block system.");
942 mfem::Array2D<const mfem::HypreParMatrix*> hypre_blocks(row_blocks, col_blocks);
944 for (
int i = 0; i < row_blocks; ++i) {
945 for (
int j = 0; j < col_blocks; ++j) {
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.");
952 hypre_blocks(i, j) = hypre_block;
954 hypre_blocks(i, j) =
nullptr;
960 return std::unique_ptr<mfem::HypreParMatrix>(mfem::HypreParMatrixFromBlocks(hypre_blocks));
966 auto* block_operator =
dynamic_cast<const mfem::BlockOperator*
>(&op);
969 if (block_operator) {
972 superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*monolithic_mat);
975 auto* matrix =
dynamic_cast<const mfem::HypreParMatrix*
>(&op);
977 SLIC_ERROR_ROOT_IF(!matrix,
"Matrix must be an assembled HypreParMatrix for use with SuperLU");
979 superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*matrix);
982 superlu_solver_.SetOperator(*superlu_mat_);
985 #ifdef MFEM_USE_STRUMPACK
987 void StrumpackSolver::Mult(
const mfem::Vector& input, mfem::Vector& output)
const
989 SLIC_ERROR_ROOT_IF(!strumpack_mat_,
"Operator must be set prior to solving with Strumpack");
992 strumpack_solver_.Mult(input, output);
995 void StrumpackSolver::SetOperator(
const mfem::Operator& op)
998 auto* block_operator =
dynamic_cast<const mfem::BlockOperator*
>(&op);
1001 if (block_operator) {
1004 strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*monolithic_mat);
1007 auto* matrix =
dynamic_cast<const mfem::HypreParMatrix*
>(&op);
1009 SLIC_ERROR_ROOT_IF(!matrix,
"Matrix must be an assembled HypreParMatrix for use with Strumpack");
1011 strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*matrix);
1014 strumpack_solver_.SetOperator(*strumpack_mat_);
1023 std::unique_ptr<mfem::NewtonSolver> nonlinear_solver;
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);
1032 "LBFGS does not support nonzero min_iterations or max_line_search_iterations");
1033 nonlinear_solver = std::make_unique<mfem::LBFGSSolver>(comm);
1035 nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts);
1037 nonlinear_solver = std::make_unique<TrustRegion>(comm, nonlinear_opts, linear_opts, prec);
1038 #ifdef SMITH_USE_PETSC
1040 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1042 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1044 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1046 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1051 #ifdef SMITH_USE_SUNDIALS
1054 "kinsol solvers do not support min_iterations or max_line_search_iterations");
1056 int kinsol_strat = KIN_NONE;
1060 kinsol_strat = KIN_NONE;
1063 kinsol_strat = KIN_LINESEARCH;
1066 kinsol_strat = KIN_PICARD;
1069 kinsol_strat = KIN_NONE;
1070 SLIC_ERROR_ROOT(
"Unknown KINSOL nonlinear solver type given.");
1072 auto kinsol_solver = std::make_unique<mfem::KINSolver>(comm, kinsol_strat,
true);
1073 nonlinear_solver = std::move(kinsol_solver);
1075 SLIC_ERROR_ROOT(
"KINSOL was not enabled when MFEM was built");
1079 nonlinear_solver->SetRelTol(nonlinear_opts.
relative_tol);
1080 nonlinear_solver->SetAbsTol(nonlinear_opts.
absolute_tol);
1082 nonlinear_solver->SetPrintLevel(nonlinear_opts.
print_level);
1087 nonlinear_solver->iterative_mode =
true;
1089 return nonlinear_solver;
1098 auto lin_solver = std::make_unique<SuperLUSolver>(linear_opts.
print_level, comm);
1099 return {std::move(lin_solver), std::move(preconditioner)};
1102 #ifdef MFEM_USE_STRUMPACK
1105 auto lin_solver = std::make_unique<StrumpackSolver>(linear_opts.
print_level, comm);
1106 return {std::move(lin_solver), std::move(preconditioner)};
1111 std::unique_ptr<mfem::IterativeSolver> iter_lin_solver;
1115 iter_lin_solver = std::make_unique<mfem::CGSolver>(comm);
1118 iter_lin_solver = std::make_unique<mfem::GMRESSolver>(comm);
1120 #ifdef SMITH_USE_PETSC
1122 iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPCG, std::string());
1125 iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPGMRES, std::string());
1130 SLIC_ERROR_ROOT(
"PETSc linear solver requested for non-PETSc build.");
1135 SLIC_ERROR_ROOT(
"Linear solver type not recognized.");
1142 iter_lin_solver->SetPrintLevel(linear_opts.
print_level);
1144 if (preconditioner) {
1145 iter_lin_solver->SetPreconditioner(*preconditioner);
1148 return {std::move(iter_lin_solver), std::move(preconditioner)};
1151 #ifdef MFEM_USE_AMGX
1152 std::unique_ptr<mfem::AmgXSolver> buildAMGX(
const AMGXOptions& options,
const MPI_Comm comm)
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";
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;
1176 static const auto solver_names = []() {
1177 std::unordered_map<AMGXSolver, std::string> names;
1195 options_node[
"solver/solver"] = solver_names.at(options.solver);
1196 options_node[
"solver/smoother"] = solver_names.at(options.smoother);
1199 amgx->ReadParameters(options_node.to_json(), mfem::AmgXSolver::INTERNAL);
1200 amgx->InitExclusiveGPU(comm);
1208 std::unique_ptr<mfem::Solver> preconditioner_solver;
1214 auto amg_preconditioner = std::make_unique<mfem::HypreBoomerAMG>();
1215 amg_preconditioner->SetPrintLevel(print_level);
1216 preconditioner_solver = std::move(amg_preconditioner);
1218 auto jac_preconditioner = std::make_unique<mfem::HypreSmoother>();
1219 jac_preconditioner->SetType(mfem::HypreSmoother::Type::Jacobi);
1220 preconditioner_solver = std::move(jac_preconditioner);
1222 auto jacl1_preconditioner = std::make_unique<mfem::HypreSmoother>();
1223 jacl1_preconditioner->SetType(mfem::HypreSmoother::Type::l1Jacobi);
1224 preconditioner_solver = std::move(jacl1_preconditioner);
1226 auto gs_preconditioner = std::make_unique<mfem::HypreSmoother>();
1227 gs_preconditioner->SetType(mfem::HypreSmoother::Type::GS);
1228 preconditioner_solver = std::move(gs_preconditioner);
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);
1235 #ifdef MFEM_USE_AMGX
1236 preconditioner_solver = buildAMGX(linear_opts.
amgx_options, comm);
1238 SLIC_ERROR_ROOT(
"AMGX requested in non-GPU build");
1241 #ifdef SMITH_USE_PETSC
1242 preconditioner_solver = mfem_ext::buildPetscPreconditioner(linear_opts.
petsc_preconditioner, comm);
1244 SLIC_ERROR_ROOT(
"PETSc preconditioner requested in non-PETSc build");
1247 SLIC_ERROR_ROOT_IF(preconditioner !=
Preconditioner::None,
"Unknown preconditioner type requested");
1250 return preconditioner_solver;
1255 auto& linear_container = container.addStruct(
"linear",
"Linear Equation Solver Parameters");
1256 linear_container.required().registerVerifier([](
const axom::inlet::Container& container_to_verify) {
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;
1266 linear_container.addString(
"type",
"The type of solver parameters to use (iterative|direct)")
1268 .validValues({
"iterative",
"direct"});
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");
1280 auto& direct_container = linear_container.addStruct(
"direct_options",
"Direct solver parameters");
1281 direct_container.addInt(
"print_level",
"Linear print level.").defaultValue(0);
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");
1301 std::string
type = base[
"type"];
1303 if (
type ==
"direct") {
1305 options.
print_level = base[
"direct_options/print_level"];
1309 auto config = base[
"iterative_options"];
1314 std::string solver_type = config[
"solver_type"];
1315 if (solver_type ==
"gmres") {
1317 }
else if (solver_type ==
"cg") {
1320 std::string msg = axom::fmt::format(
"Unknown Linear solver type given: '{0}'", solver_type);
1321 SLIC_ERROR_ROOT(msg);
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") {
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"];
1345 std::string msg = axom::fmt::format(
"Unknown preconditioner type given: '{0}'", prec_type);
1346 SLIC_ERROR_ROOT(msg);
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") {
1369 SLIC_ERROR_ROOT(axom::fmt::format(
"Unknown nonlinear solver type given: '{0}'", solver_type));
1382 std::move(linear_solver), std::move(preconditioner));
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.
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
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.
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.
constexpr SMITH_HOST_DEVICE int size(const tensor< T, n... > &)
returns the total number of stored values in a tensor
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.
@ PetscNewtonBacktracking
@ PetscNewtonCriticalPoint
@ KINBacktrackingLineSearch
Various helper functions and macros for profiling using Caliper.
#define SMITH_MARK_FUNCTION
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.