7#ifndef NEWTON_SOLVER_H_
8#define NEWTON_SOLVER_H_
17#include <AMReX_Config.H>
30template<
class Vec,
class Ops>
45 void Define ( const Vec& a_U,
46 Ops* a_ops ) override;
48 void Solve ( Vec& a_U,
52 int a_step) const override;
56 int& a_maxits )
override
63 inline void CurTime ( amrex::Real a_time )
const
83 amrex::Print() <<
"Newton linear solver: " << linsol_name <<
"\n";
197template <
class Vec,
class Ops>
204 "Newton nonlinear solver object is already defined!");
219 m_linear_solver = std::make_unique<AMReXGMRES<Vec,JacobianFunctionMF<Vec,Ops>>>();
221#ifdef AMREX_USE_PETSC
222 m_linear_solver = std::make_unique<PETScKSP<Vec,JacobianFunctionMF<Vec,Ops>>>();
224 WARPX_ABORT_WITH_MESSAGE(
"NewtonSolver::Define(): must compile with PETSc to use petsc_ksp (AMREX_USE_PETSC must be defined)");
242 std::filesystem::path
const diagnostic_dir = diagnostic_path.parent_path();
243 if (!diagnostic_dir.empty()) {
244 std::filesystem::create_directories(diagnostic_dir);
247 std::ofstream diagnostic_file{this->
m_diagnostic_file, std::ofstream::out | std::ofstream::trunc};
249 diagnostic_file <<
"#";
250 diagnostic_file <<
"[" << c++ <<
"]step()";
251 diagnostic_file <<
" ";
252 diagnostic_file <<
"[" << c++ <<
"]time(s)";
253 diagnostic_file <<
" ";
254 diagnostic_file <<
"[" << c++ <<
"]iters";
255 diagnostic_file <<
" ";
256 diagnostic_file <<
"[" << c++ <<
"]total_iters";
257 diagnostic_file <<
" ";
258 diagnostic_file <<
"[" << c++ <<
"]norm_abs";
259 diagnostic_file <<
" ";
260 diagnostic_file <<
"[" << c++ <<
"]norm_rel";
261 diagnostic_file <<
" ";
262 diagnostic_file <<
"[" << c++ <<
"]gmres_iters";
263 diagnostic_file <<
" ";
264 diagnostic_file <<
"[" << c++ <<
"]gmres_total_iters";
265 diagnostic_file <<
" ";
266 diagnostic_file <<
"[" << c++ <<
"]gmres_last_res";
267 diagnostic_file <<
"\n";
268 diagnostic_file.close();
272template <
class Vec,
class Ops>
304template <
class Vec,
class Ops>
314 "NewtonSolver::Solve() called on undefined object");
315 using namespace amrex::literals;
325 amrex::Real norm_abs = 0.;
326 amrex::Real norm0 = 1._rt;
327 amrex::Real norm_rel = 0.;
330 int linear_solver_iters = 0;
337 norm_abs =
m_F.norm2();
339 if (norm_abs > 0.) { norm0 = norm_abs; }
340 else { norm0 = 1._rt; }
342 norm_rel = norm_abs/norm0;
346 amrex::Print() <<
"Newton: iteration = " << std::setw(3) << iter <<
", norm = "
347 << std::scientific << std::setprecision(5) << norm_abs <<
" (abs.), "
348 << std::scientific << std::setprecision(5) << norm_rel <<
" (rel.)" <<
"\n";
353 amrex::Print() <<
"Newton: exiting at iteration = " << std::setw(3) << iter
354 <<
". Satisfied absolute tolerance " <<
m_atol <<
"\n";
361 amrex::Print() <<
"Newton: exiting at iteration = " << std::setw(3) << iter
362 <<
". Satisfied relative tolerance " <<
m_rtol <<
"\n";
367 if (norm_abs > 100._rt*norm0) {
368 amrex::Print() <<
"Newton: exiting at iteration = " << std::setw(3) << iter
369 <<
". SOLVER DIVERGED! relative tolerance = " <<
m_rtol <<
"\n";
370 std::stringstream convergenceMsg;
371 convergenceMsg <<
"Newton: exiting at iteration " << std::setw(3) << iter <<
372 ". SOLVER DIVERGED! absolute norm = " << norm_abs <<
373 " has increased by 100X from that after first iteration.";
388 amrex::Print() <<
"Newton: exiting at iter = " << std::setw(3) << iter
389 <<
". Maximum iteration reached: iter = " <<
m_maxits <<
"\n";
401 std::stringstream convergenceMsg;
402 convergenceMsg <<
"Newton solver failed to converge after " << iter <<
403 " iterations. Relative norm is " << norm_rel <<
404 " and the relative tolerance is " <<
m_rtol <<
405 ". Absolute norm is " << norm_abs <<
406 " and the absolute tolerance is " <<
m_atol;
416 ((a_step+1)%this->m_diagnostic_interval==0 || a_step==0)) {
417 std::ofstream diagnostic_file{this->
m_diagnostic_file, std::ofstream::out | std::ofstream::app};
418 diagnostic_file << std::setprecision(14);
419 diagnostic_file << a_step+1;
420 diagnostic_file <<
" ";
421 diagnostic_file << a_time + a_dt;
422 diagnostic_file <<
" ";
423 diagnostic_file << iter;
424 diagnostic_file <<
" ";
426 diagnostic_file <<
" ";
427 diagnostic_file << norm_abs;
428 diagnostic_file <<
" ";
429 diagnostic_file << norm_rel;
430 diagnostic_file <<
" ";
431 diagnostic_file << linear_solver_iters;
432 diagnostic_file <<
" ";
434 diagnostic_file <<
" ";
436 diagnostic_file <<
"\n";
437 diagnostic_file.close();
442template <
class Vec,
class Ops>
450 m_ops->ComputeRHS(
m_R, a_U, a_time, a_iter,
false );
LinearSolverType
struct to select the linear solver for implicit schemes
Definition LinearSolverLibrary.H:14
@ petsc_ksp
Definition LinearSolverLibrary.H:14
@ amrex_gmres
Definition LinearSolverLibrary.H:14
PreconditionerType
Types for preconditioners for field solvers.
Definition Preconditioner.H:20
@ none
Definition Preconditioner.H:20
#define WARPX_ABORT_WITH_MESSAGE(MSG)
Definition TextMsg.H:15
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:13
Vec m_dU
Intermediate Vec containers used by the solver.
Definition NewtonSolver.H:99
int m_maxits
Maximum iterations for the Newton solver.
Definition NewtonSolver.H:124
int m_linsol_verbose_int
Verbosity flag for linear solver.
Definition NewtonSolver.H:154
int m_total_linsol_iters
Total linear iterations for the diagnostic file.
Definition NewtonSolver.H:149
~NewtonSolver() override=default
amrex::Real m_linsol_atol
Absolute tolerance for linear solver.
Definition NewtonSolver.H:139
PreconditionerType m_pc_type
Preconditioner type.
Definition NewtonSolver.H:164
int m_linsol_restart_length
Restart iteration for linear solver.
Definition NewtonSolver.H:159
void Solve(Vec &a_U, const Vec &a_b, amrex::Real a_time, amrex::Real a_dt, int a_step) const override
Solve the specified nonlinear equation for U. Picard: U = b + R(U). Newton: F(U) = U - b - R(U) = 0.
Definition NewtonSolver.H:305
void Define(const Vec &a_U, Ops *a_ops) override
Read user-provided parameters that control the nonlinear solver. Allocate intermediate data container...
Definition NewtonSolver.H:198
void GetSolverParams(amrex::Real &a_rtol, amrex::Real &a_atol, int &a_maxits) override
Return the convergence parameters used by the nonlinear solver.
Definition NewtonSolver.H:54
void ParseParameters()
Definition NewtonSolver.H:273
void PrintParams() const override
Print parameters used by the nonlinear solver.
Definition NewtonSolver.H:75
std::unique_ptr< JacobianFunctionMF< Vec, Ops > > m_linear_function
The linear function used by linear solver to compute A*v. In the contect of JFNK, A = dF/dU (i....
Definition NewtonSolver.H:172
amrex::Real m_atol
Absolute tolerance for the Newton solver.
Definition NewtonSolver.H:119
bool m_require_convergence
Flag to determine whether convergence is required.
Definition NewtonSolver.H:109
NewtonSolver & operator=(const NewtonSolver &)=delete
LinearSolverType m_linear_solver_type
Choice of linear solver.
Definition NewtonSolver.H:177
NewtonSolver(const NewtonSolver &)=delete
amrex::Real m_cur_time
Definition NewtonSolver.H:166
NewtonSolver(NewtonSolver &&) noexcept=delete
std::unique_ptr< LinearSolver< Vec, JacobianFunctionMF< Vec, Ops > > > m_linear_solver
The linear solver object.
Definition NewtonSolver.H:182
Vec m_F
Definition NewtonSolver.H:99
amrex::Real m_linsol_rtol
Relative tolerance for linear solver.
Definition NewtonSolver.H:134
Vec m_R
Definition NewtonSolver.H:99
amrex::Real m_dt
Definition NewtonSolver.H:166
void CurTime(amrex::Real a_time) const
Definition NewtonSolver.H:63
int m_total_iters
Total nonlinear iterations for the diagnostic file.
Definition NewtonSolver.H:129
void CurTimeStep(amrex::Real a_dt) const
Definition NewtonSolver.H:69
int m_linsol_maxits
Maximum iterations for linear solver.
Definition NewtonSolver.H:144
void EvalResidual(Vec &a_F, const Vec &a_U, const Vec &a_b, amrex::Real a_time, int a_iter) const
Compute the nonlinear residual: F(U) = U - b - R(U).
Definition NewtonSolver.H:443
amrex::Real m_rtol
Relative tolerance for the Newton solver.
Definition NewtonSolver.H:114
Ops * m_ops
Pointer to Ops class.
Definition NewtonSolver.H:104
bool m_usePC
Definition NonlinearSolver.H:92
bool m_is_defined
Definition NonlinearSolver.H:88
std::string m_diagnostic_file
Definition NonlinearSolver.H:90
NonlinearSolver()=default
int m_diagnostic_interval
Definition NonlinearSolver.H:91
bool m_verbose
Definition NonlinearSolver.H:89
int query(std::string_view name, bool &ref, int ival=FIRST) const
void WMRecordWarning(const std::string &topic, const std::string &text, const WarnPriority &priority=WarnPriority::medium)
Helper function to abbreviate the call to WarnManager::GetInstance().RecordWarning (recording a warni...
Definition WarnManager.cpp:320
bool IOProcessor() noexcept
std::string getEnumNameString(T const &v)
bool FileExists(const std::string &filename)