WarpX
Loading...
Searching...
No Matches
NewtonSolver.H
Go to the documentation of this file.
1/* Copyright 2024 Justin Angus
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef NEWTON_SOLVER_H_
8#define NEWTON_SOLVER_H_
9
10#include "LinearSolverLibrary.H"
11#include "NonlinearSolver.H"
12#include "JacobianFunctionMF.H"
13#include "Preconditioner.H"
14#include "Utils/TextMsg.H"
15
16#include <AMReX_ParmParse.H>
17#include <AMReX_Config.H>
18
19#include <vector>
20#include <istream>
21#include <filesystem>
22
29
30template<class Vec, class Ops>
31class NewtonSolver : public NonlinearSolver<Vec,Ops>
32{
33public:
34
35 NewtonSolver() = default;
36
37 ~NewtonSolver() override = default;
38
39 // Prohibit Move and Copy operations
40 NewtonSolver(const NewtonSolver&) = delete;
42 NewtonSolver(NewtonSolver&&) noexcept = delete;
43 NewtonSolver& operator=(NewtonSolver&&) noexcept = delete;
44
45 void Define ( const Vec& a_U,
46 Ops* a_ops ) override;
47
48 void Solve ( Vec& a_U,
49 const Vec& a_b,
50 amrex::Real a_time,
51 amrex::Real a_dt,
52 int a_step) const override;
53
54 void GetSolverParams ( amrex::Real& a_rtol,
55 amrex::Real& a_atol,
56 int& a_maxits ) override
57 {
58 a_rtol = m_rtol;
59 a_atol = m_atol;
60 a_maxits = m_maxits;
61 }
62
63 inline void CurTime ( amrex::Real a_time ) const
64 {
65 m_cur_time = a_time;
66 m_linear_function->curTime( a_time );
67 }
68
69 inline void CurTimeStep ( amrex::Real a_dt ) const
70 {
71 m_dt = a_dt;
72 m_linear_function->curTimeStep( a_dt );
73 }
74
75 void PrintParams () const override
76 {
77 amrex::Print() << "Newton verbose: " << (this->m_verbose?"true":"false") << "\n";
78 amrex::Print() << "Newton max iterations: " << m_maxits << "\n";
79 amrex::Print() << "Newton relative tolerance: " << m_rtol << "\n";
80 amrex::Print() << "Newton absolute tolerance: " << m_atol << "\n";
81 amrex::Print() << "Newton require convergence: " << (m_require_convergence?"true":"false") << "\n";
83 amrex::Print() << "Newton linear solver: " << linsol_name << "\n";
84 amrex::Print() << "Linear solver (" << linsol_name << ") verbose: " << m_linsol_verbose_int << "\n";
85 amrex::Print() << "Linear solver (" << linsol_name << ") restart length: " << m_linsol_restart_length << "\n";
86 amrex::Print() << "Linear solver (" << linsol_name << ") max iterations: " << m_linsol_maxits << "\n";
87 amrex::Print() << "Linear solver (" << linsol_name << ") relative tolerance: " << m_linsol_rtol << "\n";
88 amrex::Print() << "Linear solver (" << linsol_name << ") absolute tolerance: " << m_linsol_atol << "\n";
89 amrex::Print() << "Preconditioner type: " << amrex::getEnumNameString(m_pc_type) << "\n";
90
91 m_linear_function->printParams();
92 }
93
94private:
95
99 mutable Vec m_dU, m_F, m_R;
100
104 Ops* m_ops = nullptr;
105
110
114 amrex::Real m_rtol = 1.0e-6;
115
119 amrex::Real m_atol = 0.;
120
124 int m_maxits = 100;
125
129 mutable int m_total_iters = 0;
130
134 amrex::Real m_linsol_rtol = 1.0e-4;
135
139 amrex::Real m_linsol_atol = 0.;
140
144 int m_linsol_maxits = 1000;
145
149 mutable int m_total_linsol_iters = 0;
150
155
160
165
166 mutable amrex::Real m_cur_time, m_dt;
167
172 std::unique_ptr<JacobianFunctionMF<Vec,Ops>> m_linear_function;
173
178
182 std::unique_ptr<LinearSolver<Vec,JacobianFunctionMF<Vec,Ops>>> m_linear_solver;
183
184 void ParseParameters ();
185
189 void EvalResidual ( Vec& a_F,
190 const Vec& a_U,
191 const Vec& a_b,
192 amrex::Real a_time,
193 int a_iter ) const;
194
195};
196
197template <class Vec, class Ops>
198void NewtonSolver<Vec,Ops>::Define ( const Vec& a_U,
199 Ops* a_ops )
200{
201 BL_PROFILE("NewtonSolver::Define()");
203 !this->m_is_defined,
204 "Newton nonlinear solver object is already defined!");
205
207
208 m_dU.Define(a_U);
209 m_F.Define(a_U); // residual function F(U) = U - b - R(U) = 0
210 m_R.Define(a_U); // right hand side function R(U)
211
212 m_ops = a_ops;
213
214 m_linear_function = std::make_unique<JacobianFunctionMF<Vec,Ops>>();
216 this->m_usePC = m_linear_function->usePreconditioner();
217
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>>>();
223#else
224 WARPX_ABORT_WITH_MESSAGE("NewtonSolver::Define(): must compile with PETSc to use petsc_ksp (AMREX_USE_PETSC must be defined)");
225#endif
226 } else {
227 WARPX_ABORT_WITH_MESSAGE("NewtonSolver::Define(): invalid linear solver type");
228 }
231 m_linear_solver->setRestartLength( m_linsol_restart_length );
232 m_linear_solver->setMaxIters( m_linsol_maxits );
233
234 this->m_is_defined = true;
235
236 // Create diagnostic file and write header
238 && !this->m_diagnostic_file.empty()
239 && !amrex::FileExists(this->m_diagnostic_file)) {
240
241 std::filesystem::path const diagnostic_path(this->m_diagnostic_file);
242 std::filesystem::path const diagnostic_dir = diagnostic_path.parent_path();
243 if (!diagnostic_dir.empty()) {
244 std::filesystem::create_directories(diagnostic_dir);
245 }
246
247 std::ofstream diagnostic_file{this->m_diagnostic_file, std::ofstream::out | std::ofstream::trunc};
248 int c = 0;
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();
269 }
270}
271
272template <class Vec, class Ops>
274{
275 const amrex::ParmParse pp_newton("newton");
276 pp_newton.query("verbose", this->m_verbose);
277 pp_newton.query("linear_solver", m_linear_solver_type);
278 pp_newton.query("absolute_tolerance", m_atol);
279 pp_newton.query("relative_tolerance", m_rtol);
280 pp_newton.query("max_iterations", m_maxits);
281 pp_newton.query("require_convergence", m_require_convergence);
282 pp_newton.query("diagnostic_file", this->m_diagnostic_file);
283 pp_newton.query("diagnostic_interval", this->m_diagnostic_interval);
284
286 pp_l.query("verbose_int", m_linsol_verbose_int);
287 pp_l.query("restart_length", m_linsol_restart_length);
288 pp_l.query("absolute_tolerance", m_linsol_atol);
289 pp_l.query("relative_tolerance", m_linsol_rtol);
290 pp_l.query("max_iterations", m_linsol_maxits);
291
292 /* backward compatibility */
293 const amrex::ParmParse pp_gmres("gmres");
294 pp_gmres.query("verbose_int", m_linsol_verbose_int);
295 pp_gmres.query("restart_length", m_linsol_restart_length);
296 pp_gmres.query("absolute_tolerance", m_linsol_atol);
297 pp_gmres.query("relative_tolerance", m_linsol_rtol);
298 pp_gmres.query("max_iterations", m_linsol_maxits);
299
300 const amrex::ParmParse pp_jac("jacobian");
301 pp_jac.query("pc_type", m_pc_type);
302}
303
304template <class Vec, class Ops>
306 const Vec& a_b,
307 amrex::Real a_time,
308 amrex::Real a_dt,
309 int a_step) const
310{
311 BL_PROFILE("NewtonSolver::Solve()");
313 this->m_is_defined,
314 "NewtonSolver::Solve() called on undefined object");
315 using namespace amrex::literals;
316
317 //
318 // Newton routine to solve nonlinear equation of form:
319 // F(U) = U - b - R(U) = 0
320 //
321
322 CurTime(a_time);
323 CurTimeStep(a_dt);
324
325 amrex::Real norm_abs = 0.;
326 amrex::Real norm0 = 1._rt;
327 amrex::Real norm_rel = 0.;
328
329 int iter;
330 int linear_solver_iters = 0;
331 for (iter = 0; iter < m_maxits;) {
332
333 // Compute residual: F(U) = U - b - R(U)
334 EvalResidual(m_F, a_U, a_b, a_time, iter);
335
336 // Compute norm of the residual
337 norm_abs = m_F.norm2();
338 if (iter == 0) {
339 if (norm_abs > 0.) { norm0 = norm_abs; }
340 else { norm0 = 1._rt; }
341 }
342 norm_rel = norm_abs/norm0;
343
344 // Check for convergence criteria
345 if (this->m_verbose || iter == m_maxits) {
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";
349 }
350
351 if (norm_abs < m_atol) {
352 if (this->m_verbose) {
353 amrex::Print() << "Newton: exiting at iteration = " << std::setw(3) << iter
354 << ". Satisfied absolute tolerance " << m_atol << "\n";
355 }
356 break;
357 }
358
359 if (norm_rel < m_rtol) {
360 if (this->m_verbose) {
361 amrex::Print() << "Newton: exiting at iteration = " << std::setw(3) << iter
362 << ". Satisfied relative tolerance " << m_rtol << "\n";
363 }
364 break;
365 }
366
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.";
374 WARPX_ABORT_WITH_MESSAGE(convergenceMsg.str());
375 }
376
377 // Solve linear system for Newton step [Jac]*dU = F
378 m_dU.zero();
380 linear_solver_iters += m_linear_solver->getNumIters();
381
382 // Update solution
383 a_U -= m_dU;
384
385 iter++;
386 if (iter >= m_maxits) {
387 if (this->m_verbose) {
388 amrex::Print() << "Newton: exiting at iter = " << std::setw(3) << iter
389 << ". Maximum iteration reached: iter = " << m_maxits << "\n";
390 }
391 break;
392 }
393
394 }
395
396 // Update total iteration count
397 m_total_iters += iter;
398 m_total_linsol_iters += linear_solver_iters;
399
400 if (m_rtol > 0. && iter == m_maxits) {
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;
407 if (this->m_verbose) { amrex::Print() << convergenceMsg.str() << "\n"; }
409 WARPX_ABORT_WITH_MESSAGE(convergenceMsg.str());
410 } else {
411 ablastr::warn_manager::WMRecordWarning("NewtonSolver", convergenceMsg.str());
412 }
413 }
414
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 << " ";
425 diagnostic_file << m_total_iters;
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 << " ";
433 diagnostic_file << m_total_linsol_iters;
434 diagnostic_file << " ";
435 diagnostic_file << m_linear_solver->getResidualNorm();
436 diagnostic_file << "\n";
437 diagnostic_file.close();
438 }
439
440}
441
442template <class Vec, class Ops>
444 const Vec& a_U,
445 const Vec& a_b,
446 amrex::Real a_time,
447 int a_iter ) const
448{
449 BL_PROFILE("NewtonSolver::EvalResidual()");
450 m_ops->ComputeRHS( m_R, a_U, a_time, a_iter, false );
451
452 // set base U and R(U) for matrix-free Jacobian action calculation
453 m_linear_function->setBaseSolution(a_U);
454 m_linear_function->setBaseRHS(m_R);
455
456 // update preconditioner
457 m_linear_function->updatePreCondMat(a_U);
458
459 // Compute residual: F(U) = U - b - R(U)
460 a_F.Copy(a_U);
461 a_F -= m_R;
462 a_F -= a_b;
463
464}
465
466#endif
#define BL_PROFILE(a)
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
NewtonSolver()=default
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)