7 #ifndef NEWTON_SOLVER_H_
8 #define NEWTON_SOLVER_H_
26 template<
class Vec,
class Ops>
41 void Define ( const Vec& a_U,
42 Ops* a_ops ) override;
44 void Solve ( Vec& a_U,
47 amrex::Real a_dt ) const override;
51 int& a_maxits )
override
58 inline void CurTime ( amrex::Real a_time )
const
170 template <
class Vec,
class Ops>
176 "Newton nonlinear solver object is already defined!");
186 m_linear_function = std::make_unique<JacobianFunctionMF<Vec,Ops>>();
187 m_linear_function->define(m_F, m_ops);
189 m_linear_solver = std::make_unique<amrex::GMRES<Vec,JacobianFunctionMF<Vec,Ops>>>();
190 m_linear_solver->define(*m_linear_function);
191 m_linear_solver->setVerbose( m_gmres_verbose_int );
192 m_linear_solver->setRestartLength( m_gmres_restart_length );
193 m_linear_solver->setMaxIters( m_gmres_maxits );
195 this->m_is_defined =
true;
199 template <
class Vec,
class Ops>
203 pp_newton.
query(
"verbose", this->m_verbose);
204 pp_newton.
query(
"absolute_tolerance", m_atol);
205 pp_newton.
query(
"relative_tolerance", m_rtol);
206 pp_newton.
query(
"max_iterations", m_maxits);
207 pp_newton.
query(
"require_convergence", m_require_convergence);
210 pp_gmres.
query(
"verbose_int", m_gmres_verbose_int);
211 pp_gmres.
query(
"restart_length", m_gmres_restart_length);
212 pp_gmres.
query(
"absolute_tolerance", m_gmres_atol);
213 pp_gmres.
query(
"relative_tolerance", m_gmres_rtol);
214 pp_gmres.
query(
"max_iterations", m_gmres_maxits);
217 template <
class Vec,
class Ops>
221 amrex::Real a_dt )
const
226 "NewtonSolver::Solve() called on undefined object");
227 using namespace amrex::literals;
237 amrex::Real norm_abs = 0.;
238 amrex::Real norm0 = 1._rt;
239 amrex::Real norm_rel = 0.;
242 for (iter = 0; iter < m_maxits;) {
245 EvalResidual(m_F, a_U, a_b, a_time, a_dt, iter);
248 norm_abs = m_F.norm2();
250 if (norm_abs > 0.) { norm0 = norm_abs; }
251 else { norm0 = 1._rt; }
253 norm_rel = norm_abs/norm0;
256 if (this->m_verbose || iter == m_maxits) {
257 amrex::Print() <<
"Newton: iteration = " << std::setw(3) << iter <<
", norm = "
258 << std::scientific << std::setprecision(5) << norm_abs <<
" (abs.), "
259 << std::scientific << std::setprecision(5) << norm_rel <<
" (rel.)" <<
"\n";
262 if (norm_abs < m_rtol) {
263 amrex::Print() <<
"Newton: exiting at iteration = " << std::setw(3) << iter
264 <<
". Satisfied absolute tolerance " << m_atol << std::endl;
268 if (norm_rel < m_rtol) {
269 amrex::Print() <<
"Newton: exiting at iteration = " << std::setw(3) << iter
270 <<
". Satisfied relative tolerance " << m_rtol << std::endl;
274 if (norm_abs > 100._rt*norm0) {
275 amrex::Print() <<
"Newton: exiting at iteration = " << std::setw(3) << iter
276 <<
". SOLVER DIVERGED! relative tolerance = " << m_rtol << std::endl;
277 std::stringstream convergenceMsg;
278 convergenceMsg <<
"Newton: exiting at iteration " << std::setw(3) << iter <<
279 ". SOLVER DIVERGED! absolute norm = " << norm_abs <<
280 " has increased by 100X from that after first iteration.";
286 m_linear_solver->solve( m_dU, m_F, m_gmres_rtol, m_gmres_atol );
292 if (iter >= m_maxits) {
293 amrex::Print() <<
"Newton: exiting at iter = " << std::setw(3) << iter
294 <<
". Maximum iteration reached: iter = " << m_maxits << std::endl;
300 if (m_rtol > 0. && iter == m_maxits) {
301 std::stringstream convergenceMsg;
302 convergenceMsg <<
"Newton solver failed to converge after " << iter <<
303 " iterations. Relative norm is " << norm_rel <<
304 " and the relative tolerance is " << m_rtol <<
305 ". Absolute norm is " << norm_abs <<
306 " and the absolute tolerance is " << m_atol;
307 if (this->m_verbose) {
amrex::Print() << convergenceMsg.str() << std::endl; }
308 if (m_require_convergence) {
317 template <
class Vec,
class Ops>
326 m_ops->ComputeRHS( m_R, a_U, a_time, a_dt, a_iter,
false );
329 m_linear_function->setBaseSolution(a_U);
330 m_linear_function->setBaseRHS(m_R);
333 if (m_update_pc || m_update_pc_init) {
334 m_linear_function->updatePreCondMat(a_U);
336 m_update_pc_init =
false;
#define WARPX_ABORT_WITH_MESSAGE(MSG)
Definition: TextMsg.H:15
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition: TextMsg.H:13
Newton method to solve nonlinear equation of form: F(U) = U - b - R(U) = 0. U is the solution vector,...
Definition: NewtonSolver.H:28
Vec m_dU
Intermediate Vec containers used by the solver.
Definition: NewtonSolver.H:89
int m_maxits
Maximum iterations for the Newton solver.
Definition: NewtonSolver.H:114
int m_gmres_maxits
Maximum iterations for GMRES.
Definition: NewtonSolver.H:129
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:171
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:49
void ParseParameters()
Definition: NewtonSolver.H:200
void PrintParams() const override
Print parameters used by the nonlinear solver.
Definition: NewtonSolver.H:70
std::unique_ptr< JacobianFunctionMF< Vec, Ops > > m_linear_function
The linear function used by GMRES to compute A*v. In the contect of JFNK, A = dF/dU (i....
Definition: NewtonSolver.H:149
amrex::Real m_gmres_atol
Absolute tolerance for GMRES.
Definition: NewtonSolver.H:124
amrex::Real m_atol
Absolute tolerance for the Newton solver.
Definition: NewtonSolver.H:109
amrex::Real m_gmres_rtol
Relative tolerance for GMRES.
Definition: NewtonSolver.H:119
bool m_require_convergence
Flag to determine whether convergence is required.
Definition: NewtonSolver.H:99
NewtonSolver(const NewtonSolver &)=delete
amrex::Real m_cur_time
Definition: NewtonSolver.H:141
std::unique_ptr< amrex::GMRES< Vec, JacobianFunctionMF< Vec, Ops > > > m_linear_solver
The linear solver (GMRES) object.
Definition: NewtonSolver.H:154
void EvalResidual(Vec &a_F, const Vec &a_U, const Vec &a_b, amrex::Real a_time, amrex::Real a_dt, int a_iter) const
Compute the nonlinear residual: F(U) = U - b - R(U).
Definition: NewtonSolver.H:318
int m_gmres_verbose_int
Verbosity flag for GMRES.
Definition: NewtonSolver.H:134
NewtonSolver(NewtonSolver &&) noexcept=delete
bool m_update_pc
Definition: NewtonSolver.H:142
Vec m_F
Definition: NewtonSolver.H:89
bool m_update_pc_init
Definition: NewtonSolver.H:143
NewtonSolver & operator=(const NewtonSolver &)=delete
Vec m_R
Definition: NewtonSolver.H:89
amrex::Real m_dt
Definition: NewtonSolver.H:141
void CurTime(amrex::Real a_time) const
Definition: NewtonSolver.H:58
void CurTimeStep(amrex::Real a_dt) const
Definition: NewtonSolver.H:64
int m_gmres_restart_length
Restart iteration for GMRES.
Definition: NewtonSolver.H:139
void Solve(Vec &a_U, const Vec &a_b, amrex::Real a_time, amrex::Real a_dt) 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:218
amrex::Real m_rtol
Relative tolerance for the Newton solver.
Definition: NewtonSolver.H:104
Ops * m_ops
Pointer to Ops class.
Definition: NewtonSolver.H:94
Top-level class for the nonlinear solver.
Definition: NonlinearSolver.H:28
bool m_verbose
Definition: NonlinearSolver.H:83
int query(const char *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:318