Cantera  4.0.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
25void SteadyStateSystem::setInitialGuess(span<const double> x)
26{
29 m_state->assign(x.begin(), x.end());
30}
31
32void SteadyStateSystem::getState(span<double> x) const
33{
34 copy(m_xnew.begin(), m_xnew.end(), x.begin());
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) {
54 try {
55 m_jac->updateTransient(m_rdt, m_mask);
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, m_xnew, *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, m_xnew, 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, m_xnew)));
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, span<double> x,
112 span<double> r, int loglevel)
113{
114 // set the Jacobian age parameter to the transient value
116
117 int n = 0;
118 int successiveFailures = 0;
119
120 // Only output this if nothing else under this function call will be output
121 if (loglevel == 1) {
122 writelog("\n============================");
123 writelog("\n{:<5s} {:<11s} {:<7s}\n", "step", "dt (s)", "log(ss)");
124 writelog("============================");
125 }
126 while (n < nsteps) {
127 if (loglevel == 1) { // At level 1, output concise information
128 double ss = ssnorm(x, r);
129 writelog("\n{:<5d} {:<6.4e} {:>7.4f}", n, dt, log10(ss));
130 } else if (loglevel > 1) {
131 double ss = ssnorm(x, r);
132 writelog("\nTimestep ({}) dt= {:<11.4e} log(ss)= {:<7.4f}", n, dt, log10(ss));
133 }
134
135 // set up for time stepping with stepsize dt
136 initTimeInteg(dt, x);
137
138 int j0 = m_jac->nEvals(); // Store the current number of Jacobian evaluations
139
140 // solve the transient problem
141 int status = newton().solve(x, r, *this, loglevel);
142
143 // successful time step. Copy the new solution in r to
144 // the current solution in x.
145 if (status >= 0) {
146 if (loglevel > 1) {
147 writelog("\nTimestep ({}) succeeded", n);
148 }
149 successiveFailures = 0;
150 m_nsteps++;
151 n += 1;
152 copy(r.begin(), r.end(), x.begin());
153 // No Jacobian evaluations were performed, so a larger timestep can be used
154 if (m_jac->nEvals() == j0) {
155 dt *= 1.5;
156 }
159 }
160 dt = std::min(dt, m_tmax);
161 if (m_nsteps >= m_nsteps_max) {
162 throw CanteraError("OneDim::timeStep",
163 "Took maximum number of timesteps allowed ({}) without "
164 "reaching steady-state solution.", m_nsteps_max);
165 }
166 } else {
167 // No solution could be found with this time step.
168 // Decrease the stepsize and try again.
169 successiveFailures++;
170 if (loglevel == 1) {
171 writelog("\nTimestep failed");
172 } else if (loglevel > 1) {
173 writelog("\nTimestep ({}) failed", n);
174 }
175 if (successiveFailures > 2) {
176 debuglog("--> Resetting negative species concentrations", loglevel);
178 successiveFailures = 0;
179 } else {
180 debuglog("--> Reducing timestep", loglevel);
181 dt *= m_tfactor;
182 if (dt < m_tmin) {
183 string err_msg = fmt::format(
184 "Time integration failed. Minimum timestep ({}) reached.", m_tmin);
185 throw CanteraError("OneDim::timeStep", err_msg);
186 }
187 }
188 }
189 }
190
191 // Write the final step to the log
192 if (loglevel == 1) {
193 double ss = ssnorm(x, r);
194 writelog("\n{:<5d} {:<6.4e} {:>7.4f}", n, dt, log10(ss));
195 writelog("\n============================");
196 } else if (loglevel > 1) {
197 double ss = ssnorm(x, r);
198 writelog("\nTimestep ({}) dt= {:<11.4e} log10(ss)= {:<7.4f}\n", n, dt, log10(ss));
199 }
200
201 // return the value of the last stepsize, which may be smaller
202 // than the initial stepsize
203 return dt;
204}
205
206double SteadyStateSystem::ssnorm(span<const double> x, span<double> r)
207{
208 eval(x, r, 0.0, 0);
209 double ss = 0.0;
210 for (size_t i = 0; i < m_size; i++) {
211 ss = std::max(fabs(r[i]),ss);
212 }
213 return ss;
214}
215
216void SteadyStateSystem::setTimeStep(double stepsize, span<const int> tsteps)
217{
218 m_tstep = stepsize;
219 m_steps.assign(tsteps.begin(), tsteps.end());
220}
221
223{
224 m_state->resize(size());
225 m_xnew.resize(size());
226 m_newt->resize(size());
227 m_mask.resize(size());
228 if (!m_jac) {
229 throw CanteraError("SteadyStateSystem::resize",
230 "Jacobian evaluator must be instantiated before calling resize()");
231 }
232 m_jac->initialize(size());
233 m_jac->setBandwidth(bandwidth());
234 m_jac->clearStats();
235 m_jac_ok = false;
236}
237
238void SteadyStateSystem::setLinearSolver(shared_ptr<SystemJacobian> solver)
239{
240 m_jac = solver;
241 m_jac->initialize(size());
242 m_jac->setBandwidth(bandwidth());
243 m_jac->clearStats();
244 m_jac_ok = false;
245}
246
247void SteadyStateSystem::evalSSJacobian(span<const double> x)
248{
249 double rdt_save = m_rdt;
250 m_jac_ok = false;
252 evalJacobian(x);
253 m_rdt = rdt_save;
254}
255
256void SteadyStateSystem::initTimeInteg(double dt, span<const double> x)
257{
258 double rdt_old = m_rdt;
259 m_rdt = 1.0/dt;
260
261 // if the stepsize has changed, then update the transient part of the Jacobian
262 if (fabs(rdt_old - m_rdt) > Tiny) {
263 m_jac->updateTransient(m_rdt, m_mask);
264 }
265}
266
268{
269 if (m_rdt == 0) {
270 return;
271 }
272
273 m_rdt = 0.0;
274 m_jac->updateTransient(m_rdt, m_mask);
275}
276
277void SteadyStateSystem::setJacAge(int ss_age, int ts_age)
278{
279 m_ss_jac_age = ss_age;
280 if (ts_age > 0) {
281 m_ts_jac_age = ts_age;
282 } else {
284 }
285}
286
288{
289 return *m_newt;
290}
291
292}
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
void setOptions(int maxJacAge=5)
Set options.
int solve(span< const double > x0, span< double > x1, SteadyStateSystem &r, int loglevel)
Find the solution to F(x) = 0 by damped Newton iteration.
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...
double timeStep(int nsteps, double dt, span< double > x, span< double > r, int loglevel)
Take time steps using Backward Euler.
virtual void resetBadValues(span< double > x)
Reset values such as negative species concentrations.
vector< double > m_xnew
Work array used to hold the residual or the new solution.
unique_ptr< MultiNewton > m_newt
Newton iterator.
void getState(span< double > x) const
Get the converged steady-state solution after calling solve().
size_t size() const
Total solution vector length;.
virtual void initTimeInteg(double dt, span< const double > x)
Prepare for time stepping beginning with solution x and timestep dt.
void evalSSJacobian(span< const double > x)
Evaluate the steady-state Jacobian, accessible via linearSolver()
size_t bandwidth() const
Jacobian bandwidth.
double m_rdt
Reciprocal of time step.
shared_ptr< SystemJacobian > m_jac
Jacobian evaluator.
void setInitialGuess(span< const double > x)
Set the initial guess. Should be called before solve().
shared_ptr< vector< double > > m_state
Solution vector.
virtual void eval(span< const double > x, span< double > r, double rdt=-1.0, int count=1)=0
Evaluate the residual function.
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...
double ssnorm(span< const double > x, span< double > r)
Steady-state max norm (infinity norm) of the residual evaluated using solution x.
void setTimeStep(double stepsize, span< const int > tsteps)
Set the number of time steps to try when the steady Newton solver is unsuccessful.
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.
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...
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 evalJacobian(span< const double > x0)=0
Evaluates the Jacobian at x0 using finite differences.
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...
int m_max_history
Constant that determines the maximum number of states stored in the debug log file generated by write...
Func1 * m_time_step_callback
User-supplied function called after each successful timestep.
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:176