teqp 0.23.1
Loading...
Searching...
No Matches
iteration.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <optional>
4
6#include "teqp/cpp/derivs.hpp"
7#include "teqp/exceptions.hpp"
8
9#include <boost/multiprecision/cpp_bin_float.hpp>
10
11namespace teqp {
12namespace iteration {
13
14using namespace cppinterface;
15
17
24 const int N;
25 const Eigen::Array2d& x;
26 const Eigen::Array2d& dx;
27 const Eigen::Array2d& r;
28 const std::vector<int>& nonconstant_indices;
29};
30
38public:
40 virtual std::string desc() = 0;
41 virtual ~StoppingCondition() = default;
42};
43
45private:
46 double threshold;
47public:
48 MaxAbsErrorCondition(double threshold) : threshold(threshold) {};
51 return (data.r.abs().maxCoeff() < threshold) ? s::success : s::keep_going;
52 };
53 std::string desc() override{
54 return "MaxAbsErrorCondition";
55 };
56};
57
59private:
60 double threshold;
61public:
62 MinRelStepsizeCondition(double threshold) : threshold(threshold) {};
63
69 double minval = 1e99;
70 for (auto idx : data.nonconstant_indices){
71 double val = std::abs(data.dx(idx)/data.x(idx));
72 minval = std::min(val, minval);
73 }
74 return (minval < threshold) ? s::success : s::keep_going;
75 };
76 std::string desc() override{
77 return "MinRelStepsizeCondition";
78 };
79};
80
82public:
88 bool allfinite = true;
89 for (auto idx : data.nonconstant_indices){
90 if (!std::isfinite(data.x(idx)) || !std::isfinite(data.dx(idx))){
91 allfinite = false; break;
92 }
93 }
94 return (!allfinite) ? s::fatal : s::keep_going;
95 };
96 std::string desc() override{
97 return "NanXDXErrorCondition";
98 };
99};
100
102public:
104 using s = StoppingConditionReason;
105 bool allpositive = true;
106 for (auto idx : data.nonconstant_indices){
107 if (data.x(idx) < 0){
108 allpositive = false; break;
109 }
110 }
111 return (!allpositive) ? s::fatal : s::keep_going;
112 };
113 std::string desc() override{
114 return "NegativeXErrorCondition";
115 };
116};
117
119private:
120 const std::size_t Nmax;
121public:
122 StepCountErrorCondition(int Nmax) : Nmax(Nmax){}
124 using s = StoppingConditionReason;
125 return (data.N >= Nmax) ? s::fatal : s::keep_going;
126 };
127 std::string desc() override{
128 return "StepCountErrorCondition";
129 };
130};
131
132
133
134
136public:
137 std::shared_ptr<AbstractModel> model_ideal_gas, model_residual;
138
139 template<typename Z>
140 auto get_R(const Z& z) const{
141 return model_residual->get_R(z);
142 }
143
144 auto get_R(const std::vector<double>& z) const{
145 const auto zz = Eigen::Map<const Eigen::ArrayXd>(&z[0], z.size());
146 return model_residual->get_R(zz);
147 }
148
149 template<typename Z>
150 Eigen::Array33d get_deriv_mat2(double T, double rho, const Z& z) const{
151 return model_ideal_gas->get_deriv_mat2(T, rho, z) + model_residual->get_deriv_mat2(T, rho, z);
152 }
153
154 Eigen::Array33d get_deriv_mat2(double T, double rho, const std::vector<double>& z) const{
155 const auto zz = Eigen::Map<const Eigen::ArrayXd>(&z[0], z.size());
156 return model_ideal_gas->get_deriv_mat2(T, rho, zz) + model_residual->get_deriv_mat2(T, rho, zz);
157 }
158
159 template<typename Z>
160 auto get_A00A10A01(double T, double rho, const Z& z) const{
161 auto A10 = model_ideal_gas->get_Ar10(T, rho, z) + model_residual->get_Ar10(T, rho, z);
162 Eigen::Array2d A00A01 = model_ideal_gas->get_Ar01n(T, rho, z) + model_residual->get_Ar01n(T, rho, z); // This gives [A00, A01] at the same time
163 return std::make_tuple(A00A01[0], A10, A00A01[1]);
164 }
165
176 template<typename Array>
177 auto get_vals(const std::vector<char>& vars, const double R, const double T, const double rho, const Array &z) const{
178
179 Array v(2);
180 auto Trecip = 1.0/T;
181 auto dTrecipdT = -Trecip*Trecip;
182
183 // Define some lambda functions for things that *might* be needed
184 // The lambdas are used to get a sort of lazy evaluation. The lambdas are
185 // only evaluated on an as-needed basis. If the lambda is not called,
186 // no cost is incurred.
187 //
188 // Probably the compiler will inline these functions anyhow.
189 //
190 auto [A00, A10, A01] = get_A00A10A01(T, rho, z);
191 // Derivatives of total alpha(Trecip, rho)
192 auto alpha = [&](){ return A00; };
193 auto dalphadTrecip = [&](){ return A10/Trecip; };
194 auto dalphadrho = [&](){ return A01/rho; };
195 // Derivatives of total Helmholtz energy a = alpha*R/Trecip in
196 // terms of derivatives of alpha(Trecip, rho)
197 auto a = [&](){ return alpha()*R*T; };
198 auto dadTrecip = [&](){ return R/(Trecip*Trecip)*(Trecip*dalphadTrecip()-alpha()); };
199 auto dadrho = [&](){ return R/Trecip*(dalphadrho()); };
200
201 for (auto i = 0; i < vars.size(); ++i){
202 switch(vars[i]){
203 case 'T':
204 v(i) = T; break;
205 case 'D':
206 v(i) = rho; break;
207 case 'P':
208 v(i) = rho*rho*dadrho(); break;
209 case 'S':
210 v(i) = Trecip*Trecip*dadTrecip(); break;
211 case 'H':
212 v(i) = a() + Trecip*dadTrecip() + rho*dadrho(); break;
213 case 'U':
214 v(i) = a() + Trecip*dadTrecip(); break;
215 default:
216 throw std::invalid_argument("bad var: " + std::to_string(vars[i]));
217 }
218 }
219 return v;
220 }
221};
222static_assert(std::is_copy_constructible_v<AlphaModel>);
223
227 std::string msg;
228 Eigen::Array2d Trho;
229};
230
236private:
237 const AlphaModel alphamodel;
238 const std::vector<char> vars;
239 const Eigen::Array2d vals;
240
241 const Eigen::ArrayXd z;
242 Eigen::Array2d Trho, r;
243 double R;
244
245 std::tuple<bool, bool> relative_error;
246 const std::vector<std::shared_ptr<StoppingCondition>> stopping_conditions;
247
248 int step_counter = 0;
249 std::string msg;
250 std::vector<int> nonconstant_indices;
251 bool isTD;
252
253public:
254 NRIterator(const AlphaModel& alphamodel, const std::vector<char>& vars, const Eigen::Array2d& vals, double T, double rho, const Eigen::Ref<const Eigen::ArrayXd>& z, const std::tuple<bool, bool> &relative_error, const std::vector<std::shared_ptr<StoppingCondition>> stopping_conditions) : alphamodel(alphamodel), vars(vars), vals(vals), Trho((Eigen::Array2d() << T,rho).finished()), z(z), R(alphamodel.get_R(z)), relative_error(relative_error), stopping_conditions(stopping_conditions) {
255 if(!(vars[0] == 'T' || vars[1] == 'T')){ nonconstant_indices.push_back(0); }
256 if(!(vars[0] == 'D' || vars[1] == 'D')){ nonconstant_indices.push_back(1); }
257 isTD = (nonconstant_indices.size() == 0);
258 }
259 bool verbose = false;
260
262 std::vector<char> get_vars() const { return vars; }
264 Eigen::Array2d get_vals() const { return vals; }
266 auto get_T() const { return Trho(0); }
268 auto get_rho() const { return Trho(1); }
270 Eigen::ArrayXd get_molefrac() const { return z; }
272 int get_step_count() const { return step_counter; }
274 std::string get_message() const { return msg; }
276 auto get_nonconstant_indices() const { return nonconstant_indices; }
277
278 void reset(double T, double rho){
279 Trho = (Eigen::Array2d() << T, rho).finished();
280 step_counter = 0;
281 }
282
288 auto calc_matrices(double T, double rho) const{
289 auto A = alphamodel.get_deriv_mat2(T, rho, z);
290 return build_iteration_Jv(vars, A, R, T, rho, z);
291 }
292 auto calc_step(double T, double rho) const{
293 auto im = calc_matrices(T, rho);
294 return std::make_tuple((im.J.matrix().fullPivLu().solve((-(im.v-vals)).matrix())).eval(), im);
295 }
296 auto calc_just_step(double T, double rho) const {
297 return std::get<0>(calc_step(T, rho));
298 }
299 auto calc_vals(double T, double rho) const{
300 return calc_matrices(T, rho).v;
301 }
302 auto calc_r(double T, double rho) const{
303 auto im = calc_matrices(T, rho);
304 Eigen::Array2d r = im.v-vals;
305 if (std::get<0>(relative_error)){ r(0) /= vals(0);}
306 if (std::get<1>(relative_error)){ r(1) /= vals(1);}
307 return r;
308 }
309 auto calc_J(double T, double rho) const{
310 auto im = calc_matrices(T, rho);
311 if (std::get<0>(relative_error)){ im.J.row(0) /= vals(0);}
312 if (std::get<1>(relative_error)){ im.J.row(1) /= vals(1);}
313 return im.J;
314 }
315
316 auto calc_maxabsr(double T, double rho) const{
317 Eigen::Array2d r = alphamodel.get_vals(vars, R, T, rho, z)-vals;
318 if (std::get<0>(relative_error)){ r(0) /= vals(0);}
319 if (std::get<1>(relative_error)){ r(1) /= vals(1);}
320 return r.abs().maxCoeff();
321 }
322
323 auto get_maxabsr() const{
324 return r.abs().maxCoeff();
325 }
326
331 auto take_steps(int N, bool apply_stopping=true){
332 if (N <= 0){
333 throw teqp::InvalidArgument("N must be greater than 0");
334 }
336 if(isTD){
338 auto im = calc_matrices(Trho(0), Trho(1));
339 r = im.v-vals;
340
341 if (std::get<0>(relative_error)){ r(0) /= vals(0); im.J.row(0) /= vals(0); }
342 if (std::get<1>(relative_error)){ r(1) /= vals(1); im.J.row(1) /= vals(1); }
343
344 Eigen::Array2d dTrho = im.J.matrix().fullPivLu().solve((-r).matrix());
345 Trho += dTrho;
346 step_counter++;
347 reason = StoppingConditionReason::success; msg = "Only one step is needed for DT inputs";
348 return reason;
349 }
350
351 for (auto K = 0; K < N; ++K){
352 auto im = calc_matrices(Trho(0), Trho(1));
353 r = im.v-vals;
354
355 if (std::get<0>(relative_error)){ r(0) /= vals(0); im.J.row(0) /= vals(0); }
356 if (std::get<1>(relative_error)){ r(1) /= vals(1); im.J.row(1) /= vals(1); }
357
358 Eigen::Array2d dTrho = im.J.matrix().fullPivLu().solve((-r).matrix());
359
360 bool stop = false;
361 if (apply_stopping){
362 // Check whether a stopping condition (either good[complete] or bad[error])
363 const StoppingData data{K, Trho, dTrho, r, nonconstant_indices};
364 for (auto& condition : stopping_conditions){
365 using s = StoppingConditionReason;
366 auto this_reason = condition->stop(data);
367 if (this_reason != s::keep_going){
368 stop = true; reason = this_reason; msg = condition->desc(); break;
369 }
370 }
371 }
372
373 Trho += dTrho;
374 step_counter++;
375 if (stop){
376 break;
377 }
378 }
379 return reason;
380 }
381
382
386 auto take_steps_logrho(int N){
387 if (N <= 0){
388 throw teqp::InvalidArgument("N must be greater than 0");
389 }
390 auto tic = std::chrono::steady_clock::now();
392
393 for (auto K = 0; K < N; ++K){
394 auto im = calc_matrices(Trho(0), Trho(1));
395 r = im.v-vals;
396
397 im.J.col(1) *= Trho(1); // This will make the step be [dT, dln(rho)] instead of [dT, drho]
398 if (std::get<0>(relative_error)){ r(0) /= vals(0); im.J.row(0) /= vals(0); }
399 if (std::get<1>(relative_error)){ r(1) /= vals(1); im.J.row(1) /= vals(1); }
400
401// Eigen::JacobiSVD<Eigen::Matrix2d> svd(im.J);
402// double cond = svd.singularValues()(0) / svd.singularValues()(svd.singularValues().size()-1);
403
404 Eigen::Array2d step;
405 if (false){//cond > 1){
406 // The idea in this block was to use extended precision to obtain the step
407 // but it seems to make no difference at all, because the condition number is a
408 // function of the problem, but the loss in precision is only meaningful relative to
409 // the number of digits of working precision available via the lost digits approximated
410 // by log10(kappa), where kappa is the condition number
411 using my_float = boost::multiprecision::number<boost::multiprecision::cpp_bin_float<164U>>;
412 Eigen::Matrix<my_float, 2, 2> Jmp = im.J.cast<my_float>();
413 Eigen::Vector<my_float, 2> rmp = r.cast<my_float>();
414 step = Jmp.fullPivLu().solve(-rmp).cast<double>();
415 Eigen::Array<my_float, 2, -1> newvals = (Trho.cast<my_float>() + Jmp.fullPivLu().solve(-rmp).array());
416
417 Eigen::Array2d new_double = (Trho + im.J.matrix().fullPivLu().solve((-r).matrix()).array());
418 std::cout << "-- " << ((newvals-new_double.cast<my_float>())/Trho.cast<my_float>()).cast<double>() << std::endl;
419 }
420 else{
421 step = im.J.matrix().fullPivLu().solve((-r).matrix());
422 }
423 Eigen::Array2d dTrho = (Eigen::Array2d() << step(0), Trho(1)*(exp(step(1))-1)).finished();
424
425 // Check whether a stopping condition (either good[complete] or bad[error])
426 bool stop = false;
427 const StoppingData data{K, Trho, dTrho, r, nonconstant_indices};
428 for (auto& condition : stopping_conditions){
429 using s = StoppingConditionReason;
430 auto this_reason = condition->stop(data);
431 if (this_reason != s::keep_going){
432 stop = true; reason = this_reason; msg = condition->desc(); break;
433 }
434 }
435
436 Trho += dTrho;
437 step_counter++;
438
439 if(nonconstant_indices.size() == 0){
440 stop = true; reason = StoppingConditionReason::success; msg = "Only one step is needed for DT inputs";
441 }
442// if(nonconstant_indices.size() == 2){
443// // If the step size gets too small on a relative basis then you
444// // need to stop
445// // There is a loss in precision caused by the Jacobian having a large condition number
446// // and use that to determine whether you need to stop stepping
447// //
448// // The condition number only is meaningful for 2x2 systems where both variables
449// // are being iterated for.
450// if ((dTrho/Trho)(nonconstant_indices).abs().minCoeff() < 1e-15*cond){
451// std::stringstream ss; ss << im.J;
452// stop = true; reason = StoppingConditionReason::success; msg = "MinRelStepSize ~= cond of " + std::to_string(cond) + ": J" + ss.str();
453// }
454// }
455 if (verbose){
456// std::cout << step_counter << "," << r.abs().maxCoeff() << "," << Trho << "|" << (dTrho/Trho).abs().minCoeff() << "cond: " << cond << std::endl;
457// std::cout << step_counter << "," << r.abs().maxCoeff() << "," << Trho << "|" << (dTrho/Trho).abs().minCoeff() << im.J.matrix() << im.v-vals << "cond: " << cond << std::endl;
458 }
459 if (stop){
460 break;
461 }
462 }
463 auto toc = std::chrono::steady_clock::now();
464// std::cout << "Time difference = " << std::chrono::duration_cast<std::chrono::nanoseconds>(toc - tic).count()/1e3/step_counter << "[µs/call]" << std::endl;
465 return reason;
466 }
467 auto path_integration(double T, double rho, std::size_t N){
468 // Start at the given value of T, rho
469
470 // Calculate the given variables
471 auto im = calc_matrices(T, rho);
472
473 Eigen::Array2d vals_current = im.v;
474 double x_current = vals_current(0);
475 double y_current = vals_current(1);
476
477 // Assume (for now) a linear path between starting point and the target point
478 double x_target = vals(0);
479 double y_target = vals(1);
480 auto dxdt = (x_target-x_current)/(1.0-0.0); // The 1.0 to make clear what is going on; we are integrating from 0 -> 1
481 auto dydt = (y_target-y_current)/(1.0-0.0);
482 double dt = 1.0/N;
483
484 // Determine the value of the residual function at the end point
485 for (auto i = 0U; i < N; ++i){
486
487 auto dxdT__rho = im.J(0, 0);
488 auto dxdrho__T = im.J(0, 1);
489 auto dydT__rho = im.J(1, 0);
490 auto dydrho__T = im.J(1, 1);
491
492 auto dxdT__y = dxdT__rho - dxdrho__T*dydT__rho/dydrho__T;
493 auto dxdrho__y = dxdrho__T - dxdT__rho*dydrho__T/dydT__rho;
494
495 auto dydT__x = dydT__rho - dydrho__T*dxdT__rho/dxdrho__T;
496 auto dydrho__x = dydrho__T - dydT__rho*dxdrho__T/dxdT__rho;
497
498 // Calculate the increment in temperature and density along the integration path (needs to be adaptive eventually)
499 double dTdt = dxdt/dxdT__y + dydt/dydT__x;
500 double drhodt = dxdt/dxdrho__y + dydt/dydrho__x;
501
502 T += dTdt*dt;
503 rho += drhodt*dt;
504 im = calc_matrices(T, rho);
505 if (verbose){
506 std::cout << i << "," << (im.v-vals).abs().maxCoeff() << std::endl;
507 }
508 }
509 return std::make_tuple(T, rho, im.v(0), im.v(1));
510 }
511};
512
513
514}
515}
std::shared_ptr< AbstractModel > model_residual
auto get_A00A10A01(double T, double rho, const Z &z) const
std::shared_ptr< AbstractModel > model_ideal_gas
auto get_R(const Z &z) const
auto get_R(const std::vector< double > &z) const
auto get_vals(const std::vector< char > &vars, const double R, const double T, const double rho, const Array &z) const
A convenience function for calculation of Jacobian terms of the form and where is one of the therm...
Eigen::Array33d get_deriv_mat2(double T, double rho, const std::vector< double > &z) const
Eigen::Array33d get_deriv_mat2(double T, double rho, const Z &z) const
StoppingConditionReason stop(const StoppingData &data) override
Definition iteration.hpp:49
StoppingConditionReason stop(const StoppingData &data) override
Definition iteration.hpp:67
Eigen::ArrayXd get_molefrac() const
Return the current mole fractions.
auto calc_maxabsr(double T, double rho) const
Calculate the maximum absolute value of the values in r.
std::vector< char > get_vars() const
Return the variables that are being used in the iteration.
auto calc_J(double T, double rho) const
auto get_T() const
Return the current temperature.
auto path_integration(double T, double rho, std::size_t N)
auto get_rho() const
Return the current molar density.
auto calc_vals(double T, double rho) const
std::string get_message() const
Return the current message relaying information about success or failure of the iteration.
auto take_steps(int N, bool apply_stopping=true)
auto calc_just_step(double T, double rho) const
auto get_nonconstant_indices() const
Get the indices of the nonconstant variables.
auto calc_step(double T, double rho) const
Eigen::Array2d get_vals() const
Return the target values to be obtained.
auto calc_matrices(double T, double rho) const
void reset(double T, double rho)
NRIterator(const AlphaModel &alphamodel, const std::vector< char > &vars, const Eigen::Array2d &vals, double T, double rho, const Eigen::Ref< const Eigen::ArrayXd > &z, const std::tuple< bool, bool > &relative_error, const std::vector< std::shared_ptr< StoppingCondition > > stopping_conditions)
auto get_maxabsr() const
Get the maximum absolute value of residual vector, using the cached value for r.
int get_step_count() const
Return the current step counter.
auto calc_r(double T, double rho) const
StoppingConditionReason stop(const StoppingData &data) override
Definition iteration.hpp:86
StoppingConditionReason stop(const StoppingData &data) override
StoppingConditionReason stop(const StoppingData &data) override
virtual StoppingConditionReason stop(const StoppingData &)=0
virtual ~StoppingCondition()=default
virtual std::string desc()=0
auto build_iteration_Jv(const std::vector< char > &vars, const Eigen::Array< double, 3, 3 > &A, const double R, const double T, const double rho, const Array &z)
A convenience function for calculation of Jacobian terms of the form and where is one of the therm...
Definition derivs.hpp:26
StoppingConditionReason reason
const Eigen::Array2d & r
The set of residual functions.
Definition iteration.hpp:27
const std::vector< int > & nonconstant_indices
The indices where the independent variable can be varied (are not T or rho)
Definition iteration.hpp:28
const Eigen::Array2d & dx
The complete set of steps in independent variables.
Definition iteration.hpp:26
const int N
The number of steps that have been taken so far.
Definition iteration.hpp:24
const Eigen::Array2d & x
The complete set of independent variables.
Definition iteration.hpp:25