teqp 0.19.1
Loading...
Searching...
No Matches
derivs.hpp
Go to the documentation of this file.
1#pragma once
2
3#if defined(TEQP_COMPLEXSTEP_ENABLED)
4#include <complex>
5#endif
6#include <map>
7#include <tuple>
8#include <numeric>
9#include <concepts>
10
11#include "teqp/types.hpp"
12#include "teqp/exceptions.hpp"
13
14#if defined(TEQP_MULTICOMPLEX_ENABLED)
15#include "MultiComplex/MultiComplex.hpp"
16#endif
17
18// autodiff include
19#include <autodiff/forward/dual.hpp>
20#include <autodiff/forward/dual/eigen.hpp>
21#include <autodiff/forward/real.hpp>
22#include <autodiff/forward/real/eigen.hpp>
23
24using namespace autodiff;
25
26namespace teqp {
27
32template <typename TType, typename ContainerType, typename FuncType>
33typename ContainerType::value_type derivT(const FuncType& f, TType T, const ContainerType& rho) {
34 double h = 1e-100;
35 return f(std::complex<TType>(T, h), rho).imag() / h;
36}
37
38#if defined(TEQP_MULTICOMPLEX_ENABLED)
43template <typename TType, typename ContainerType, typename FuncType>
44typename ContainerType::value_type derivTmcx(const FuncType& f, TType T, const ContainerType& rho) {
45 using fcn_t = std::function<mcx::MultiComplex<double>(const mcx::MultiComplex<double>&)>;
46 fcn_t wrapper = [&rho, &f](const auto& T_) {return f(T_, rho); };
47 auto ders = diff_mcx1(wrapper, T, 1);
48 return ders[0];
49}
50#endif
51
56template <typename TType, typename ContainerType, typename FuncType, typename Integer>
57typename ContainerType::value_type derivrhoi(const FuncType& f, TType T, const ContainerType& rho, Integer i) {
58 double h = 1e-100;
59 using comtype = std::complex<typename ContainerType::value_type>;
60 Eigen::ArrayX<comtype> rhocom(rho.size());
61 for (auto j = 0; j < rho.size(); ++j) {
62 rhocom[j] = comtype(rho[j], 0.0);
63 }
64 rhocom[i] = comtype(rho[i], h);
65 return f(T, rhocom).imag() / h;
66}
67
70template<typename T, size_t ... I>
71auto build_duplicated_tuple_impl(const T& val, std::index_sequence<I...>) {
72 return std::make_tuple((static_cast<void>(I), val)...); // The comma operator (a,b) evaluates the first argument and discards it, and keeps the second argument
73}
74
76template<int N, typename T>
77auto build_duplicated_tuple(const T& val) {
78 return build_duplicated_tuple_impl(val, std::make_index_sequence<N>());
79}
80
84struct wrt_helper {
85 template<typename... Args>
86 auto operator()(Args&&... args) const
87 {
88 return Wrt<Args&&...>{ std::forward_as_tuple(std::forward<Args>(args)...) };
89 }
90};
91
92enum class ADBackends { autodiff
93#if defined(TEQP_MULTICOMPLEX_ENABLED)
94 ,multicomplex
95#endif
96#if defined(TEQP_COMPLEXSTEP_ENABLED)
97 ,complex_step
98#endif
99};
100
101template<typename T, typename U, typename V, typename W>
102concept CallableAlpha = requires(T t, U u, V v, W w) {
103 { t.alpha(u,v,w) };
104};
105
106template<typename T, typename U, typename V, typename W>
107concept CallableAlphar = requires(T t, U u, V v, W w) {
108 { t.alphar(u,v,w) };
109};
110
111template<typename T, typename U, typename V, typename W>
112concept CallableAlpharTauDelta = requires(T t, U u, V v, W w) {
113 { t.alphar_taudelta(u,v,w) };
114};
115
116template<typename Model, typename Scalar = double, typename VectorType = Eigen::ArrayXd>
118
119 static auto get_Ar00(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
120 return model.alphar(T, rho, molefrac);
121 }
122
123 template<typename AlphaWrapper, typename S1, typename S2, typename Vec>
124 static auto AlphaCaller(const AlphaWrapper& w, const S1& T, const S2& rho, const Vec& molefrac) requires CallableAlpha<AlphaWrapper, S1, S2, Vec>{
125 return w.alpha(T, rho, molefrac);
126 }
127 template<typename AlphaWrapper, typename S1, typename S2, typename Vec>
128 static auto AlphaCaller(const AlphaWrapper& w, const S1& T, const S2& rho, const Vec& molefrac) requires CallableAlphar<AlphaWrapper, S1, S2, Vec>{
129 return w.alphar(T, rho, molefrac);
130 }
131
132 template<typename AlphaWrapper, typename S1, typename S2, typename Vec>
133 static auto AlpharTauDeltaCaller(const AlphaWrapper& w, const S1& T, const S2& rho, const Vec& molefrac) requires CallableAlpharTauDelta<AlphaWrapper, S1, S2, Vec>{
134 return w.alphar_taudelta(T, rho, molefrac);
135 }
136 template<typename AlphaWrapper, typename S1, typename S2, typename Vec>
137 static auto AlpharTauDeltaCaller(const AlphaWrapper& , const S1&, const S2&, const Vec& molefrac){
138 throw teqp::NotImplementedError("Cannot take derivatives of a class that doesn't define the alphar_taudelta method");
139 return std::common_type_t<S1, S2, decltype(molefrac[0])>(1e99);
140 }
141
150 template<int iT, int iD, ADBackends be = ADBackends::autodiff, class AlphaWrapper>
151 static auto get_Agenxy(const AlphaWrapper& w, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
152
153 if constexpr (iT == 0 && iD == 0){
154 return AlphaCaller(w, T, rho, molefrac);
155 }
156 else if constexpr (iT == 0 && iD > 0) {
157 if constexpr (be == ADBackends::autodiff) {
158 // If a pure derivative, then we can use autodiff::Real for that variable and Scalar for other variable
159 autodiff::Real<iD, Scalar> rho_ = rho;
160 auto f = [&w, &T, &molefrac](const auto& rho__) { return AlphaCaller(w, T, rho__, molefrac); };
161 return powi(rho, iD)*derivatives(f, along(1), at(rho_))[iD];
162 }
163#if defined(TEQP_COMPLEXSTEP_ENABLED)
164 else if constexpr (iD == 1 && be == ADBackends::complex_step) {
165 double h = 1e-100;
166 auto rho_ = std::complex<Scalar>(rho, h);
167 return powi(rho, iD) * AlphaCaller(w, T, rho_, molefrac).imag() / h;
168 }
169#endif
170#if defined(TEQP_MULTICOMPLEX_ENABLED)
171 else if constexpr (be == ADBackends::multicomplex) {
172 using fcn_t = std::function<mcx::MultiComplex<Scalar>(const mcx::MultiComplex<Scalar>&)>;
173 fcn_t f = [&](const auto& rhomcx) { return AlphaCaller(w, T, rhomcx, molefrac); };
174 auto ders = diff_mcx1(f, rho, iD, true /* and_val */);
175 return powi(rho, iD)*ders[iD];
176 }
177#endif
178 else {
179 throw std::invalid_argument("algorithmic differentiation backend is invalid in get_Agenxy for iT == 0");
180 }
181 }
182 else if constexpr (iT > 0 && iD == 0) {
183 Scalar Trecip = 1.0 / T;
184 if constexpr (be == ADBackends::autodiff) {
185 // If a pure derivative, then we can use autodiff::Real for that variable and Scalar for other variable
186 autodiff::Real<iT, Scalar> Trecipad = Trecip;
187 auto f = [&w, &rho, &molefrac](const auto& Trecip__) {return AlphaCaller(w, forceeval(1.0/Trecip__), rho, molefrac); };
188 return powi(Trecip, iT)*derivatives(f, along(1), at(Trecipad))[iT];
189 }
190#if defined(TEQP_COMPLEXSTEP_ENABLED)
191 else if constexpr (iT == 1 && be == ADBackends::complex_step) {
192 double h = 1e-100;
193 auto Trecipcsd = std::complex<Scalar>(Trecip, h);
194 return powi(Trecip, iT)* w.alpha(1/Trecipcsd, rho, molefrac).imag()/h;
195 }
196#endif
197#if defined(TEQP_MULTICOMPLEX_ENABLED)
198 else if constexpr (be == ADBackends::multicomplex) {
199 using fcn_t = std::function<mcx::MultiComplex<Scalar>(const mcx::MultiComplex<Scalar>&)>;
200 fcn_t f = [&](const auto& Trecipmcx) { return AlphaCaller(w, 1.0/Trecipmcx, rho, molefrac); };
201 auto ders = diff_mcx1(f, Trecip, iT, true /* and_val */);
202 return powi(Trecip, iT)*ders[iT];
203 }
204#endif
205 else {
206 throw std::invalid_argument("algorithmic differentiation backend is invalid in get_Agenxy for iD == 0");
207 }
208 }
209 else { // iT > 0 and iD > 0
210 if constexpr (be == ADBackends::autodiff) {
211 using adtype = autodiff::HigherOrderDual<iT + iD, double>;
212 adtype Trecipad = 1.0 / T, rhoad = rho;
213 auto f = [&w, &molefrac](const adtype& Trecip, const adtype& rho_) {
214 adtype T_ = 1.0/Trecip;
215 return eval(AlphaCaller(w, T_, rho_, molefrac)); };
216 auto wrts = std::tuple_cat(build_duplicated_tuple<iT>(std::ref(Trecipad)), build_duplicated_tuple<iD>(std::ref(rhoad)));
217 auto der = derivatives(f, std::apply(wrt_helper(), wrts), at(Trecipad, rhoad));
218 return powi(forceeval(1.0 / T), iT) * powi(rho, iD) * der[der.size() - 1];
219 }
220#if defined(TEQP_MULTICOMPLEX_ENABLED)
221 else if constexpr (be == ADBackends::multicomplex) {
222 using fcn_t = std::function< mcx::MultiComplex<double>(const std::valarray<mcx::MultiComplex<double>>&)>;
223 const fcn_t func = [&w, &molefrac](const auto& zs) {
224 auto Trecip = zs[0], rhomolar = zs[1];
225 return AlphaCaller(w, 1.0 / Trecip, rhomolar, molefrac);
226 };
227 std::vector<double> xs = { 1.0 / T, rho};
228 std::vector<int> order = { iT, iD };
229 auto der = mcx::diff_mcxN(func, xs, order);
230 return powi(1.0 / T, iT)*powi(rho, iD)*der;
231 }
232#endif
233 else {
234 throw std::invalid_argument("algorithmic differentiation backend is invalid in get_Agenxy for iD > 0 and iT > 0");
235 }
236 }
237 // return static_cast<Scalar>(-999999999*T); // This will never hit, only to make compiler happy because it doesn't know the return type
238 }
239
247 template<int iT, int iD, int iXi, typename AlphaWrapper>
248 static auto get_ATrhoXi(const AlphaWrapper& w, const Scalar& T, const Scalar& rho, const VectorType& molefrac, int i){
249 using adtype = autodiff::HigherOrderDual<iT + iD + iXi, double>;
250 adtype Trecipad = 1.0 / T, rhoad = rho, xi = molefrac[i];
251 auto f = [&w, &molefrac, &i](const adtype& Trecip, const adtype& rho_, const adtype& xi_) {
252 adtype T_ = 1.0/Trecip;
253 Eigen::ArrayX<adtype> molefracdual = molefrac.template cast<adtype>();
254 molefracdual[i] = xi_;
255 return eval(AlphaCaller(w, T_, rho_, molefracdual)); };
256 auto wrts = std::tuple_cat(build_duplicated_tuple<iT>(std::ref(Trecipad)), build_duplicated_tuple<iD>(std::ref(rhoad)), build_duplicated_tuple<iXi>(std::ref(xi)));
257 auto der = derivatives(f, std::apply(wrt_helper(), wrts), at(Trecipad, rhoad, xi));
258 return powi(forceeval(1.0 / T), iT) * powi(rho, iD) * der[der.size() - 1];
259 }
260
261 #define get_ATrhoXi_runtime_combinations \
262 X(0,0,1) \
263 X(0,0,2) \
264 X(0,0,3) \
265 X(1,0,1) \
266 X(1,0,2) \
267 X(1,0,3) \
268 X(0,1,1) \
269 X(0,1,2) \
270 X(0,1,3) \
271 X(2,0,1) \
272 X(2,0,2) \
273 X(2,0,3) \
274 X(1,1,1) \
275 X(1,1,2) \
276 X(1,1,3) \
277 X(0,2,1) \
278 X(0,2,2) \
279 X(0,2,3)
280
281 template<typename AlphaWrapper>
282 static auto get_ATrhoXi_runtime(const AlphaWrapper& w, const Scalar& T, int iT, const Scalar& rho, int iD, const VectorType& molefrac, int i, int iXi){
283 #define X(a,b,c) if (iT == a && iD == b && iXi == c) { return get_ATrhoXi<a,b,c>(w, T, rho, molefrac, i); }
285 #undef X
286 throw teqp::InvalidArgument("Can't match these derivative counts");
287 }
288
289 #define get_ATrhoXiXj_runtime_combinations \
290 X(0,0,1,0) \
291 X(0,0,2,0) \
292 X(0,0,0,1) \
293 X(0,0,0,2) \
294 X(0,0,1,1) \
295 X(1,0,1,0) \
296 X(1,0,2,0) \
297 X(1,0,0,1) \
298 X(1,0,0,2) \
299 X(1,0,1,1) \
300 X(0,1,1,0) \
301 X(0,1,2,0) \
302 X(0,1,0,1) \
303 X(0,1,0,2) \
304 X(0,1,1,1)
305
306 template<typename AlphaWrapper>
307 static auto get_ATrhoXiXj_runtime(const AlphaWrapper& w, const Scalar& T, int iT, const Scalar& rho, int iD, const VectorType& molefrac, int i, int iXi, int j, int iXj){
308 #define X(a,b,c,d) if (iT == a && iD == b && iXi == c && iXj == d) { return get_ATrhoXiXj<a,b,c,d>(w, T, rho, molefrac, i, j); }
310 #undef X
311 throw teqp::InvalidArgument("Can't match these derivative counts");
312 }
313
314 #define get_ATrhoXiXjXk_runtime_combinations \
315 X(0,0,0,1,1) \
316 X(0,0,1,0,1) \
317 X(0,0,1,1,0) \
318 X(0,0,1,1,1) \
319 X(1,0,0,1,1) \
320 X(1,0,1,0,1) \
321 X(1,0,1,1,0) \
322 X(1,0,1,1,1) \
323 X(0,1,0,1,1) \
324 X(0,1,1,0,1) \
325 X(0,1,1,1,0) \
326 X(0,1,1,1,1)
327
328 template<typename AlphaWrapper>
329 static auto get_ATrhoXiXjXk_runtime(const AlphaWrapper& w, const Scalar& T, int iT, const Scalar& rho, int iD, const VectorType& molefrac, int i, int iXi, int j, int iXj, int k, int iXk){
330 #define X(a,b,c,d,e) if (iT == a && iD == b && iXi == c && iXj == d && iXk == e) { return get_ATrhoXiXjXk<a,b,c,d,e>(w, T, rho, molefrac, i, j, k); }
332 #undef X
333 throw teqp::InvalidArgument("Can't match these derivative counts");
334 }
335
336 template<int iT, int iD, int iXi, typename AlphaWrapper>
337 static auto get_AtaudeltaXi(const AlphaWrapper& w, const Scalar& tau, const Scalar& delta, const VectorType& molefrac, const int i) {
338 using adtype = autodiff::HigherOrderDual<iT + iD + iXi, double>;
339 adtype tauad = tau, deltaad = delta, xi = molefrac[i];
340 auto f = [&w, &molefrac, &i](const adtype& tau_, const adtype& delta_, const adtype& xi_) {
341 Eigen::ArrayX<adtype> molefracdual = molefrac.template cast<adtype>();
342 molefracdual[i] = xi_;
343 return eval(AlpharTauDeltaCaller(w, tau_, delta_, molefracdual)); };
344 auto wrts = std::tuple_cat(build_duplicated_tuple<iT>(std::ref(tauad)), build_duplicated_tuple<iD>(std::ref(deltaad)), build_duplicated_tuple<iXi>(std::ref(xi)));
345 auto der = derivatives(f, std::apply(wrt_helper(), wrts), at(tauad, deltaad, xi));
346 return powi(tau, iT) * powi(delta, iD) * der[der.size() - 1];
347 }
348
349 template<int iT, int iD, int iXi, int iXj, typename AlphaWrapper>
350 static auto get_AtaudeltaXiXj(const AlphaWrapper& w, const Scalar& tau, const Scalar& delta, const VectorType& molefrac, const int i, const int j) {
351 using adtype = autodiff::HigherOrderDual<iT + iD + iXi + iXj, double>;
352 if (i == j){
353 throw teqp::InvalidArgument("i cannot equal j");
354 }
355 adtype tauad = tau, deltaad = delta, xi = molefrac[i], xj = molefrac[j];
356 auto f = [&w, &molefrac, i, j](const adtype& tau_, const adtype& delta_, const adtype& xi_, const adtype& xj_) {
357 Eigen::ArrayX<adtype> molefracdual = molefrac.template cast<adtype>();
358 molefracdual[i] = xi_;
359 molefracdual[j] = xj_;
360 return eval(AlpharTauDeltaCaller(w, tau_, delta_, molefracdual)); };
361 auto wrts = std::tuple_cat(build_duplicated_tuple<iT>(std::ref(tauad)), build_duplicated_tuple<iD>(std::ref(deltaad)), build_duplicated_tuple<iXi>(std::ref(xi)), build_duplicated_tuple<iXj>(std::ref(xj)));
362 auto der = derivatives(f, std::apply(wrt_helper(), wrts), at(tauad, deltaad, xi, xj));
363 return powi(tau, iT) * powi(delta, iD) * der[der.size() - 1];
364 }
365
366 template<int iT, int iD, int iXi, int iXj, int iXk, typename AlphaWrapper>
367 static auto get_AtaudeltaXiXjXk(const AlphaWrapper& w, const Scalar& tau, const Scalar& delta, const VectorType& molefrac, const int i, const int j, const int k) {
368 using adtype = autodiff::HigherOrderDual<iT + iD + iXi + iXj + iXk, double>;
369 if (i == j || j == k || i == k){
370 throw teqp::InvalidArgument("i, j, and k must all be unique");
371 }
372 adtype tauad = tau, deltaad = delta, xi = molefrac[i], xj = molefrac[j], xk = molefrac[k];
373 auto f = [&w, &molefrac, i, j, k](const adtype& tau_, const adtype& delta_, const adtype& xi_, const adtype& xj_, const adtype& xk_) {
374 Eigen::ArrayX<adtype> molefracdual = molefrac.template cast<adtype>();
375 molefracdual[i] = xi_;
376 molefracdual[j] = xj_;
377 molefracdual[k] = xk_;
378 return eval(AlpharTauDeltaCaller(w, tau_, delta_, molefracdual)); };
379 auto wrts = std::tuple_cat(build_duplicated_tuple<iT>(std::ref(tauad)), build_duplicated_tuple<iD>(std::ref(deltaad)), build_duplicated_tuple<iXi>(std::ref(xi)), build_duplicated_tuple<iXj>(std::ref(xj)), build_duplicated_tuple<iXk>(std::ref(xk)));
380 auto der = derivatives(f, std::apply(wrt_helper(), wrts), at(tauad, deltaad, xi, xj, xk));
381 return powi(tau, iT) * powi(delta, iD) * der[der.size() - 1];
382 }
383
384 template<typename AlphaWrapper>
385 static auto get_AtaudeltaXi_runtime(const AlphaWrapper& w, const Scalar& tau, const int iT, const Scalar& delta, const int iD, const VectorType& molefrac, const int i, const int iXi){
386 #define X(a,b,c) if (iT == a && iD == b && iXi == c) { return get_AtaudeltaXi<a,b,c>(w, tau, delta, molefrac, i); }
388 #undef X
389 throw teqp::InvalidArgument("Can't match these derivative counts");
390 }
391
392 template<typename AlphaWrapper>
393 static auto get_AtaudeltaXiXj_runtime(const AlphaWrapper& w, const Scalar& tau, const int iT, const Scalar& delta, const int iD, const VectorType& molefrac, const int i, const int iXi, const int j, const int iXj){
394 #define X(a,b,c,d) if (iT == a && iD == b && iXi == c && iXj == d) { return get_AtaudeltaXiXj<a,b,c,d>(w, tau, delta, molefrac, i, j); }
396 #undef X
397 throw teqp::InvalidArgument("Can't match these derivative counts");
398 }
399
400 template<typename AlphaWrapper>
401 static auto get_AtaudeltaXiXjXk_runtime(const AlphaWrapper& w, const Scalar& tau, const int iT, const Scalar& delta, const int iD, const VectorType& molefrac, const int i, int iXi, const int j, const int iXj, const int k, const int iXk){
402 #define X(a,b,c,d,e) if (iT == a && iD == b && iXi == c && iXj == d && iXk == e) { return get_AtaudeltaXiXjXk<a,b,c,d,e>(w, tau, delta, molefrac, i, j, k); }
404 #undef X
405 throw teqp::InvalidArgument("Can't match these derivative counts");
406 }
407
415 template<int iT, int iD, int iXi, int iXj, typename AlphaWrapper>
416 static auto get_ATrhoXiXj(const AlphaWrapper& w, const Scalar& T, const Scalar& rho, const VectorType& molefrac, int i, int j){
417 if (i == j){
418 throw teqp::InvalidArgument("i cannot equal j");
419 }
420 using adtype = autodiff::HigherOrderDual<iT + iD + iXi + iXj, double>;
421 adtype Trecipad = 1.0 / T, rhoad = rho, xi = molefrac[i], xj = molefrac[j];
422 auto f = [&w, &molefrac, i, j](const adtype& Trecip, const adtype& rho_, const adtype& xi_, const adtype& xj_) {
423 adtype T_ = 1.0/Trecip;
424 Eigen::ArrayX<adtype> molefracdual = molefrac.template cast<adtype>();
425 molefracdual[i] = xi_;
426 molefracdual[j] = xj_;
427 return eval(AlphaCaller(w, T_, rho_, molefracdual)); };
428 auto wrts = std::tuple_cat(build_duplicated_tuple<iT>(std::ref(Trecipad)), build_duplicated_tuple<iD>(std::ref(rhoad)), build_duplicated_tuple<iXi>(std::ref(xi)), build_duplicated_tuple<iXj>(std::ref(xj)));
429 auto der = derivatives(f, std::apply(wrt_helper(), wrts), at(Trecipad, rhoad, xi, xj));
430 return powi(forceeval(1.0 / T), iT) * powi(rho, iD) * der[der.size() - 1];
431 }
432
440 template<int iT, int iD, int iXi, int iXj, int iXk, typename AlphaWrapper>
441 static auto get_ATrhoXiXjXk(const AlphaWrapper& w, const Scalar& T, const Scalar& rho, const VectorType& molefrac, int i, int j, int k){
442 if (i == j || j == k || i == k){
443 throw teqp::InvalidArgument("i, j, and k must all be unique");
444 }
445 using adtype = autodiff::HigherOrderDual<iT + iD + iXi + iXj + iXk, double>;
446 adtype Trecipad = 1.0 / T, rhoad = rho, xi = molefrac[i], xj = molefrac[j], xk = molefrac[k];
447 auto f = [&w, &molefrac, i, j, k](const adtype& Trecip, const adtype& rho_, const adtype& xi_, const adtype& xj_, const adtype& xk_) {
448 adtype T_ = 1.0/Trecip;
449 Eigen::ArrayX<adtype> molefracdual = molefrac.template cast<adtype>();
450 molefracdual[i] = xi_;
451 molefracdual[j] = xj_;
452 molefracdual[k] = xk_;
453 return eval(AlphaCaller(w, T_, rho_, molefracdual)); };
454 auto wrts = std::tuple_cat(build_duplicated_tuple<iT>(std::ref(Trecipad)), build_duplicated_tuple<iD>(std::ref(rhoad)), build_duplicated_tuple<iXi>(std::ref(xi)), build_duplicated_tuple<iXj>(std::ref(xj)), build_duplicated_tuple<iXk>(std::ref(xk)));
455 auto der = derivatives(f, std::apply(wrt_helper(), wrts), at(Trecipad, rhoad, xi, xj, xk));
456 return powi(forceeval(1.0 / T), iT) * powi(rho, iD) * der[der.size() - 1];
457 }
458
467 template<int iT, int iD, ADBackends be = ADBackends::autodiff>
468 static auto get_Arxy(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
469 return get_Agenxy<iT, iD, be>(model, T, rho, molefrac);
470 }
471
480 template<int iT, int iD, ADBackends be = ADBackends::autodiff>
481 static auto get_Aigxy(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
482 return get_Agenxy<iT, iD, be>(model, T, rho, molefrac);
483 }
484
485 template<ADBackends be = ADBackends::autodiff>
486 static auto get_Ar10(const Model& model, const Scalar &T, const Scalar &rho, const VectorType& molefrac) {
487 return get_Arxy<1, 0, be>(model, T, rho, molefrac);
488 }
489
490 template<ADBackends be = ADBackends::autodiff>
491 static Scalar get_Ar01(const Model& model, const Scalar&T, const Scalar &rho, const VectorType& molefrac){
492 return get_Arxy<0, 1, be>(model, T, rho, molefrac);
493 }
494
495 template<ADBackends be = ADBackends::autodiff>
496 static auto get_Ar02(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
497 return get_Arxy<0, 2, be>(model, T, rho, molefrac);
498 }
499
500 template<ADBackends be = ADBackends::autodiff>
501 static auto get_Ar03(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
502 return get_Arxy<0, 3, be>(model, T, rho, molefrac);
503 }
504
505 template<ADBackends be = ADBackends::autodiff>
506 static auto get_Ar20(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
507 return get_Arxy<2, 0, be>(model, T, rho, molefrac);
508 }
509
510 template<ADBackends be = ADBackends::autodiff>
511 static auto get_Ar30(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
512 return get_Arxy<3, 0, be>(model, T, rho, molefrac);
513 }
514
515 template<ADBackends be = ADBackends::autodiff>
516 static auto get_Ar21(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
517 return get_Arxy<2, 1, be>(model, T, rho, molefrac);
518 }
519
520 template<ADBackends be = ADBackends::autodiff>
521 static auto get_Ar12(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
522 return get_Arxy<1, 2, be>(model, T, rho, molefrac);
523 }
524
525 template<ADBackends be = ADBackends::autodiff>
526 static auto get_Ar11(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
527 return get_Arxy<1, 1, be>(model, T, rho, molefrac);
528 }
529
530 template<ADBackends be = ADBackends::autodiff>
531 static auto get_Aig11(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
532 return get_Aigxy<1, 1, be>(model, T, rho, molefrac);
533 }
534
535 template<int Nderiv, ADBackends be = ADBackends::autodiff, class AlphaWrapper>
536 static auto get_Agen0n(const AlphaWrapper& w, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
537 std::valarray<Scalar> o(Nderiv+1);
538 if constexpr (be == ADBackends::autodiff) {
539 // If a pure derivative, then we can use autodiff::Real for that variable and Scalar for other variable
540 autodiff::Real<Nderiv, Scalar> rho_ = rho;
541 auto f = [&w, &T, &molefrac](const auto& rho__) { return AlphaCaller(w, T, rho__, molefrac); };
542 auto ders = derivatives(f, along(1), at(rho_));
543 for (auto n = 0; n <= Nderiv; ++n) {
544 o[n] = forceeval(powi(rho, n) * ders[n]);
545 }
546 return o;
547 }
548#if defined(TEQP_MULTICOMPLEX_ENABLED)
549 else {
550 using fcn_t = std::function<mcx::MultiComplex<Scalar>(const mcx::MultiComplex<Scalar>&)>;
551 bool and_val = true;
552 fcn_t f = [&w, &T, &molefrac](const auto& rhomcx) { return AlphaCaller(w, T, rhomcx, molefrac); };
553 auto ders = diff_mcx1(f, rho, Nderiv, and_val);
554 for (auto n = 0; n <= Nderiv; ++n) {
555 o[n] = powi(rho, n) * ders[n];
556 }
557 return o;
558 }
559#endif
560// throw std::invalid_argument("algorithmic differentiation backend is invalid in get_Ar0n");
561 }
562
563 template<int Nderiv, ADBackends be = ADBackends::autodiff, class AlphaWrapper>
564 static auto get_Agenn0(const AlphaWrapper& w, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
565 std::valarray<Scalar> o(Nderiv+1);
566 Scalar Trecip = 1.0 / T;
567 if constexpr (be == ADBackends::autodiff) {
568 // If a pure derivative, then we can use autodiff::Real for that variable and Scalar for other variable
569 autodiff::Real<Nderiv, Scalar> Trecipad = Trecip;
570 auto f = [&w, &rho, &molefrac](const auto& Trecip__) {return AlphaCaller(w, forceeval(1.0/Trecip__), rho, molefrac); };
571 auto ders = derivatives(f, along(1), at(Trecipad));
572 for (auto n = 0; n <= Nderiv; ++n) {
573 o[n] = powi(Trecip, n) * ders[n];
574 }
575 }
576#if defined(TEQP_MULTICOMPLEX_ENABLED)
577 else if constexpr (be == ADBackends::multicomplex) {
578 using fcn_t = std::function<mcx::MultiComplex<Scalar>(const mcx::MultiComplex<Scalar>&)>;
579 fcn_t f = [&](const auto& Trecipmcx) { return AlphaCaller(w, 1.0/Trecipmcx, rho, molefrac); };
580 auto ders = diff_mcx1(f, Trecip, Nderiv+1, true /* and_val */);
581 for (auto n = 0; n <= Nderiv; ++n) {
582 o[n] = powi(Trecip, n) * ders[n];
583 }
584 }
585#endif
586 else {
587 throw std::invalid_argument("algorithmic differentiation backend is invalid in get_Agenn0");
588 }
589 return o;
590 }
591
600 template<int iT, ADBackends be = ADBackends::autodiff>
601 static auto get_Arn0(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
602 return get_Agenn0<iT, be>(model, T, rho, molefrac);
603 }
604
613 template<int iD, ADBackends be = ADBackends::autodiff>
614 static auto get_Ar0n(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
615 return get_Agen0n<iD, be>(model, T, rho, molefrac);
616 }
617
618
619 template<ADBackends be = ADBackends::autodiff>
620 static auto get_Ar(const int itau, const int idelta, const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
621 if (itau == 0) {
622 if (idelta == 0) {
623 return get_Ar00(model, T, rho, molefrac);
624 }
625 else if (idelta == 1) {
626 return get_Ar01(model, T, rho, molefrac);
627 }
628 else if (idelta == 2) {
629 return get_Ar02(model, T, rho, molefrac);
630 }
631 else if (idelta == 3) {
632 return get_Ar03(model, T, rho, molefrac);
633 }
634 else {
635 throw std::invalid_argument("Invalid value for idelta");
636 }
637 }
638 else if (itau == 1){
639 if (idelta == 0) {
640 return get_Ar10(model, T, rho, molefrac);
641 }
642 else if (idelta == 1) {
643 return get_Ar11(model, T, rho, molefrac);
644 }
645 else if (idelta == 2) {
646 return get_Ar12(model, T, rho, molefrac);
647 }
648 else {
649 throw std::invalid_argument("Invalid value for idelta");
650 }
651 }
652 else if (itau == 2) {
653 if (idelta == 0) {
654 return get_Ar20(model, T, rho, molefrac);
655 }
656 else if (idelta == 1) {
657 return get_Ar21(model, T, rho, molefrac);
658 }
659 else {
660 throw std::invalid_argument("Invalid value for idelta");
661 }
662 }
663 else if (itau == 3) {
664 if (idelta == 0) {
665 return get_Ar30(model, T, rho, molefrac);
666 }
667 else {
668 throw std::invalid_argument("Invalid value for idelta");
669 }
670 }
671 else {
672 throw std::invalid_argument("Invalid value for itau");
673 }
674 }
675
676 template<ADBackends be = ADBackends::autodiff>
677 static auto get_neff(const Model& model, const Scalar& T, const Scalar& rho, const VectorType& molefrac) {
678 auto Ar01 = get_Ar01<be>(model, T, rho, molefrac);
679 auto Ar11 = get_Ar11<be>(model, T, rho, molefrac);
680 auto Ar20 = get_Ar20<be>(model, T, rho, molefrac);
681 return -3*(Ar01-Ar11)/Ar20;
682 }
683};
684
685template<typename Model, typename Scalar = double, typename VectorType = Eigen::ArrayXd>
687
688 static auto get_B2vir(const Model& model, const Scalar &T, const VectorType& molefrac) {
689 Scalar h = 1e-100;
690 // B_2 = dalphar/drho|T,z at rho=0
691 auto B2 = model.alphar(T, std::complex<Scalar>(0.0, h), molefrac).imag()/h;
692 return B2;
693 }
694
705 template <int Nderiv, ADBackends be = ADBackends::autodiff>
706 static auto get_Bnvir(const Model& model, const Scalar &T, const VectorType& molefrac)
707 {
708 std::map<int, double> dnalphardrhon;
709 if constexpr(be == ADBackends::autodiff){
710 auto f = [&model, &T, &molefrac](const auto& rho_) { return model.alphar(T, rho_, molefrac); };
711 autodiff::Real<Nderiv, Scalar> rhoreal = 0.0;
712 auto derivs = derivatives(f, along(1), at(rhoreal));
713
714 for (auto n = 1; n < Nderiv; ++n){
715 dnalphardrhon[n] = derivs[n];
716 }
717 }
718#if defined(TEQP_MULTICOMPLEX_ENABLED)
719 else if constexpr(be == ADBackends::multicomplex){
720 using namespace mcx;
721 using fcn_t = std::function<MultiComplex<double>(const MultiComplex<double>&)>;
722 fcn_t f = [&model, &T, &molefrac](const auto& rho_) { return model.alphar(T, rho_, molefrac); };
723 auto derivs = diff_mcx1(f, 0.0, Nderiv, true /* and_val */);
724 for (auto n = 1; n < Nderiv; ++n){
725 dnalphardrhon[n] = derivs[n];
726 }
727 }
728#endif
729 else{
730 //static_assert(false, "algorithmic differentiation backend is invalid");
731 throw std::invalid_argument("algorithmic differentiation backend is invalid in get_Bnvir");
732 }
733
734 std::map<int, Scalar> o;
735 for (int n = 2; n <= Nderiv; ++n) {
736 o[n] = dnalphardrhon[n-1];
737 // 0! = 1, 1! = 1, so only n>3 terms need factorial correction
738 if (n > 3) {
739 auto factorial = [](int N) {return tgamma(N + 1); };
740 o[n] /= factorial(n-2);
741 }
742 }
743 return o;
744 }
745
748 template <ADBackends be = ADBackends::autodiff>
749 static auto get_Bnvir_runtime(const int Nderiv, const Model& model, const Scalar &T, const VectorType& molefrac) {
750 switch(Nderiv){
751 case 2: return get_Bnvir<2,be>(model, T, molefrac);
752 case 3: return get_Bnvir<3,be>(model, T, molefrac);
753 case 4: return get_Bnvir<4,be>(model, T, molefrac);
754 case 5: return get_Bnvir<5,be>(model, T, molefrac);
755 case 6: return get_Bnvir<6,be>(model, T, molefrac);
756 case 7: return get_Bnvir<7,be>(model, T, molefrac);
757 case 8: return get_Bnvir<8,be>(model, T, molefrac);
758 case 9: return get_Bnvir<9,be>(model, T, molefrac);
759 default: throw std::invalid_argument("Only Nderiv up to 9 is supported, get_Bnvir templated function allows more");
760 }
761 }
762
775 template <int Nderiv, int NTderiv, ADBackends be = ADBackends::autodiff>
776 static auto get_dmBnvirdTm(const Model& model, const Scalar& T, const VectorType& molefrac)
777 {
778 std::map<int, Scalar> o;
779 auto factorial = [](int N) {return tgamma(N + 1); };
780 if constexpr (be == ADBackends::autodiff) {
781 autodiff::HigherOrderDual<NTderiv + Nderiv-1, double> rhodual = 0.0, Tdual = T;
782 auto f = [&model, &molefrac](const auto& T_, const auto& rho_) { return model.alphar(T_, rho_, molefrac); };
783 auto wrts = std::tuple_cat(build_duplicated_tuple<NTderiv>(std::ref(Tdual)), build_duplicated_tuple<Nderiv-1>(std::ref(rhodual)));
784 auto derivs = derivatives(f, std::apply(wrt_helper(), wrts), at(Tdual, rhodual));
785 return derivs.back() / factorial(Nderiv - 2);
786 }
787#if defined(TEQP_MULTICOMPLEX_ENABLED)
788 else if constexpr (be == ADBackends::multicomplex) {
789 using namespace mcx;
790 using fcn_t = std::function<MultiComplex<double>(const std::valarray<MultiComplex<double>>&)>;
791 fcn_t f = [&model, &molefrac](const auto& zs) {
792 auto T_ = zs[0], rho_ = zs[1];
793 return model.alphar(T_, rho_, molefrac);
794 };
795 std::valarray<double> at = { T, 0.0 };
796 auto deriv = diff_mcxN(f, at, { NTderiv, Nderiv-1});
797 return deriv / factorial(Nderiv - 2);
798 }
799#endif
800 else {
801 //static_assert(false, "algorithmic differentiation backend is invalid");
802 throw std::invalid_argument("algorithmic differentiation backend is invalid in get_Bnvir");
803 }
804 }
805
808 template <ADBackends be = ADBackends::autodiff>
809 static auto get_dmBnvirdTm_runtime(const int Nderiv, const int NTderiv, const Model& model, const Scalar& T, const VectorType& molefrac) {
810 if (Nderiv == 2) { // B_2
811 switch (NTderiv) {
812 case 0: return get_Bnvir<2, be>(model, T, molefrac)[2];
813 case 1: return get_dmBnvirdTm<2, 1, be>(model, T, molefrac);
814 case 2: return get_dmBnvirdTm<2, 2, be>(model, T, molefrac);
815 case 3: return get_dmBnvirdTm<2, 3, be>(model, T, molefrac);
816 default: throw std::invalid_argument("NTderiv is invalid in get_dmBnvirdTm_runtime");
817 }
818 }
819 else if (Nderiv == 3) { // B_3
820 switch (NTderiv) {
821 case 0: return get_Bnvir<3, be>(model, T, molefrac)[3];
822 case 1: return get_dmBnvirdTm<3, 1, be>(model, T, molefrac);
823 case 2: return get_dmBnvirdTm<3, 2, be>(model, T, molefrac);
824 case 3: return get_dmBnvirdTm<3, 3, be>(model, T, molefrac);
825 default: throw std::invalid_argument("NTderiv is invalid in get_dmBnvirdTm_runtime");
826 }
827 }
828 else if (Nderiv == 4) { // B_4
829 switch (NTderiv) {
830 case 0: return get_Bnvir<4, be>(model, T, molefrac)[4];
831 case 1: return get_dmBnvirdTm<4, 1, be>(model, T, molefrac);
832 case 2: return get_dmBnvirdTm<4, 2, be>(model, T, molefrac);
833 case 3: return get_dmBnvirdTm<4, 3, be>(model, T, molefrac);
834 default: throw std::invalid_argument("NTderiv is invalid in get_dmBnvirdTm_runtime");
835 }
836 }
837 else {
838 throw std::invalid_argument("Nderiv is invalid in get_dmBnvirdTm_runtime");
839 }
840 }
841
848 static auto get_B12vir(const Model& model, const Scalar &T, const VectorType& molefrac) {
849 if (molefrac.size() != 2) { throw std::invalid_argument("length of mole fraction vector must be 2 in get_B12vir"); }
850 auto B2 = get_B2vir(model, T, molefrac); // Overall B2 for mixture
851 const auto xpure0 = (Eigen::ArrayXd(2) << 1,0).finished();
852 const auto xpure1 = (Eigen::ArrayXd(2) << 0,1).finished();
853 auto B20 = get_B2vir(model, T, xpure0); // Pure first component with index 0
854 auto B21 = get_B2vir(model, T, xpure1); // Pure second component with index 1
855 auto z0 = molefrac[0];
856 auto B12 = (B2 - z0*z0*B20 - (1-z0)*(1-z0)*B21)/(2*z0*(1-z0));
857 return B12;
858 }
859};
860
948template<typename Model, typename Scalar = double, typename VectorType = Eigen::ArrayXd>
950
954 static auto get_splus(const Model& model, const Scalar &T, const VectorType& rhovec) {
955 auto rhotot = rhovec.sum();
956 auto molefrac = rhovec / rhotot;
957 return model.alphar(T, rhotot, molefrac) - get_Ar10(model, T, rhovec);
958 }
959
963 static auto get_pr(const Model& model, const Scalar &T, const VectorType& rhovec)
964 {
965 auto rhotot_ = rhovec.sum();
966 auto molefrac = (rhovec / rhotot_).eval();
967 auto h = 1e-100;
968 auto Ar01 = model.alphar(T, std::complex<double>(rhotot_, h), molefrac).imag() / h * rhotot_;
969 return Ar01*rhotot_*model.R(molefrac)*T;
970 }
971
972 static auto get_Ar00(const Model& model, const Scalar& T, const VectorType& rhovec) {
973 auto rhotot = rhovec.sum();
974 auto molefrac = (rhovec / rhotot).eval();
975 return model.alphar(T, rhotot, molefrac);
976 }
977
978 static auto get_Ar10(const Model& model, const Scalar& T, const VectorType& rhovec) {
979 auto rhotot = rhovec.sum();
980 auto molefrac = (rhovec / rhotot).eval();
981 return -T * derivT([&model, &rhotot, &molefrac](const auto& T, const auto& /*rhovec*/) { return model.alphar(T, rhotot, molefrac); }, T, rhovec);
982 }
983
984 static auto get_Ar01(const Model& model, const Scalar &T, const VectorType& rhovec) {
985 auto rhotot_ = std::accumulate(std::begin(rhovec), std::end(rhovec), (decltype(rhovec[0]))0.0);
986 decltype(rhovec[0] * T) Ar01 = 0.0;
987 for (auto i = 0; i < rhovec.size(); ++i) {
988 auto Ar00 = [&model](const auto &T, const auto&rhovec) {
989 auto rhotot = rhovec.sum();
990 auto molefrac = rhovec / rhotot;
991 return model.alphar(T, rhotot, molefrac);
992 };
993 Ar01 += rhovec[i] * derivrhoi(Ar00, T, rhovec, i);
994 }
995 return Ar01;
996 }
997
1001 static auto get_Psir(const Model& model, const Scalar &T, const VectorType& rhovec) {
1002 auto rhotot_ = rhovec.sum();
1003 auto molefrac = rhovec / rhotot_;
1004 return model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_;
1005 }
1006
1010 static auto get_dPsirdT_constrhovec(const Model& model, const Scalar& T, const VectorType& rhovec) {
1011 auto rhotot_ = rhovec.sum();
1012 auto molefrac = (rhovec / rhotot_).eval();
1013 autodiff::Real<1, Scalar> Tad = T;
1014 auto f = [&model, &rhotot_, &molefrac](const auto& T_) {return rhotot_*model.R(molefrac)*T_*model.alphar(T_, rhotot_, molefrac); };
1015 return derivatives(f, along(1), at(Tad))[1];
1016 }
1017
1023 static auto build_Psir_Hessian_autodiff(const Model& model, const Scalar& T, const VectorType& rho) {
1024 // Double derivatives in each component's concentration
1025 // N^N matrix (symmetric)
1026
1027 dual2nd u; // the output scalar u = f(x), evaluated together with Hessian below
1028 ArrayXdual2nd g;
1029 ArrayXdual2nd rhovecc(rho.size()); for (auto i = 0; i < rho.size(); ++i) { rhovecc[i] = rho[i]; }
1030 auto hfunc = [&model, &T](const ArrayXdual2nd& rho_) {
1031 auto rhotot_ = rho_.sum();
1032 auto molefrac = (rho_ / rhotot_).eval();
1033 return eval(model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_);
1034 };
1035 return autodiff::hessian(hfunc, wrt(rhovecc), at(rhovecc), u, g).eval(); // evaluate the function value u, its gradient, and its Hessian matrix H
1036 }
1037
1043 static auto build_Psir_fgradHessian_autodiff(const Model& model, const Scalar& T, const VectorType& rho) {
1044 // Double derivatives in each component's concentration
1045 // N^N matrix (symmetric)
1046
1047 dual2nd u; // the output scalar u = f(x), evaluated together with Hessian below
1048 ArrayXdual g;
1049 ArrayXdual2nd rhovecc(rho.size()); for (auto i = 0; i < rho.size(); ++i) { rhovecc[i] = rho[i]; }
1050 auto hfunc = [&model, &T](const ArrayXdual2nd& rho_) {
1051 auto rhotot_ = rho_.sum();
1052 auto molefrac = (rho_ / rhotot_).eval();
1053 return eval(model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_);
1054 };
1055 // Evaluate the function value u, its gradient, and its Hessian matrix H
1056 Eigen::MatrixXd H = autodiff::hessian(hfunc, wrt(rhovecc), at(rhovecc), u, g);
1057 // Remove autodiff stuff from the numerical values
1058 auto f = getbaseval(u);
1059 auto gg = g.cast<double>().eval();
1060 return std::make_tuple(f, gg, H);
1061 }
1062
1068 static auto build_Psi_Hessian_autodiff(const Model& model, const Scalar& T, const VectorType& rho) {
1069 auto rhotot_ = rho.sum();
1070 auto molefrac = (rho / rhotot_).eval();
1071 auto H = build_Psir_Hessian_autodiff(model, T, rho).eval();
1072 for (auto i = 0; i < 2; ++i) {
1073 H(i, i) += model.R(molefrac) * T / rho[i];
1074 }
1075 return H;
1076 }
1077
1078#if defined(TEQP_MULTICOMPLEX_ENABLED)
1084 static auto build_Psir_Hessian_mcx(const Model& model, const Scalar& T, const VectorType& rho) {
1085 // Double derivatives in each component's concentration
1086 // N^N matrix (symmetric)
1087 using namespace mcx;
1088
1089 // Lambda function for getting Psir with multicomplex concentrations
1090 using fcn_t = std::function< MultiComplex<double>(const Eigen::ArrayX<MultiComplex<double>>&)>;
1091 fcn_t func = [&model, &T](const auto& rhovec) {
1092 auto rhotot_ = rhovec.sum();
1093 auto molefrac = (rhovec / rhotot_).eval();
1094 return model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_;
1095 };
1096 using mattype = Eigen::ArrayXXd;
1097 auto H = get_Hessian<mattype, fcn_t, VectorType, HessianMethods::Multiple>(func, rho);
1098 return H;
1099 }
1100#endif
1101
1107 static auto build_Psir_gradient_autodiff(const Model& model, const Scalar& T, const VectorType& rho) {
1108 ArrayXdual rhovecc(rho.size()); for (auto i = 0; i < rho.size(); ++i) { rhovecc[i] = rho[i]; }
1109 auto psirfunc = [&model, &T](const ArrayXdual& rho_) {
1110 auto rhotot_ = rho_.sum();
1111 auto molefrac = (rho_ / rhotot_).eval();
1112 return eval(model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_);
1113 };
1114 auto val = autodiff::gradient(psirfunc, wrt(rhovecc), at(rhovecc)).eval(); // evaluate the gradient
1115 return val;
1116 }
1117
1118#if defined(TEQP_MULTICOMPLEX_ENABLED)
1124 static auto build_Psir_gradient_multicomplex(const Model& model, const Scalar& T, const VectorType& rho) {
1125 using rho_type = typename VectorType::value_type;
1126 Eigen::ArrayX<mcx::MultiComplex<rho_type>> rhovecc(rho.size()); //for (auto i = 0; i < rho.size(); ++i) { rhovecc[i] = rho[i]; }
1127 auto psirfunc = [&model](const auto &T, const auto& rho_) {
1128 auto rhotot_ = rho_.sum();
1129 auto molefrac = (rho_ / rhotot_).eval();
1130 return eval(model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_);
1131 };
1132 VectorType out(rho.size());
1133 for (int i = 0; i < rho.size(); ++i) {
1134 out[i] = derivrhoi(psirfunc, T, rho, i);
1135 }
1136 return out;
1137 }
1138#endif
1144 static auto build_Psir_gradient_complex_step(const Model& model, const Scalar& T, const VectorType& rho) {
1145 using rho_type = typename VectorType::value_type;
1146 Eigen::ArrayX<std::complex<rho_type>> rhovecc(rho.size()); for (auto i = 0; i < rho.size(); ++i) { rhovecc[i] = rho[i]; }
1147 auto psirfunc = [&model](const auto& T, const auto& rho_) {
1148 auto rhotot_ = rho_.sum();
1149 auto molefrac = (rho_ / rhotot_).eval();
1150 return forceeval(model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_);
1151 };
1152 VectorType out(rho.size());
1153 for (int i = 0; i < rho.size(); ++i) {
1154 auto rhocopy = rhovecc;
1155 rho_type h = 1e-100;
1156 rhocopy[i] = rhocopy[i] + std::complex<rho_type>(0,h);
1157 auto calc = psirfunc(T, rhocopy);
1158 out[i] = calc.imag / static_cast<double>(h);
1159 }
1160 return out;
1161 }
1162
1163 /* Convenience function to select the correct implementation at compile-time */
1164 template<ADBackends be = ADBackends::autodiff>
1165 static auto build_Psir_gradient(const Model& model, const Scalar& T, const VectorType& rho) {
1166 if constexpr (be == ADBackends::autodiff) {
1167 return build_Psir_gradient_autodiff(model, T, rho);
1168 }
1169#if defined(TEQP_MULTICOMPLEX_ENABLED)
1170 else if constexpr (be == ADBackends::multicomplex) {
1171 return build_Psir_gradient_multicomplex(model, T, rho);
1172 }
1173#endif
1174#if defined(TEQP_COMPLEXSTEP_ENABLED)
1175 else if constexpr (be == ADBackends::complex_step) {
1176 return build_Psir_gradient_complex_step(model, T, rho);
1177 }
1178#endif
1179 }
1180
1188 static auto get_chempotVLE_autodiff(const Model& model, const Scalar& T, const VectorType& rho) {
1189 typename VectorType::value_type rhotot = rho.sum();
1190 auto molefrac = (rho / rhotot).eval();
1191 auto rhorefideal = 1.0;
1192 return (build_Psir_gradient_autodiff(model, T, rho).array() + model.R(molefrac)*T*(rhorefideal + log(rho / rhorefideal))).eval();
1193 }
1194
1200 template<ADBackends be = ADBackends::autodiff>
1201 static auto get_fugacity_coefficients(const Model& model, const Scalar& T, const VectorType& rhovec) {
1202 VectorType lnphi = get_ln_fugacity_coefficients<be>(model, T, rhovec);
1203 return exp(lnphi).eval();
1204 }
1205
1211 template<ADBackends be = ADBackends::autodiff>
1212 static auto get_ln_fugacity_coefficients(const Model& model, const Scalar& T, const VectorType& rhovec) {
1213 auto rhotot = forceeval(rhovec.sum());
1214 auto molefrac = (rhovec / rhotot).eval();
1215 auto R = model.R(molefrac);
1217 auto Z = 1.0 + tdx::template get_Ar01<be>(model, T, rhotot, molefrac);
1218 auto grad = build_Psir_gradient<be>(model, T, rhovec).eval();
1219 auto RT = R * T;
1220 auto lnphi = ((grad / RT).array() - log(Z)).eval();
1221 return forceeval(lnphi.eval());
1222 }
1223
1224 template<ADBackends be = ADBackends::autodiff>
1225 static auto get_ln_fugacity_coefficients_Trhomolefracs(const Model& model, const Scalar& T, const Scalar& rhotot, const VectorType& molefrac) {
1226 auto R = model.R(molefrac);
1228 auto Z = 1.0 + tdx::template get_Ar01<be>(model, T, rhotot, molefrac);
1229 auto rhovec = (rhotot*molefrac).eval();
1230 auto grad = build_Psir_gradient<be>(model, T, rhovec).eval();
1231 auto RT = R * T;
1232 auto lnphi = ((grad / RT).array() - log(Z)).eval();
1233 return forceeval(lnphi.eval());
1234 }
1235
1236 template<ADBackends be = ADBackends::autodiff>
1237 static auto get_ln_fugacity_coefficients1(const Model& model, const Scalar& T, const VectorType& rhovec) {
1238 auto rhotot = forceeval(rhovec.sum());
1239 auto molefrac = (rhovec / rhotot).eval();
1240 auto R = model.R(molefrac);
1241 auto grad = build_Psir_gradient<be>(model, T, rhovec).eval();
1242 auto RT = R * T;
1243 auto lnphi = ((grad / RT).array()).eval();
1244 return forceeval(lnphi);
1245 }
1246
1247 template<ADBackends be = ADBackends::autodiff>
1248 static auto get_ln_fugacity_coefficients2(const Model& model, const Scalar& T, const VectorType& rhovec) {
1249 auto rhotot = forceeval(rhovec.sum());
1250 auto molefrac = (rhovec / rhotot).eval();
1252 auto Z = 1.0 + tdx::template get_Ar01<be>(model, T, rhotot, molefrac);
1253 return forceeval(-log(Z));
1254 }
1255
1256 static Eigen::ArrayXd build_d2PsirdTdrhoi_autodiff(const Model& model, const Scalar& T, const VectorType& rho) {
1257 Eigen::ArrayXd deriv(rho.size());
1258 // d^2psir/dTdrho_i
1259 for (auto i = 0; i < rho.size(); ++i) {
1260 auto psirfunc = [&model, &rho, i](const auto& T, const auto& rhoi) {
1261 ArrayXdual2nd rhovecc(rho.size()); for (auto j = 0; j < rho.size(); ++j) { rhovecc[j] = rho[j]; }
1262 rhovecc[i] = rhoi;
1263 auto rhotot_ = rhovecc.sum();
1264 auto molefrac = (rhovecc / rhotot_).eval();
1265 return eval(model.alphar(T, rhotot_, molefrac) * model.R(molefrac) * T * rhotot_);
1266 };
1267 dual2nd Tdual = T, rhoidual = rho[i];
1268 auto [u00, u10, u11] = derivatives(psirfunc, wrt(Tdual, rhoidual), at(Tdual, rhoidual));
1269 deriv[i] = u11;
1270 }
1271 return deriv;
1272 }
1273
1274 static Eigen::ArrayXd build_d2alphardrhodxi_constT(const Model& model, const Scalar& T, const Scalar& rhomolar, const VectorType& molefrac) {
1275 Eigen::ArrayXd deriv(molefrac.size());
1276 // d^2alphar/drhodx_i|T
1277 for (auto i = 0; i < molefrac.size(); ++i) {
1278 auto alpharfunc = [&model, &T, &molefrac, i](const auto& rho, const auto& xi) {
1279 ArrayXdual2nd molefracdual = molefrac.template cast<autodiff::dual2nd>();
1280 molefracdual[i] = xi;
1281 return forceeval(model.alphar(T, rho, molefracdual));
1282 };
1283 dual2nd rhodual = rhomolar, xidual = molefrac[i];
1284 auto [u00, u10, u11] = derivatives(alpharfunc, wrt(rhodual, xidual), at(rhodual, xidual));
1285 deriv[i] = u11;
1286 }
1287 return deriv;
1288 }
1289
1298 template<ADBackends be = ADBackends::autodiff>
1299 static auto get_d_ln_fugacity_coefficients_dT_constrhovec(const Model& model, const Scalar& T, const VectorType& rhovec) {
1300 auto rhotot = forceeval(rhovec.sum());
1301 auto molefrac = (rhovec / rhotot).eval();
1302 auto R = model.R(molefrac);
1304 auto Z = 1.0 + tdx::template get_Ar01<be>(model, T, rhotot, molefrac);
1305 auto dZdT_Z = tdx::template get_Ar11<be>(model, T, rhotot, molefrac)/(-T)/Z; // Note: (1/T)dX/d(1/T) = -TdX/dT, the deriv in RHS is what we want, the left is what we get, so divide by -T
1306 VectorType grad = build_Psir_gradient<be>(model, T, rhovec).eval();
1307 VectorType Tgrad = build_d2PsirdTdrhoi_autodiff(model, T, rhovec);
1308 return forceeval((1/(R*T)*(Tgrad - 1.0/T*grad)-dZdT_Z).eval());
1309 }
1310
1316 template<ADBackends be = ADBackends::autodiff>
1317 static auto get_lnZ_Z_dZdrho(const Model& model, const Scalar& T, const VectorType& rhovec) {
1318 auto rhotot = forceeval(rhovec.sum());
1319 auto molefrac = (rhovec / rhotot).eval();
1321 auto Ar01 = tdx::template get_Ar01<be>(model, T, rhotot, molefrac);
1322 auto Ar02 = tdx::template get_Ar02<be>(model, T, rhotot, molefrac);
1323 auto Z = 1.0 + Ar01;
1324 auto dZdrho = (Ar01 + Ar02)/rhotot; // (dZ/rhotot)_{T,x}
1325 return std::make_tuple(log(Z), Z, dZdrho);
1326 }
1327
1335 template<ADBackends be = ADBackends::autodiff>
1336 static auto get_d_ln_fugacity_coefficients_drho_constTmolefracs(const Model& model, const Scalar& T, const VectorType& rhovec) {
1337 auto rhotot = forceeval(rhovec.sum());
1338 auto molefrac = (rhovec / rhotot).eval();
1339 auto R = model.R(molefrac);
1340 auto [lnZ, Z, dZdrho] = get_lnZ_Z_dZdrho(model, T, rhovec);
1341 auto hessian = build_Psir_Hessian_autodiff(model, T, rhovec);
1342 return forceeval((1/(R*T)*(hessian*molefrac.matrix()).array() - dZdrho/Z).eval());
1343 }
1344
1352 template<ADBackends be = ADBackends::autodiff>
1353 static auto get_d_ln_fugacity_coefficients_dv_constTmolefracs(const Model& model, const Scalar& T, const VectorType& rhovec) {
1354 auto rhotot = forceeval(rhovec.sum());
1355 auto drhodv = -rhotot*rhotot; // rho = 1/v; drho/dv = -1/v^2 = -rho^2
1356 return get_d_ln_fugacity_coefficients_drho_constTmolefracs(model, T, rhovec)*drhodv;
1357 }
1358
1366 template<ADBackends be = ADBackends::autodiff>
1367 static auto get_d_ln_fugacity_coefficients_dmolefracs_constTrho(const Model& model, const Scalar& T, const VectorType& rhovec) {
1368 auto rhotot = forceeval(rhovec.sum());
1369 auto molefrac = (rhovec / rhotot).eval();
1370 auto R = model.R(molefrac);
1371
1373 auto Ar01 = tdx::template get_Ar01<be>(model, T, rhotot, molefrac);
1374 auto Z = 1.0 + Ar01;
1375 VectorType dZdx = rhotot*build_d2alphardrhodxi_constT(model, T, rhotot, molefrac);
1376 Eigen::RowVector<decltype(rhotot), Eigen::Dynamic> dZdx_Z = dZdx/Z;
1377
1378 // Starting matrix is from the first term
1379 auto hessian = build_Psir_Hessian_autodiff(model, T, rhovec);
1380 Eigen::ArrayXXd out = rhotot/(R*T)*hessian;
1381
1382 // Then each row gets the second part
1383 out.rowwise() -= dZdx_Z.array();
1384 return out;
1385 }
1386
1387 template<ADBackends be = ADBackends::autodiff>
1388 static auto get_d_ln_fugacity_coefficients_dmolefracs_constTrho1(const Model& model, const Scalar& T, const VectorType& rhovec) {
1389 auto rhotot = forceeval(rhovec.sum());
1390 auto molefrac = (rhovec / rhotot).eval();
1391 auto R = model.R(molefrac);
1392
1393 auto hessian = build_Psir_Hessian_autodiff(model, T, rhovec);
1394 // Starting matrix is from the first term
1395 Eigen::ArrayXXd out = 1/(R*T)*rhotot*hessian;
1396 return out;
1397 }
1398
1403 static auto get_dchempotdT_autodiff(const Model& model, const Scalar& T, const VectorType& rhovec) {
1404 auto rhotot = rhovec.sum();
1405 auto molefrac = (rhovec / rhotot).eval();
1406 auto rhorefideal = 1.0;
1407 return (build_d2PsirdTdrhoi_autodiff(model, T, rhovec) + model.R(molefrac)*(rhorefideal + log(rhovec/rhorefideal))).eval();
1408 }
1409
1413 static auto get_dpdT_constrhovec(const Model& model, const Scalar& T, const VectorType& rhovec) {
1414 auto rhotot = rhovec.sum();
1415 auto molefrac = (rhovec / rhotot).eval();
1416 auto dPsirdT = get_dPsirdT_constrhovec(model, T, rhovec);
1417 return rhotot*model.R(molefrac) - dPsirdT + rhovec.matrix().dot(build_d2PsirdTdrhoi_autodiff(model, T, rhovec).matrix());
1418 }
1419
1423 static auto get_dpdrhovec_constT(const Model& model, const Scalar& T, const VectorType& rhovec) {
1424 auto rhotot = rhovec.sum();
1425 auto molefrac = (rhovec / rhotot).eval();
1426 auto RT = model.R(molefrac)*T;
1427 auto [func, grad, hessian] = build_Psir_fgradHessian_autodiff(model, T, rhovec); // The hessian matrix
1428 return (RT + (hessian*rhovec.matrix()).array()).eval(); // at constant temperature
1429 }
1430
1474 static auto get_partial_molar_volumes(const Model& model, const Scalar& T, const VectorType& rhovec) {
1475 auto rhotot = rhovec.sum();
1476 auto molefrac = (rhovec / rhotot).eval();
1478 auto RT = model.R(molefrac)*T;
1479
1480 auto dpdrhovec = get_dpdrhovec_constT(model, T, rhovec);
1481 auto numerator = -rhotot*dpdrhovec;
1482
1483 auto ders = tdx::template get_Ar0n<2>(model, T, rhotot, molefrac);
1484 auto denominator = -pow2(rhotot)*RT*(1 + 2*ders[1] + ders[2]);
1485 return (numerator/denominator).eval();
1486 }
1487
1488 static VectorType get_Psir_sigma_derivs(const Model& model, const Scalar& T, const VectorType& rhovec, const VectorType& v) {
1489 autodiff::Real<4, double> sigma = 0.0;
1490 auto rhovecad = rhovec.template cast<decltype(sigma)>(), vad = v.template cast<decltype(sigma)>();
1491 auto wrapper = [&rhovecad, &vad, &T, &model](const auto& sigma_1) {
1492 auto rhovecused = (rhovecad + sigma_1 * vad).eval();
1493 auto rhotot = rhovecused.sum();
1494 auto molefrac = (rhovecused / rhotot).eval();
1495 return forceeval(model.alphar(T, rhotot, molefrac) * model.R(molefrac) * T * rhotot);
1496 };
1497 auto der = derivatives(wrapper, along(1), at(sigma));
1498 VectorType ret(der.size());
1499 for (auto i = 0; i < ret.size(); ++i){ ret[i] = der[i];}
1500 return ret;
1501 }
1502};
1503
1504template<int Nderivsmax>
1506
1507public:
1508 Eigen::Array<double, Nderivsmax+1, Nderivsmax+1> derivs;
1509
1510 template<typename Model, typename Scalar, typename VecType>
1511 DerivativeHolderSquare(const Model& model, const Scalar& T, const Scalar& rho, const VecType& z) {
1512 using tdx = TDXDerivatives<decltype(model), Scalar, VecType>;
1513 static_assert(Nderivsmax == 2, "It's gotta be 2 for now");
1514
1515
1516 auto AX02 = tdx::template get_Agen0n<2>(model, T, rho, z);
1517 derivs(0, 0) = AX02[0];
1518 derivs(0, 1) = AX02[1];
1519 derivs(0, 2) = AX02[2];
1520
1521 auto AX20 = tdx::template get_Agenn0<2>(model, T, rho, z);
1522 derivs(0, 0) = AX20[0];
1523 derivs(1, 0) = AX20[1];
1524 derivs(2, 0) = AX20[2];
1525
1526 derivs(1, 1) = tdx::template get_Agenxy<1,1>(model, T, rho, z);
1527 }
1528};
1529
1530
1531
1532}; // namespace teqp
DerivativeHolderSquare(const Model &model, const Scalar &T, const Scalar &rho, const VecType &z)
Definition derivs.hpp:1511
Eigen::Array< double, Nderivsmax+1, Nderivsmax+1 > derivs
Definition derivs.hpp:1508
#define get_ATrhoXiXjXk_runtime_combinations
Definition derivs.hpp:314
#define get_ATrhoXi_runtime_combinations
Definition derivs.hpp:261
#define get_ATrhoXiXj_runtime_combinations
Definition derivs.hpp:289
ContainerType::value_type derivT(const FuncType &f, TType T, const ContainerType &rho)
Given a function, use complex step derivatives to calculate the derivative with respect to the first ...
Definition derivs.hpp:33
auto getbaseval(const T &expr)
Definition types.hpp:84
T powi(const T &x, int n)
From Ulrich Deiters.
Definition types.hpp:133
T pow2(const T &x)
Definition types.hpp:7
ADBackends
Definition derivs.hpp:92
auto build_duplicated_tuple_impl(const T &val, std::index_sequence< I... >)
Definition derivs.hpp:71
auto build_duplicated_tuple(const T &val)
A function to generate a tuple of N repeated copies of argument val at compile-time.
Definition derivs.hpp:77
auto forceeval(T &&expr)
Definition types.hpp:46
ContainerType::value_type derivrhoi(const FuncType &f, TType T, const ContainerType &rho, Integer i)
Given a function, use complex step derivatives to calculate the derivative with respect to the given ...
Definition derivs.hpp:57
std::valarray< double > vad
Definition fwd.hpp:34
static auto get_dPsirdT_constrhovec(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate derivative w.r.t. T at constant molar concentrations.
Definition derivs.hpp:1010
static auto get_partial_molar_volumes(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the partial molar volumes of each component.
Definition derivs.hpp:1474
static auto build_Psir_Hessian_autodiff(const Model &model, const Scalar &T, const VectorType &rho)
Calculate the Hessian of w.r.t. the molar concentrations.
Definition derivs.hpp:1023
static auto get_d_ln_fugacity_coefficients_dv_constTmolefracs(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the derivative of the natural logarithm of fugacity coefficient of each component w....
Definition derivs.hpp:1353
static Eigen::ArrayXd build_d2alphardrhodxi_constT(const Model &model, const Scalar &T, const Scalar &rhomolar, const VectorType &molefrac)
Definition derivs.hpp:1274
static auto get_d_ln_fugacity_coefficients_dmolefracs_constTrho1(const Model &model, const Scalar &T, const VectorType &rhovec)
Definition derivs.hpp:1388
static auto get_ln_fugacity_coefficients(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the natural logarithm of fugacity coefficient of each component.
Definition derivs.hpp:1212
static auto get_splus(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the residual entropy ( ) from derivatives of alphar.
Definition derivs.hpp:954
static Eigen::ArrayXd build_d2PsirdTdrhoi_autodiff(const Model &model, const Scalar &T, const VectorType &rho)
Definition derivs.hpp:1256
static auto build_Psi_Hessian_autodiff(const Model &model, const Scalar &T, const VectorType &rho)
Calculate the Hessian of w.r.t. the molar concentrations.
Definition derivs.hpp:1068
static auto get_d_ln_fugacity_coefficients_drho_constTmolefracs(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the derivative of the natural logarithm of fugacity coefficient of each component w....
Definition derivs.hpp:1336
static auto build_Psir_fgradHessian_autodiff(const Model &model, const Scalar &T, const VectorType &rho)
Calculate the function value, gradient, and Hessian of w.r.t. the molar concentrations.
Definition derivs.hpp:1043
static auto get_Ar10(const Model &model, const Scalar &T, const VectorType &rhovec)
Definition derivs.hpp:978
static auto build_Psir_gradient_autodiff(const Model &model, const Scalar &T, const VectorType &rho)
Gradient of Psir = ar*rho w.r.t. the molar concentrations.
Definition derivs.hpp:1107
static auto get_d_ln_fugacity_coefficients_dmolefracs_constTrho(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the derivative of the natural logarithm of fugacity coefficient of each component w....
Definition derivs.hpp:1367
static auto get_pr(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the residual pressure from derivatives of alphar.
Definition derivs.hpp:963
static auto get_lnZ_Z_dZdrho(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate ln(Z), Z, and dZ/drho at constant temperature and mole fractions.
Definition derivs.hpp:1317
static auto get_Psir(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate .
Definition derivs.hpp:1001
static auto build_Psir_gradient_complex_step(const Model &model, const Scalar &T, const VectorType &rho)
Gradient of Psir = ar*rho w.r.t. the molar concentrations.
Definition derivs.hpp:1144
static auto get_d_ln_fugacity_coefficients_dT_constrhovec(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the derivative of the natural logarithm of fugacity coefficient of each component w....
Definition derivs.hpp:1299
static auto get_Ar01(const Model &model, const Scalar &T, const VectorType &rhovec)
Definition derivs.hpp:984
static auto get_dpdT_constrhovec(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the temperature derivative of the pressure at constant molar concentrations.
Definition derivs.hpp:1413
static auto get_ln_fugacity_coefficients2(const Model &model, const Scalar &T, const VectorType &rhovec)
Definition derivs.hpp:1248
static auto get_dpdrhovec_constT(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the molar concentration derivatives of the pressure at constant temperature.
Definition derivs.hpp:1423
static auto get_fugacity_coefficients(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the fugacity coefficient of each component.
Definition derivs.hpp:1201
static auto get_Ar00(const Model &model, const Scalar &T, const VectorType &rhovec)
Definition derivs.hpp:972
static auto get_ln_fugacity_coefficients_Trhomolefracs(const Model &model, const Scalar &T, const Scalar &rhotot, const VectorType &molefrac)
Definition derivs.hpp:1225
static auto get_ln_fugacity_coefficients1(const Model &model, const Scalar &T, const VectorType &rhovec)
Definition derivs.hpp:1237
static auto get_chempotVLE_autodiff(const Model &model, const Scalar &T, const VectorType &rho)
Calculate the chemical potential of each component.
Definition derivs.hpp:1188
static auto get_dchempotdT_autodiff(const Model &model, const Scalar &T, const VectorType &rhovec)
Calculate the temperature derivative of the chemical potential of each component.
Definition derivs.hpp:1403
static VectorType get_Psir_sigma_derivs(const Model &model, const Scalar &T, const VectorType &rhovec, const VectorType &v)
Definition derivs.hpp:1488
static auto build_Psir_gradient(const Model &model, const Scalar &T, const VectorType &rho)
Definition derivs.hpp:1165
static auto get_ATrhoXiXj_runtime(const AlphaWrapper &w, const Scalar &T, int iT, const Scalar &rho, int iD, const VectorType &molefrac, int i, int iXi, int j, int iXj)
Definition derivs.hpp:307
static auto get_Agen0n(const AlphaWrapper &w, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:536
static auto get_Ar11(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:526
static auto get_ATrhoXiXj(const AlphaWrapper &w, const Scalar &T, const Scalar &rho, const VectorType &molefrac, int i, int j)
Definition derivs.hpp:416
static auto get_Agenn0(const AlphaWrapper &w, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:564
static auto get_Aig11(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:531
static auto get_Aigxy(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:481
static auto get_ATrhoXiXjXk(const AlphaWrapper &w, const Scalar &T, const Scalar &rho, const VectorType &molefrac, int i, int j, int k)
Definition derivs.hpp:441
static auto get_Ar03(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:501
static auto get_AtaudeltaXiXj_runtime(const AlphaWrapper &w, const Scalar &tau, const int iT, const Scalar &delta, const int iD, const VectorType &molefrac, const int i, const int iXi, const int j, const int iXj)
Definition derivs.hpp:393
static auto get_Ar20(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:506
static auto get_ATrhoXi_runtime(const AlphaWrapper &w, const Scalar &T, int iT, const Scalar &rho, int iD, const VectorType &molefrac, int i, int iXi)
Definition derivs.hpp:282
static auto AlphaCaller(const AlphaWrapper &w, const S1 &T, const S2 &rho, const Vec &molefrac)
Definition derivs.hpp:124
static auto get_Ar(const int itau, const int idelta, const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:620
static auto get_Ar12(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:521
static auto get_AtaudeltaXiXjXk_runtime(const AlphaWrapper &w, const Scalar &tau, const int iT, const Scalar &delta, const int iD, const VectorType &molefrac, const int i, int iXi, const int j, const int iXj, const int k, const int iXk)
Definition derivs.hpp:401
static auto get_AtaudeltaXi_runtime(const AlphaWrapper &w, const Scalar &tau, const int iT, const Scalar &delta, const int iD, const VectorType &molefrac, const int i, const int iXi)
Definition derivs.hpp:385
static auto AlphaCaller(const AlphaWrapper &w, const S1 &T, const S2 &rho, const Vec &molefrac)
Definition derivs.hpp:128
static auto AlpharTauDeltaCaller(const AlphaWrapper &, const S1 &, const S2 &, const Vec &molefrac)
Definition derivs.hpp:137
static Scalar get_Ar01(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:491
static auto get_neff(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:677
static auto get_Arn0(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:601
static auto get_AtaudeltaXi(const AlphaWrapper &w, const Scalar &tau, const Scalar &delta, const VectorType &molefrac, const int i)
Definition derivs.hpp:337
static auto get_Ar02(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:496
static auto AlpharTauDeltaCaller(const AlphaWrapper &w, const S1 &T, const S2 &rho, const Vec &molefrac)
Definition derivs.hpp:133
static auto get_AtaudeltaXiXj(const AlphaWrapper &w, const Scalar &tau, const Scalar &delta, const VectorType &molefrac, const int i, const int j)
Definition derivs.hpp:350
static auto get_ATrhoXi(const AlphaWrapper &w, const Scalar &T, const Scalar &rho, const VectorType &molefrac, int i)
Definition derivs.hpp:248
static auto get_Ar0n(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:614
static auto get_Ar10(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:486
static auto get_Arxy(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:468
static auto get_Ar00(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:119
static auto get_AtaudeltaXiXjXk(const AlphaWrapper &w, const Scalar &tau, const Scalar &delta, const VectorType &molefrac, const int i, const int j, const int k)
Definition derivs.hpp:367
static auto get_ATrhoXiXjXk_runtime(const AlphaWrapper &w, const Scalar &T, int iT, const Scalar &rho, int iD, const VectorType &molefrac, int i, int iXi, int j, int iXj, int k, int iXk)
Definition derivs.hpp:329
static auto get_Ar30(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:511
static auto get_Ar21(const Model &model, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:516
static auto get_Agenxy(const AlphaWrapper &w, const Scalar &T, const Scalar &rho, const VectorType &molefrac)
Definition derivs.hpp:151
static auto get_dmBnvirdTm(const Model &model, const Scalar &T, const VectorType &molefrac)
Temperature derivatives of a virial coefficient.
Definition derivs.hpp:776
static auto get_B12vir(const Model &model, const Scalar &T, const VectorType &molefrac)
Calculate the cross-virial coefficient .
Definition derivs.hpp:848
static auto get_dmBnvirdTm_runtime(const int Nderiv, const int NTderiv, const Model &model, const Scalar &T, const VectorType &molefrac)
Definition derivs.hpp:809
static auto get_Bnvir(const Model &model, const Scalar &T, const VectorType &molefrac)
Definition derivs.hpp:706
static auto get_B2vir(const Model &model, const Scalar &T, const VectorType &molefrac)
Definition derivs.hpp:688
static auto get_Bnvir_runtime(const int Nderiv, const Model &model, const Scalar &T, const VectorType &molefrac)
Definition derivs.hpp:749
auto operator()(Args &&... args) const
Definition derivs.hpp:86