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;
102 if (
norm == 0.0)
return;
105 mfem::out <<
"Newton iteration " << std::setw(3) << 0 <<
" : ||r|| = " << std::setw(13) <<
norm <<
"\n";
108 norm_goal =
std::max(rel_tol * initial_norm, abs_tol);
109 prec->iterative_mode =
false;
113 MFEM_ASSERT(mfem::IsFinite(
norm),
"norm = " <<
norm);
115 mfem::out <<
"Newton iteration " << std::setw(3) << it <<
" : ||r|| = " << std::setw(13) <<
norm;
117 mfem::out <<
", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ?
norm / initial_norm :
norm);
123 mfem::out <<
"Initial residual for Newton iteration is undefined/nan.\n";
124 mfem::out <<
"Newton: No convergence!\n";
131 }
else if (it >= max_iter) {
136 real_t norm_nm1 =
norm;
143 x0.SetSize(x.Size());
147 real_t stepScale = 1.0;
148 add(
x0, -stepScale, c, x);
152 static constexpr real_t reduction = 0.5;
154 const double sufficientDecreaseParam = 0.0;
155 const double cMagnitudeInR = sufficientDecreaseParam != 0.0 ?
std::abs(Dot(c, r)) / norm_nm1 : 0.0;
157 auto is_improved = [=](real_t currentNorm, real_t c_scale) {
158 return currentNorm < norm_nm1 - sufficientDecreaseParam * c_scale * cMagnitudeInR;
164 for (; !is_improved(
norm, stepScale) && ls_iter < max_ls_iters; ++ls_iter, ++ls_iter_sum) {
165 stepScale *= reduction;
166 add(
x0, -stepScale, c, x);
171 if (max_ls_iters > 0 && ls_iter == max_ls_iters && !is_improved(
norm, stepScale)) {
173 add(
x0, stepScale, c, x);
177 for (; !is_improved(
norm, stepScale) && ls_iter < max_ls_iters; ++ls_iter, ++ls_iter_sum) {
178 stepScale *= reduction;
179 add(
x0, stepScale, c, x);
185 if (ls_iter == max_ls_iters && !is_improved(
norm, stepScale)) {
187 stepScale *= reduction;
188 add(
x0, -stepScale, c, x);
195 mfem::out <<
"Number of line search steps taken = " << ls_iter_sum << std::endl;
197 if (
print_level >= 2 && (ls_iter_sum == 2 * max_ls_iters + 1)) {
198 mfem::out <<
"The maximum number of line search cut back have occurred, the resulting residual may not have "
209 mfem::out <<
"Newton iteration " << std::setw(3) << final_iter <<
" : ||r|| = " << std::setw(13) <<
norm <<
'\n';
212 mfem::out <<
"Newton: No convergence!\n";
300 void printTrustRegionInfo(
double realObjective,
double modelObjective,
size_t cgIters,
double trSize,
bool willAccept)
302 mfem::out <<
"real energy = " << std::setw(13) << realObjective <<
", model energy = " << std::setw(13)
303 << modelObjective <<
", cg iter = " << std::setw(7) << cgIters <<
", next tr size = " << std::setw(8)
304 << trSize <<
", accepting = " << willAccept << std::endl;
325 mutable std::vector<std::shared_ptr<mfem::Vector>>
left_mosts;
367 double deltadelta_m_zz = delta * delta - zz;
368 if (deltadelta_m_zz == 0)
return;
369 double tau = (
std::sqrt(deltadelta_m_zz * dd + zd * zd) - zd) / dd;
374 template <
typename HessVecFunc>
376 [[maybe_unused]]
const std::vector<const mfem::Vector*> ds,
377 [[maybe_unused]]
const std::vector<const mfem::Vector*> Hds,
378 [[maybe_unused]]
const mfem::Vector& g, [[maybe_unused]]
double delta,
379 [[maybe_unused]]
int num_leftmost)
const
381 #ifdef SMITH_USE_SLEPC
385 std::vector<const mfem::Vector*> directions;
387 directions.emplace_back(d);
390 directions.emplace_back(left.get());
393 std::vector<const mfem::Vector*> H_directions;
394 for (
auto& Hd : Hds) {
395 H_directions.emplace_back(Hd);
398 H_directions.emplace_back(H_left.get());
402 std::tie(directions, H_directions) = removeDependentDirections(directions, H_directions);
403 }
catch (
const std::exception& e) {
405 mfem::out <<
"remove dependent directions failed with " << e.what() << std::endl;
414 std::vector<std::shared_ptr<mfem::Vector>> leftvecs;
415 std::vector<double> leftvals;
416 double energy_change;
419 std::tie(sol, leftvecs, leftvals, energy_change) =
420 solveSubspaceProblem(directions, H_directions, b, delta, num_leftmost);
421 }
catch (
const std::exception& e) {
423 mfem::out <<
"subspace solve failed with " << e.what() << std::endl;
429 for (
auto& lv : leftvecs) {
434 double subspace_energy =
computeEnergy(g, hess_vec_func, sol);
437 double leftval = leftvals.size() ? leftvals[0] : 1.0;
438 mfem::out <<
"Energy using subspace solver from: " << base_energy <<
", to: " << subspace_energy <<
" / "
439 << energy_change <<
". Min eig: " << leftval << std::endl;
442 if (subspace_energy < base_energy) {
452 double dd = yy - 2 * zy + zz;
454 double tau = (
std::sqrt((trSize * trSize - zz) * dd + zd * zd) - zd) / dd;
460 void doglegStep(
const mfem::Vector& cp,
const mfem::Vector& newtonP,
double trSize, mfem::Vector& s)
const
464 double cc = Dot(cp, cp);
465 double nn = Dot(newtonP, newtonP);
466 double tt = trSize * trSize;
471 }
else if (cc > nn) {
473 mfem::out <<
"cp outside newton, preconditioner likely inaccurate\n";
476 }
else if (nn > tt) {
478 double cn = Dot(cp, newtonP);
486 template <
typename HessVecFunc>
487 double computeEnergy(
const mfem::Vector& r_local,
const HessVecFunc& H,
const mfem::Vector& z)
const
490 double rz = Dot(r_local, z);
491 mfem::Vector tmp(r_local);
494 return rz + 0.5 * Dot(z, tmp);
498 template <
typename HessVecFunc,
typename PrecondFunc>
511 auto& Pr = results.
Pr;
512 auto& Hd = results.
H_d;
514 const double cg_tol_squared = settings.
cg_tol * settings.
cg_tol;
518 mfem::out <<
"Trust region solution state within tolerance on first iteration."
533 double rPr = Dot(rCurrent, Pr);
535 double dd = Dot(d, d);
541 if (Dot(d, rCurrent) > 0) {
543 results.
interior_status = TrustRegionResults::Status::NonDescentDirection;
546 hess_vec_func(d, Hd);
547 const double curvature = Dot(d, Hd);
548 const double alphaCg = curvature != 0.0 ? rPr / curvature : 0.0;
553 add(z, alphaCg, d, zPred);
554 double zzNp1 = Dot(zPred, zPred);
556 const bool go_to_boundary = curvature <= 0 || zzNp1 >= trSize * trSize;
557 if (go_to_boundary) {
559 if (curvature <= 0) {
560 results.
interior_status = TrustRegionResults::Status::NegativeCurvature;
569 if (results.
interior_status == TrustRegionResults::Status::NonDescentDirection) {
571 mfem::out <<
"Found a non descent direction\n";
576 add(rCurrent, alphaCg, Hd, rCurrent);
579 double rPrNp1 = Dot(rCurrent, Pr);
581 if (Dot(rCurrent, rCurrent) <= cg_tol_squared && cgIter >= settings.
min_cg_iterations) {
585 double beta = rPrNp1 / rPr;
587 add(-1.0, Pr, beta, d, d);
601 grad = &oper->GetGradient(x);
603 auto* grad_blocked =
dynamic_cast<mfem::BlockOperator*
>(grad);
618 void hessVec(
const mfem::Vector& x_, mfem::Vector& v_)
const
626 void precond(
const mfem::Vector& x_, mfem::Vector& v_)
const
634 void Mult(
const mfem::Vector&, mfem::Vector& X)
const
636 MFEM_ASSERT(oper != NULL,
"the Operator is not set (use SetOperator).");
637 MFEM_ASSERT(prec != NULL,
"the Solver is not set (use SetSolver).");
642 using real_t = mfem::real_t;
650 real_t
norm, norm_goal = 0.0;
652 if (
norm == 0.0)
return;
654 norm_goal =
std::max(rel_tol * initial_norm, abs_tol);
657 mfem::out <<
"TrustRegion iteration " << std::setw(3) << 0 <<
" : ||r|| = " << std::setw(13) <<
norm <<
"\n";
660 prec->iterative_mode =
false;
675 settings.
cg_tol = 0.5 * norm_goal;
682 size_t cumulative_cg_iters_from_last_precond_update = 0;
686 MFEM_ASSERT(mfem::IsFinite(
norm),
"norm = " <<
norm);
688 mfem::out <<
"TrustRegion iteration " << std::setw(3) << it <<
" : ||r|| = " << std::setw(13) <<
norm;
690 mfem::out <<
", ||r||/||r_0|| = " << std::setw(13) << (initial_norm != 0.0 ?
norm / initial_norm :
norm);
691 mfem::out <<
", x_incr = " << std::setw(13) << trResults.
d.Norml2();
693 mfem::out <<
", norm goal = " << std::setw(13) << norm_goal;
699 mfem::out <<
"Initial residual for trust-region iteration is undefined/nan." << std::endl;
700 mfem::out <<
"TrustRegion: No convergence!\n";
707 }
else if (it >= max_iter) {
717 cumulative_cg_iters_from_last_precond_update = 0;
720 auto hess_vec_func = [&](
const mfem::Vector& x_, mfem::Vector& v_) {
hessVec(x_, v_); };
721 auto precond_func = [&](
const mfem::Vector& x_, mfem::Vector& v_) {
precond(x_, v_); };
723 double cauchyPointNormSquared = tr_size * tr_size;
726 hess_vec_func(r, trResults.
H_d);
727 const double gKg = Dot(r, trResults.
H_d);
729 const double alphaCp = -Dot(r, r) / gKg;
733 const double alphaTr = -tr_size /
std::sqrt(Dot(r, r));
736 mfem::out <<
"Negative curvature un-preconditioned cauchy point direction found."
741 if (cauchyPointNormSquared >= tr_size * tr_size) {
743 mfem::out <<
"Un-preconditioned gradient cauchy point outside trust region, step size = "
744 <<
std::sqrt(cauchyPointNormSquared) <<
"\n";
757 bool have_computed_Hvs =
false;
759 int lineSearchIter = 0;
765 bool use_with_option1 =
766 (subspace_option >= 1) && (trResults.
interior_status == TrustRegionResults::Status::NonDescentDirection ||
767 trResults.
interior_status == TrustRegionResults::Status::NegativeCurvature ||
768 ((Norm(trResults.
d) > (1.0 - 1.0e-6) * tr_size) && lineSearchIter > 1));
769 bool use_with_option2 = (subspace_option >= 2) && (Norm(trResults.
d) > (1.0 - 1.0e-6) * tr_size);
770 bool use_with_option3 = (subspace_option >= 3);
772 if (use_with_option1 || use_with_option2 || use_with_option3) {
773 if (!have_computed_Hvs) {
774 have_computed_Hvs =
true;
775 hess_vec_func(trResults.
z, trResults.
H_z);
782 H_left_mosts.emplace_back(std::make_shared<mfem::Vector>(*left));
786 std::vector<const mfem::Vector*> ds{&trResults.
z, &trResults.
d_old, &trResults.
cauchy_point};
791 static constexpr
double roundOffTol = 0.0;
793 hess_vec_func(trResults.
d, trResults.
H_d);
794 double dHd = Dot(trResults.
d, trResults.
H_d);
795 double modelObjective = Dot(r, trResults.
d) + 0.5 * dHd - roundOffTol;
803 double obj1 = 0.5 * (Dot(r, trResults.
d) + Dot(
r_pred, trResults.
d)) - roundOffTol;
804 realObjective = obj1;
805 }
catch (
const std::exception&) {
810 if (normPred <= norm_goal) {
811 trResults.
d_old = trResults.
d;
823 double modelImprove = -modelObjective;
824 double realImprove = -realObjective;
826 double rho = realImprove / modelImprove;
827 if (modelObjective > 0) {
829 mfem::out <<
"Found a positive model objective increase. Debug if you see this.\n";
831 rho = realImprove / -modelImprove;
837 if (!(rho >= settings.
eta2) ||
838 rho > settings.
eta4) {
839 tr_size *= settings.
t1;
840 }
else if ((rho > settings.
eta3 && rho <= settings.
eta4 &&
841 trResults.
interior_status == TrustRegionResults::Status::OnBoundary) ||
842 (rho > 0.95 && rho < 1.05 &&
844 TrustRegionResults::Status::NegativeCurvature)) {
846 tr_size *= settings.
t2;
853 bool willAccept = rho >= settings.
eta1 && rho <= settings.
eta4;
862 trResults.
d_old = trResults.
d;
875 mfem::out <<
"TrustRegion iteration " << std::setw(3) << final_iter <<
" : ||r|| = " << std::setw(13) <<
norm
879 mfem::out <<
"TrustRegion: No convergence!\n";
896 lin_solver_ = std::move(lin_solver);
902 std::unique_ptr<mfem::Solver> linear_solver,
903 std::unique_ptr<mfem::Solver> preconditioner)
905 SLIC_ERROR_ROOT_IF(!nonlinear_solver,
"Nonlinear solvers must be given to construct an EquationSolver");
906 SLIC_ERROR_ROOT_IF(!linear_solver,
"Linear solvers must be given to construct an EquationSolver");
908 nonlin_solver_ = std::move(nonlinear_solver);
909 lin_solver_ = std::move(linear_solver);
915 nonlin_solver_->SetOperator(op);
918 if (!nonlin_solver_set_solver_called_) {
920 nonlin_solver_set_solver_called_ =
true;
926 mfem::Vector
zero(x);
930 nonlin_solver_->Mult(
zero, x);
935 SLIC_ERROR_ROOT_IF(!superlu_mat_,
"Operator must be set prior to solving with SuperLU");
938 superlu_solver_.Mult(input, output);
943 int row_blocks = block_operator.NumRowBlocks();
944 int col_blocks = block_operator.NumColBlocks();
946 SLIC_ERROR_ROOT_IF(row_blocks != col_blocks,
"Attempted to use a direct solver on a non-square block system.");
948 mfem::Array2D<const mfem::HypreParMatrix*> hypre_blocks(row_blocks, col_blocks);
950 for (
int i = 0; i < row_blocks; ++i) {
951 for (
int j = 0; j < col_blocks; ++j) {
953 if (!block_operator.IsZeroBlock(i, j)) {
954 auto* hypre_block =
dynamic_cast<const mfem::HypreParMatrix*
>(&block_operator.GetBlock(i, j));
955 SLIC_ERROR_ROOT_IF(!hypre_block,
956 "Trying to use SuperLU on a block operator that does not contain HypreParMatrix blocks.");
958 hypre_blocks(i, j) = hypre_block;
960 hypre_blocks(i, j) =
nullptr;
966 return std::unique_ptr<mfem::HypreParMatrix>(mfem::HypreParMatrixFromBlocks(hypre_blocks));
972 auto* block_operator =
dynamic_cast<const mfem::BlockOperator*
>(&op);
975 if (block_operator) {
978 superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*monolithic_mat);
981 auto* matrix =
dynamic_cast<const mfem::HypreParMatrix*
>(&op);
983 SLIC_ERROR_ROOT_IF(!matrix,
"Matrix must be an assembled HypreParMatrix for use with SuperLU");
985 superlu_mat_ = std::make_unique<mfem::SuperLURowLocMatrix>(*matrix);
988 superlu_solver_.SetOperator(*superlu_mat_);
991 #ifdef MFEM_USE_STRUMPACK
993 void StrumpackSolver::Mult(
const mfem::Vector& input, mfem::Vector& output)
const
995 SLIC_ERROR_ROOT_IF(!strumpack_mat_,
"Operator must be set prior to solving with Strumpack");
998 strumpack_solver_.Mult(input, output);
1001 void StrumpackSolver::SetOperator(
const mfem::Operator& op)
1004 auto* block_operator =
dynamic_cast<const mfem::BlockOperator*
>(&op);
1007 if (block_operator) {
1010 strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*monolithic_mat);
1013 auto* matrix =
dynamic_cast<const mfem::HypreParMatrix*
>(&op);
1015 SLIC_ERROR_ROOT_IF(!matrix,
"Matrix must be an assembled HypreParMatrix for use with Strumpack");
1017 strumpack_mat_ = std::make_unique<mfem::STRUMPACKRowLocMatrix>(*matrix);
1019 height = op.Height();
1021 strumpack_solver_.SetOperator(*strumpack_mat_);
1030 std::unique_ptr<mfem::NewtonSolver> nonlinear_solver;
1034 SLIC_ERROR_ROOT_IF(nonlinear_opts.
min_iterations != 0,
"Newton's method does not support nonzero min_iterations");
1035 nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts);
1038 SLIC_ERROR_ROOT_IF(nonlinear_opts.
min_iterations != 0,
"LBFGS does not support nonzero min_iterations");
1039 nonlinear_solver = std::make_unique<mfem::LBFGSSolver>(comm);
1041 nonlinear_solver = std::make_unique<NewtonSolver>(comm, nonlinear_opts);
1043 nonlinear_solver = std::make_unique<TrustRegion>(comm, nonlinear_opts, linear_opts, prec);
1044 #ifdef SMITH_USE_PETSC
1046 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1048 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1050 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1052 nonlinear_solver = std::make_unique<mfem_ext::PetscNewtonSolver>(comm, nonlinear_opts);
1057 #ifdef SMITH_USE_SUNDIALS
1059 SLIC_ERROR_ROOT_IF(nonlinear_opts.
min_iterations != 0,
"kinsol solvers do not support min_iterations");
1061 int kinsol_strat = KIN_NONE;
1065 kinsol_strat = KIN_NONE;
1068 kinsol_strat = KIN_LINESEARCH;
1071 kinsol_strat = KIN_PICARD;
1074 kinsol_strat = KIN_NONE;
1075 SLIC_ERROR_ROOT(
"Unknown KINSOL nonlinear solver type given.");
1077 auto kinsol_solver = std::make_unique<mfem::KINSolver>(comm, kinsol_strat,
true);
1078 nonlinear_solver = std::move(kinsol_solver);
1080 SLIC_ERROR_ROOT(
"KINSOL was not enabled when MFEM was built");
1084 nonlinear_solver->SetRelTol(nonlinear_opts.
relative_tol);
1085 nonlinear_solver->SetAbsTol(nonlinear_opts.
absolute_tol);
1087 nonlinear_solver->SetPrintLevel(nonlinear_opts.
print_level);
1092 nonlinear_solver->iterative_mode =
true;
1094 return nonlinear_solver;
1103 auto lin_solver = std::make_unique<SuperLUSolver>(linear_opts.
print_level, comm);
1104 return {std::move(lin_solver), std::move(preconditioner)};
1107 #ifdef MFEM_USE_STRUMPACK
1110 auto lin_solver = std::make_unique<StrumpackSolver>(linear_opts.
print_level, comm);
1111 return {std::move(lin_solver), std::move(preconditioner)};
1116 std::unique_ptr<mfem::IterativeSolver> iter_lin_solver;
1120 iter_lin_solver = std::make_unique<mfem::CGSolver>(comm);
1123 iter_lin_solver = std::make_unique<mfem::GMRESSolver>(comm);
1125 #ifdef SMITH_USE_PETSC
1127 iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPCG, std::string());
1130 iter_lin_solver = std::make_unique<smith::mfem_ext::PetscKSPSolver>(comm, KSPGMRES, std::string());
1135 SLIC_ERROR_ROOT(
"PETSc linear solver requested for non-PETSc build.");
1140 SLIC_ERROR_ROOT(
"Linear solver type not recognized.");
1147 iter_lin_solver->SetPrintLevel(linear_opts.
print_level);
1149 if (preconditioner) {
1150 iter_lin_solver->SetPreconditioner(*preconditioner);
1153 return {std::move(iter_lin_solver), std::move(preconditioner)};
1156 #ifdef MFEM_USE_AMGX
1157 std::unique_ptr<mfem::AmgXSolver> buildAMGX(
const AMGXOptions& options,
const MPI_Comm comm)
1159 auto amgx = std::make_unique<mfem::AmgXSolver>();
1160 conduit::Node options_node;
1161 options_node[
"config_version"] = 2;
1162 auto& solver_options = options_node[
"solver"];
1163 solver_options[
"solver"] =
"AMG";
1164 solver_options[
"presweeps"] = 1;
1165 solver_options[
"postsweeps"] = 2;
1166 solver_options[
"interpolator"] =
"D2";
1167 solver_options[
"max_iters"] = 2;
1168 solver_options[
"convergence"] =
"ABSOLUTE";
1169 solver_options[
"cycle"] =
"V";
1171 if (options.verbose) {
1172 options_node[
"solver/obtain_timings"] = 1;
1173 options_node[
"solver/monitor_residual"] = 1;
1174 options_node[
"solver/print_solve_stats"] = 1;
1181 static const auto solver_names = []() {
1182 std::unordered_map<AMGXSolver, std::string> names;
1200 options_node[
"solver/solver"] = solver_names.at(options.solver);
1201 options_node[
"solver/smoother"] = solver_names.at(options.smoother);
1204 amgx->ReadParameters(options_node.to_json(), mfem::AmgXSolver::INTERNAL);
1205 amgx->InitExclusiveGPU(comm);
1213 std::unique_ptr<mfem::Solver> preconditioner_solver;
1219 auto amg_preconditioner = std::make_unique<mfem::HypreBoomerAMG>();
1220 amg_preconditioner->SetPrintLevel(preconditioner_print_level);
1221 preconditioner_solver = std::move(amg_preconditioner);
1223 auto jac_preconditioner = std::make_unique<mfem::HypreSmoother>();
1224 jac_preconditioner->SetType(mfem::HypreSmoother::Type::Jacobi);
1225 preconditioner_solver = std::move(jac_preconditioner);
1227 auto jacl1_preconditioner = std::make_unique<mfem::HypreSmoother>();
1228 jacl1_preconditioner->SetType(mfem::HypreSmoother::Type::l1Jacobi);
1229 preconditioner_solver = std::move(jacl1_preconditioner);
1231 auto gs_preconditioner = std::make_unique<mfem::HypreSmoother>();
1232 gs_preconditioner->SetType(mfem::HypreSmoother::Type::GS);
1233 preconditioner_solver = std::move(gs_preconditioner);
1235 auto ilu_preconditioner = std::make_unique<mfem::HypreILU>();
1236 ilu_preconditioner->SetLevelOfFill(1);
1237 ilu_preconditioner->SetPrintLevel(preconditioner_print_level);
1238 preconditioner_solver = std::move(ilu_preconditioner);
1240 #ifdef MFEM_USE_AMGX
1241 preconditioner_solver = buildAMGX(linear_opts.
amgx_options, comm);
1243 SLIC_ERROR_ROOT(
"AMGX requested in non-GPU build");
1246 #ifdef SMITH_USE_PETSC
1247 preconditioner_solver = mfem_ext::buildPetscPreconditioner(linear_opts.
petsc_preconditioner, comm);
1249 SLIC_ERROR_ROOT(
"PETSc preconditioner requested in non-PETSc build");
1252 auto amgfcontact_preconditioner = std::make_unique<mfem::AMGFSolver>();
1254 amgfcontact_preconditioner->GetAMG().SetPrintLevel(preconditioner_print_level);
1255 amgfcontact_preconditioner->GetAMG().SetSystemsOptions(amgfcontact_opts.dim_systems_options);
1256 amgfcontact_preconditioner->GetAMG().SetRelaxType(amgfcontact_opts.relax_type);
1257 preconditioner_solver = std::move(amgfcontact_preconditioner);
1259 SLIC_ERROR_ROOT_IF(preconditioner !=
Preconditioner::None,
"Unknown preconditioner type requested");
1262 return preconditioner_solver;
1267 auto& linear_container = container.addStruct(
"linear",
"Linear Equation Solver Parameters");
1268 linear_container.required().registerVerifier([](
const axom::inlet::Container& container_to_verify) {
1270 const bool is_iterative = (container_to_verify[
"type"].get<std::string>() ==
"iterative") &&
1271 container_to_verify.contains(
"iterative_options");
1272 const bool is_direct =
1273 (container_to_verify[
"type"].get<std::string>() ==
"direct") && container_to_verify.contains(
"direct_options");
1274 return is_iterative || is_direct;
1278 linear_container.addString(
"type",
"The type of solver parameters to use (iterative|direct)")
1280 .validValues({
"iterative",
"direct"});
1282 auto& iterative_container = linear_container.addStruct(
"iterative_options",
"Iterative solver parameters");
1283 iterative_container.addDouble(
"rel_tol",
"Relative tolerance for the linear solve.").defaultValue(1.0e-6);
1284 iterative_container.addDouble(
"abs_tol",
"Absolute tolerance for the linear solve.").defaultValue(1.0e-8);
1285 iterative_container.addInt(
"max_iter",
"Maximum iterations for the linear solve.").defaultValue(5000);
1286 iterative_container.addInt(
"print_level",
"Linear print level.").defaultValue(0);
1287 iterative_container.addString(
"solver_type",
"Solver type (gmres|minres|cg).").defaultValue(
"gmres");
1288 iterative_container.addString(
"prec_type",
"Preconditioner type (JacobiSmoother|L1JacobiSmoother|AMG|ILU|Petsc).")
1289 .defaultValue(
"JacobiSmoother");
1290 iterative_container.addString(
"petsc_prec_type",
"Type of PETSc preconditioner to use.").defaultValue(
"jacobi");
1292 auto& direct_container = linear_container.addStruct(
"direct_options",
"Direct solver parameters");
1293 direct_container.addInt(
"print_level",
"Linear print level.").defaultValue(0);
1296 auto& nonlinear_container = container.addStruct(
"nonlinear",
"Newton Equation Solver Parameters").required(
false);
1297 nonlinear_container.addDouble(
"rel_tol",
"Relative tolerance for the Newton solve.").defaultValue(1.0e-2);
1298 nonlinear_container.addDouble(
"abs_tol",
"Absolute tolerance for the Newton solve.").defaultValue(1.0e-4);
1299 nonlinear_container.addInt(
"max_iter",
"Maximum iterations for the Newton solve.").defaultValue(500);
1300 nonlinear_container.addInt(
"print_level",
"Nonlinear print level.").defaultValue(0);
1301 nonlinear_container.addString(
"solver_type",
"Solver type (Newton|KINFullStep|KINLineSearch)").defaultValue(
"Newton");
1313 std::string
type = base[
"type"];
1315 if (
type ==
"direct") {
1317 options.
print_level = base[
"direct_options/print_level"];
1321 auto config = base[
"iterative_options"];
1326 std::string solver_type = config[
"solver_type"];
1327 if (solver_type ==
"gmres") {
1329 }
else if (solver_type ==
"cg") {
1332 std::string msg = std::format(
"Unknown Linear solver type given: '{0}'", solver_type);
1333 SLIC_ERROR_ROOT(msg);
1335 const std::string prec_type = config[
"prec_type"];
1336 if (prec_type ==
"JacobiSmoother") {
1338 }
else if (prec_type ==
"L1JacobiSmoother") {
1340 }
else if (prec_type ==
"HypreAMG") {
1342 }
else if (prec_type ==
"ILU") {
1344 #ifdef MFEM_USE_AMGX
1345 }
else if (prec_type ==
"AMGX") {
1348 }
else if (prec_type ==
"GaussSeidel") {
1350 #ifdef SMITH_USE_PETSC
1351 }
else if (prec_type ==
"Petsc") {
1352 const std::string petsc_prec = config[
"petsc_prec_type"];
1356 }
else if (prec_type ==
"AMGFContact") {
1359 std::string msg = std::format(
"Unknown preconditioner type given: '{0}'", prec_type);
1360 SLIC_ERROR_ROOT(msg);
1373 const std::string solver_type = base[
"solver_type"];
1374 if (solver_type ==
"Newton") {
1376 }
else if (solver_type ==
"KINFullStep") {
1378 }
else if (solver_type ==
"KINLineSearch") {
1380 }
else if (solver_type ==
"KINPicard") {
1383 SLIC_ERROR_ROOT(std::format(
"Unknown nonlinear solver type given: '{0}'", solver_type));
1396 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(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 preconditioner_print_level
Debugging print level for the preconditioner.
int max_iterations
Maximum number of iterations.
double relative_tol
Relative tolerance.
int print_level
Debugging print level for the linear solver.
AMGFContactOptions amgfcontact_options
AMGFContact Options, used for Preconditioner::AMGFContact.
Preconditioner preconditioner
PreconditionerOptions selection.
PetscPCType petsc_preconditioner
PETSc preconditioner type.
double absolute_tol
Absolute tolerance.
AMGXOptions amgx_options
AMGX Options, used for Preconditioner::AMGX.
LinearSolver linear_solver
Linear solver selection.
Nonlinear solution scheme parameters.
int min_iterations
Minimum number of iterations.
int max_iterations
Maximum number of iterations.
double trust_region_scaling
Scaling for the initial trust region size.
int max_line_search_iterations
Maximum line search cutbacks.
double relative_tol
Relative tolerance.
int print_level
Debug print level.
int num_leftmost
Number of extra leftmost eigenvector to be stored between solves.
double absolute_tol
Absolute tolerance.
NonlinearSolver nonlin_solver
Nonlinear solver selection.
SubSpaceOptions subspace_option
Option for how when the subspace solver should be utilized within trust-region solver.
bool force_monolithic
Should the gradient be converted to a monolithic matrix.
Internal structure for storing trust region stateful data.
size_t cg_iterations_count
iteration counter
mfem::Vector cauchy_point
cauchy point
mfem::Vector H_z
action of hessian on current step z
mfem::Vector z
step direction
mfem::Vector d
incrementalCG direction
mfem::Vector H_cauchy_point
action of hessian on direction of cauchy point
mfem::Vector H_d_old
action of hessian on previous step z_old
void reset()
resets trust region results for a new outer iteration
mfem::Vector H_d
action of hessian on direction d
TrustRegionResults(int size)
Constructor takes the size of the solution vector.
mfem::Vector Pr
preconditioned residual
Status
enumerates the possible final status of the trust region steps
Status interior_status
specifies if step is interior, exterior, negative curvature, etc.
mfem::Vector d_old
old step direction
Internal structure for storing trust region settings.
double eta2
non-ideal energy drop ratio. trust region decreases if energy drop is worse than this.
size_t max_cg_iterations
max cg iters should be around # of system dofs
double eta1
worse case energy drop ratio. trust region accepted if energy drop is better than this.
double eta3
ideal energy drop ratio. trust region increases if energy drop is better than this.
size_t min_cg_iterations
min cg iters
size_t max_cumulative_iteration
max cumulative iterations
double min_tr_size
minimum trust region size
double t1
trust region decrease factor
double eta4
parameter limiting how fast the energy can drop relative to the prediction (in case the energy surrog...
double t2
trust region increase factor
A sentinel struct for eliding no-op tensor operations.