WarpX
Loading...
Searching...
No Matches
PicardSolver.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 PICARD_SOLVER_H_
8#define PICARD_SOLVER_H_
9
10#include "NonlinearSolver.H"
11
12#include <AMReX_ParmParse.H>
13#include "Utils/TextMsg.H"
14
15#include <vector>
16#include <istream>
17#include <filesystem>
18
25
26template<class Vec, class Ops>
27class PicardSolver : public NonlinearSolver<Vec,Ops>
28{
29public:
30
31 PicardSolver() = default;
32
33 ~PicardSolver() override = default;
34
35 // Prohibit Move and Copy operations
36 PicardSolver(const PicardSolver&) = delete;
38 PicardSolver(PicardSolver&&) noexcept = delete;
39 PicardSolver& operator=(PicardSolver&&) noexcept = delete;
40
41 void Define ( const Vec& a_U,
42 Ops* a_ops ) override;
43
44 void Solve ( Vec& a_U,
45 const Vec& a_b,
46 amrex::Real a_time,
47 amrex::Real a_dt,
48 int a_step) const override;
49
50 void GetSolverParams ( amrex::Real& a_rtol,
51 amrex::Real& a_atol,
52 int& a_maxits ) override
53 {
54 a_rtol = m_rtol;
55 a_atol = m_atol;
56 a_maxits = m_maxits;
57 }
58
59 void PrintParams () const override
60 {
61 amrex::Print() << "Picard max iterations: " << m_maxits << "\n";
62 amrex::Print() << "Picard relative tolerance: " << m_rtol << "\n";
63 amrex::Print() << "Picard absolute tolerance: " << m_atol << "\n";
64 amrex::Print() << "Picard require convergence: " << (m_require_convergence?"true":"false") << "\n";
65 }
66
67private:
68
72 mutable Vec m_Usave, m_R;
73
77 Ops* m_ops = nullptr;
78
83
87 amrex::Real m_rtol = 1.0e-6;
88
92 amrex::Real m_atol = 0.;
93
97 int m_maxits = 100;
98
102 mutable int m_total_iters = 0;
103
104 void ParseParameters( );
105
106};
107
108template <class Vec, class Ops>
109void PicardSolver<Vec,Ops>::Define ( const Vec& a_U,
110 Ops* a_ops )
111{
113 !this->m_is_defined,
114 "Picard nonlinear solver object is already defined!");
115
117
118 m_Usave.Define(a_U);
119 m_R.Define(a_U);
120
121 m_ops = a_ops;
122
123 this->m_is_defined = true;
124
125 // Create diagnostic file and write header
127 && !this->m_diagnostic_file.empty()
128 && !amrex::FileExists(this->m_diagnostic_file)) {
129
130 std::filesystem::path const diagnostic_path(this->m_diagnostic_file);
131 std::filesystem::path const diagnostic_dir = diagnostic_path.parent_path();
132 if (!diagnostic_dir.empty()) {
133 std::filesystem::create_directories(diagnostic_dir);
134 }
135
136 std::ofstream diagnostic_file{this->m_diagnostic_file, std::ofstream::out | std::ofstream::trunc};
137 int c = 0;
138 diagnostic_file << "#";
139 diagnostic_file << "[" << c++ << "]step()";
140 diagnostic_file << " ";
141 diagnostic_file << "[" << c++ << "]time(s)";
142 diagnostic_file << " ";
143 diagnostic_file << "[" << c++ << "]iters";
144 diagnostic_file << " ";
145 diagnostic_file << "[" << c++ << "]total_iters";
146 diagnostic_file << " ";
147 diagnostic_file << "[" << c++ << "]norm_abs";
148 diagnostic_file << " ";
149 diagnostic_file << "[" << c++ << "]norm_rel";
150 diagnostic_file << "\n";
151 diagnostic_file.close();
152 }
153
154}
155
156template <class Vec, class Ops>
158{
159 const amrex::ParmParse pp_picard("picard");
160 pp_picard.query("verbose", this->m_verbose);
161 pp_picard.query("absolute_tolerance", m_atol);
162 pp_picard.query("relative_tolerance", m_rtol);
163 pp_picard.query("max_iterations", m_maxits);
164 pp_picard.query("require_convergence", m_require_convergence);
165 pp_picard.query("diagnostic_file", this->m_diagnostic_file);
166 pp_picard.query("diagnostic_interval", this->m_diagnostic_interval);
167}
168
169template <class Vec, class Ops>
171 const Vec& a_b,
172 amrex::Real a_time,
173 amrex::Real a_dt,
174 int a_step) const
175{
176 BL_PROFILE("PicardSolver::Solve()");
178 this->m_is_defined,
179 "PicardSolver::Solve() called on undefined object");
180 using namespace amrex::literals;
181
182 //
183 // Picard fixed-point iteration method to solve nonlinear
184 // equation of form: U = b + R(U)
185 //
186
187 amrex::Real norm_abs = 0.;
188 amrex::Real norm0 = 1._rt;
189 amrex::Real norm_rel = 0.;
190
191 int iter;
192 for (iter = 0; iter < m_maxits;) {
193
194 // Save previous state for norm calculation
195 m_Usave.Copy(a_U);
196
197 // Update the solver state (a_U = a_b + m_R)
198 m_ops->ComputeRHS( m_R, a_U, a_time, iter, false );
199 a_U.Copy(a_b);
200 a_U += m_R;
201
202 // Compute the step norm and update iter
203 m_Usave -= a_U;
204 norm_abs = m_Usave.norm2();
205 if (iter == 0) {
206 if (norm_abs > 0.) { norm0 = norm_abs; }
207 else { norm0 = 1._rt; }
208 }
209 norm_rel = norm_abs/norm0;
210 iter++;
211
212 // Check for convergence criteria
213 if (this->m_verbose || iter == m_maxits) {
214 amrex::Print() << "Picard: iter = " << std::setw(3) << iter << ", norm = "
215 << std::scientific << std::setprecision(5) << norm_abs << " (abs.), "
216 << std::scientific << std::setprecision(5) << norm_rel << " (rel.)" << "\n";
217 }
218
219 if (norm_abs < m_atol) {
220 if (this->m_verbose) {
221 amrex::Print() << "Picard: exiting at iter = " << std::setw(3) << iter
222 << ". Satisfied absolute tolerance " << m_atol << "\n";
223 }
224 break;
225 }
226
227 if (norm_rel < m_rtol) {
228 if (this->m_verbose) {
229 amrex::Print() << "Picard: exiting at iter = " << std::setw(3) << iter
230 << ". Satisfied relative tolerance " << m_rtol << "\n";
231 }
232 break;
233 }
234
235 if (iter >= m_maxits) {
236 if (this->m_verbose) {
237 amrex::Print() << "Picard: exiting at iter = " << std::setw(3) << iter
238 << ". Maximum iteration reached: iter = " << m_maxits << "\n";
239 }
240 break;
241 }
242
243 }
244
245 // Update total iteration count
246 m_total_iters += iter;
247
248 if (m_rtol > 0. && iter == m_maxits) {
249 std::stringstream convergenceMsg;
250 convergenceMsg << "Picard solver failed to converge after " << iter <<
251 " iterations. Relative norm is " << norm_rel <<
252 " and the relative tolerance is " << m_rtol <<
253 ". Absolute norm is " << norm_abs <<
254 " and the absolute tolerance is " << m_atol;
255 if (this->m_verbose) { amrex::Print() << convergenceMsg.str() << "\n"; }
257 WARPX_ABORT_WITH_MESSAGE(convergenceMsg.str());
258 } else {
259 ablastr::warn_manager::WMRecordWarning("PicardSolver", convergenceMsg.str());
260 }
261 }
262
264 ((a_step+1)%this->m_diagnostic_interval==0 || a_step==0)) {
265 std::ofstream diagnostic_file{this->m_diagnostic_file, std::ofstream::out | std::ofstream::app};
266 diagnostic_file << std::setprecision(14);
267 diagnostic_file << a_step+1;
268 diagnostic_file << " ";
269 diagnostic_file << a_time + a_dt;
270 diagnostic_file << " ";
271 diagnostic_file << iter;
272 diagnostic_file << " ";
273 diagnostic_file << m_total_iters;
274 diagnostic_file << " ";
275 diagnostic_file << norm_abs;
276 diagnostic_file << " ";
277 diagnostic_file << norm_rel;
278 diagnostic_file << "\n";
279 diagnostic_file.close();
280 }
281
282}
283
284#endif
#define BL_PROFILE(a)
#define WARPX_ABORT_WITH_MESSAGE(MSG)
Definition TextMsg.H:15
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:13
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
PicardSolver & operator=(const PicardSolver &)=delete
void ParseParameters()
Definition PicardSolver.H:157
PicardSolver()=default
PicardSolver(const PicardSolver &)=delete
bool m_require_convergence
Flag to determine whether convergence is required.
Definition PicardSolver.H:82
Vec m_R
Definition PicardSolver.H:72
void GetSolverParams(amrex::Real &a_rtol, amrex::Real &a_atol, int &a_maxits) override
Return the convergence parameters used by the nonlinear solver.
Definition PicardSolver.H:50
~PicardSolver() override=default
amrex::Real m_rtol
Relative tolerance for the Picard nonlinear solver.
Definition PicardSolver.H:87
Ops * m_ops
Pointer to Ops class.
Definition PicardSolver.H:77
int m_maxits
Maximum iterations for the Picard nonlinear solver.
Definition PicardSolver.H:97
amrex::Real m_atol
Absolute tolerance for the Picard nonlinear solver.
Definition PicardSolver.H:92
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 PicardSolver.H:170
Vec m_Usave
Intermediate Vec containers used by the solver.
Definition PicardSolver.H:72
PicardSolver(PicardSolver &&) noexcept=delete
void Define(const Vec &a_U, Ops *a_ops) override
Read user-provided parameters that control the nonlinear solver. Allocate intermediate data container...
Definition PicardSolver.H:109
int m_total_iters
Total iterations for the diagnostic file.
Definition PicardSolver.H:102
void PrintParams() const override
Print parameters used by the nonlinear solver.
Definition PicardSolver.H:59
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
bool FileExists(const std::string &filename)