Cantera  3.3.0a1
Loading...
Searching...
No Matches
SteadyStateSystem.cpp
Go to the documentation of this file.
1//! @file SteadyStateSystem.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
9
10using namespace std;
11
12namespace Cantera
13{
14
15SteadyStateSystem::SteadyStateSystem()
16{
17 m_state = make_shared<vector<double>>();
18 m_newt = make_unique<MultiNewton>(1);
19}
20
21SteadyStateSystem::~SteadyStateSystem()
22{
23}
24
26{
29 m_state->assign(x, x + size());
30}
31
32void SteadyStateSystem::getState(double* x) const
33{
34 copy(m_xnew.begin(), m_xnew.end(), x);
35}
36
37void SteadyStateSystem::solve(int loglevel)
38{
39 size_t istep = 0;
40 int nsteps = m_steps[istep];
41 m_nsteps = 0;
42 double dt = m_tstep;
43
44 while (true) {
45 // Keep the attempt_counter in the range of [1, max_history]
47
48 // Attempt to solve the steady problem
51 debuglog("\nAttempt Newton solution of steady-state problem.", loglevel);
52 if (!m_jac_ok) {
53 evalJacobian(m_state->data());
54 try {
55 m_jac->updateTransient(m_rdt, m_mask.data());
56 m_jac_ok = true;
57 } catch (CanteraError& err) {
58 // Allow solver to continue after failure to factorize the steady-state
59 // Jacobian by returning to time stepping mode
60 if (m_rdt == 0.0) {
61 if (loglevel > 1) {
62 writelog("\nSteady Jacobian factorization failed:"
63 "\n {}:\n {}",
64 err.getMethod(), err.getMessage());
65 }
66 } else {
67 throw;
68 }
69 }
70 }
71 int m = -1;
72 if (m_jac_ok) {
73 m = newton().solve(m_state->data(), m_xnew.data(), *this, loglevel);
74 }
75 if (m == 1) {
76 *m_state = m_xnew;
77 writeDebugInfo("NewtonSuccess", "After successful Newton solve",
78 loglevel, m_attempt_counter);
79
80 return;
81 } else {
82 debuglog("\nNewton steady-state solve failed.\n", loglevel);
83 writeDebugInfo("NewtonFail", "After unsuccessful Newton solve",
84 loglevel, m_attempt_counter);
85
86 if (loglevel > 0) {
87 writelog("\nAttempt {} timesteps.", nsteps);
88 }
89
90 dt = timeStep(nsteps, dt, m_state->data(), m_xnew.data(), loglevel-1);
92 writeDebugInfo("Timestepping", "After timestepping", loglevel,
94
95 // Repeat the last timestep's data for logging purposes
96 if (loglevel == 1) {
97 writelog("\nFinal timestep info: dt= {:<10.4g} log(ss)= {:<10.4g}\n", dt,
98 log10(ssnorm(m_state->data(), m_xnew.data())));
99 }
100 istep++;
101 if (istep >= m_steps.size()) {
102 nsteps = m_steps.back();
103 } else {
104 nsteps = m_steps[istep];
105 }
106 dt = std::min(dt, m_tmax);
107 }
108 }
109}
110
111double SteadyStateSystem::timeStep(int nsteps, double dt, double* x, double* r, int loglevel)
112{
113 // set the Jacobian age parameter to the transient value
115
116 int n = 0;
117 int successiveFailures = 0;
118
119 // Only output this if nothing else under this function call will be output
120 if (loglevel == 1) {
121 writelog("\n============================");
122 writelog("\n{:<5s} {:<11s} {:<7s}\n", "step", "dt (s)", "log(ss)");
123 writelog("============================");
124 }
125 while (n < nsteps) {
126 if (loglevel == 1) { // At level 1, output concise information
127 double ss = ssnorm(x, r);
128 writelog("\n{:<5d} {:<6.4e} {:>7.4f}", n, dt, log10(ss));
129 } else if (loglevel > 1) {
130 double ss = ssnorm(x, r);
131 writelog("\nTimestep ({}) dt= {:<11.4e} log(ss)= {:<7.4f}", n, dt, log10(ss));
132 }
133
134 // set up for time stepping with stepsize dt
135 initTimeInteg(dt,x);
136
137 int j0 = m_jac->nEvals(); // Store the current number of Jacobian evaluations
138
139 // solve the transient problem
140 int status = newton().solve(x, r, *this, loglevel);
141
142 // successful time step. Copy the new solution in r to
143 // the current solution in x.
144 if (status >= 0) {
145 if (loglevel > 1) {
146 writelog("\nTimestep ({}) succeeded", n);
147 }
148 successiveFailures = 0;
149 m_nsteps++;
150 n += 1;
151 copy(r, r + m_size, x);
152 // No Jacobian evaluations were performed, so a larger timestep can be used
153 if (m_jac->nEvals() == j0) {
154 dt *= 1.5;
155 }
158 }
159 dt = std::min(dt, m_tmax);
160 if (m_nsteps >= m_nsteps_max) {
161 throw CanteraError("OneDim::timeStep",
162 "Took maximum number of timesteps allowed ({}) without "
163 "reaching steady-state solution.", m_nsteps_max);
164 }
165 } else {
166 // No solution could be found with this time step.
167 // Decrease the stepsize and try again.
168 successiveFailures++;
169 if (loglevel == 1) {
170 writelog("\nTimestep failed");
171 } else if (loglevel > 1) {
172 writelog("\nTimestep ({}) failed", n);
173 }
174 if (successiveFailures > 2) {
175 debuglog("--> Resetting negative species concentrations", loglevel);
177 successiveFailures = 0;
178 } else {
179 debuglog("--> Reducing timestep", loglevel);
180 dt *= m_tfactor;
181 if (dt < m_tmin) {
182 string err_msg = fmt::format(
183 "Time integration failed. Minimum timestep ({}) reached.", m_tmin);
184 throw CanteraError("OneDim::timeStep", err_msg);
185 }
186 }
187 }
188 }
189
190 // Write the final step to the log
191 if (loglevel == 1) {
192 double ss = ssnorm(x, r);
193 writelog("\n{:<5d} {:<6.4e} {:>7.4f}", n, dt, log10(ss));
194 writelog("\n============================");
195 } else if (loglevel > 1) {
196 double ss = ssnorm(x, r);
197 writelog("\nTimestep ({}) dt= {:<11.4e} log10(ss)= {:<7.4f}\n", n, dt, log10(ss));
198 }
199
200 // return the value of the last stepsize, which may be smaller
201 // than the initial stepsize
202 return dt;
203}
204
205double SteadyStateSystem::ssnorm(double* x, double* r)
206{
207 eval(x, r, 0.0, 0);
208 double ss = 0.0;
209 for (size_t i = 0; i < m_size; i++) {
210 ss = std::max(fabs(r[i]),ss);
211 }
212 return ss;
213}
214
215void SteadyStateSystem::setTimeStep(double stepsize, size_t n, const int* tsteps)
216{
217 m_tstep = stepsize;
218 m_steps.resize(n);
219 for (size_t i = 0; i < n; i++) {
220 m_steps[i] = tsteps[i];
221 }
222}
223
225{
226 m_state->resize(size());
227 m_xnew.resize(size());
228 m_newt->resize(size());
229 m_mask.resize(size());
230 if (!m_jac) {
231 throw CanteraError("SteadyStateSystem::resize",
232 "Jacobian evaluator must be instantiated before calling resize()");
233 }
234 m_jac->initialize(size());
235 m_jac->setBandwidth(bandwidth());
236 m_jac->clearStats();
237 m_jac_ok = false;
238}
239
240void SteadyStateSystem::setLinearSolver(shared_ptr<SystemJacobian> solver)
241{
242 m_jac = solver;
243 m_jac->initialize(size());
244 m_jac->setBandwidth(bandwidth());
245 m_jac->clearStats();
246 m_jac_ok = false;
247}
248
249void SteadyStateSystem::evalSSJacobian(double* x, double* rsd)
250{
251 double rdt_save = m_rdt;
252 m_jac_ok = false;
254 evalJacobian(x);
255 m_rdt = rdt_save;
256}
257
258void SteadyStateSystem::initTimeInteg(double dt, double* x)
259{
260 double rdt_old = m_rdt;
261 m_rdt = 1.0/dt;
262
263 // if the stepsize has changed, then update the transient part of the Jacobian
264 if (fabs(rdt_old - m_rdt) > Tiny) {
265 m_jac->updateTransient(m_rdt, m_mask.data());
266 }
267}
268
270{
271 if (m_rdt == 0) {
272 return;
273 }
274
275 m_rdt = 0.0;
276 m_jac->updateTransient(m_rdt, m_mask.data());
277}
278
279void SteadyStateSystem::setJacAge(int ss_age, int ts_age)
280{
281 m_ss_jac_age = ss_age;
282 if (ts_age > 0) {
283 m_ts_jac_age = ts_age;
284 } else {
286 }
287}
288
290{
291 return *m_newt;
292}
293
294}
Base class for exceptions thrown by Cantera classes.
virtual string getMessage() const
Method overridden by derived classes to format the error message.
virtual string getMethod() const
Get the name of the method that threw the exception.
virtual double eval(double t) const
Evaluate the function.
Definition Func1.cpp:28
Newton iterator for multi-domain, one-dimensional problems.
Definition MultiNewton.h:24
int solve(double *x0, double *x1, SteadyStateSystem &r, int loglevel)
Find the solution to F(x) = 0 by damped Newton iteration.
void setOptions(int maxJacAge=5)
Set options.
int m_nsteps
Number of time steps taken in the current call to solve()
size_t m_size
Solution vector size
int m_nsteps_max
Maximum number of timesteps allowed per call to solve()
virtual void resize()
Call to set the size of internal data structures after first defining the system or if the problem si...
vector< double > m_xnew
Work array used to hold the residual or the new solution.
unique_ptr< MultiNewton > m_newt
Newton iterator.
virtual void resetBadValues(double *x)
Reset values such as negative species concentrations.
size_t size() const
Total solution vector length;.
double ssnorm(double *x, double *r)
Steady-state max norm (infinity norm) of the residual evaluated using solution x.
virtual void initTimeInteg(double dt, double *x)
Prepare for time stepping beginning with solution x and timestep dt.
size_t bandwidth() const
Jacobian bandwidth.
double m_rdt
Reciprocal of time step.
shared_ptr< SystemJacobian > m_jac
Jacobian evaluator.
shared_ptr< vector< double > > m_state
Solution vector.
vector< int > m_mask
Transient mask.
void solve(int loglevel=0)
Solve the steady-state problem, taking internal timesteps as necessary until the Newton solver can co...
vector< int > m_steps
Array of number of steps to take after each unsuccessful steady-state solve before re-attempting the ...
double m_tfactor
Factor time step is multiplied by if time stepping fails ( < 1 )
bool m_jac_ok
If true, Jacobian is current.
double m_tstep
Initial timestep.
int m_ts_jac_age
Maximum age of the Jacobian in time-stepping mode.
virtual void eval(double *x, double *r, double rdt=-1.0, int count=1)=0
Evaluate the residual function.
double timeStep(int nsteps, double dt, double *x, double *r, int loglevel)
Take time steps using Backward Euler.
void setJacAge(int ss_age, int ts_age=-1)
Set the maximum number of steps that can be taken using the same Jacobian before it must be re-evalua...
void evalSSJacobian(double *x, double *rsd)
Evaluate the steady-state Jacobian, accessible via linearSolver()
virtual void writeDebugInfo(const string &header_suffix, const string &message, int loglevel, int attempt_counter)
Write solver debugging based on the specified log level.
double m_tmin
Minimum timestep size.
double m_tmax
Maximum timestep size.
int m_ss_jac_age
Maximum age of the Jacobian in steady-state mode.
virtual void setSteadyMode()
Prepare to solve the steady-state problem.
virtual void clearDebugFile()
Delete debug output file that may be created by writeDebugInfo() when solving with high loglevel.
int m_attempt_counter
Counter used to manage the number of states stored in the debug log file generated by writeDebugInfo(...
void setLinearSolver(shared_ptr< SystemJacobian > solver)
Set the linear solver used to hold the Jacobian matrix and solve linear systems as part of each Newto...
void setTimeStep(double stepsize, size_t n, const int *tsteps)
Set the number of time steps to try when the steady Newton solver is unsuccessful.
int m_max_history
Constant that determines the maximum number of states stored in the debug log file generated by write...
void setInitialGuess(const double *x)
Set the initial guess. Should be called before solve().
void getState(double *x) const
Get the converged steady-state solution after calling solve().
Func1 * m_time_step_callback
User-supplied function called after each successful timestep.
virtual void evalJacobian(double *x0)=0
Evaluates the Jacobian at x0 using finite differences.
vector< double > m_xlast_ts
State vector after the last successful set of time steps.
MultiNewton & newton()
Return a reference to the Newton iterator.
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
Definition global.h:154
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition global.h:171
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const double Tiny
Small number to compare differences of mole fractions against.
Definition ct_defs.h:173