9 inline auto get_BIPdep(
const nlohmann::json& collection,
const std::vector<std::string>& identifiers,
const nlohmann::json& flags) {
10 if (!collection.is_array()){
13 if (collection.size() > 0 && !collection[0].is_object()){
18 if (flags.contains(
"force-estimate")) {
19 std::string scheme = flags.at(
"estimate");
20 if (scheme ==
"Lorentz-Berthelot") {
21 return std::make_tuple(nlohmann::json({
22 {
"betaT", 1.0}, {
"gammaT", 1.0}, {
"betaV", 1.0}, {
"gammaV", 1.0}, {
"F", 0.0}
26 throw std::invalid_argument(
"estimation scheme is not understood:" + scheme);
31 auto toupper = [](
const std::string s) {
auto data = s; std::for_each(data.begin(), data.end(), [](
char& c) { c = ::toupper(c); });
return data; };
33 std::string comp0 = toupper(identifiers[0]);
34 std::string comp1 = toupper(identifiers[1]);
36 for (
auto& el : collection) {
37 if (!el.contains(
"hash1")){
continue; }
38 std::string name1 = toupper(el.at(
"hash1"));
39 std::string name2 = toupper(el.at(
"hash2"));
40 if (comp0 == name1 && comp1 == name2) {
41 return std::make_tuple(el,
false);
43 if (comp0 == name2 && comp1 == name1) {
44 return std::make_tuple(el,
true);
48 for (
auto& el : collection) {
49 std::string name1 = toupper(el.at(
"Name1"));
50 std::string name2 = toupper(el.at(
"Name2"));
51 if (comp0 == name1 && comp1 == name2) {
52 return std::make_tuple(el,
false);
54 if (comp0 == name2 && comp1 == name1) {
55 return std::make_tuple(el,
true);
59 for (
auto& el : collection) {
60 std::string CAS1 = el.at(
"CAS1");
61 std::string CAS2 = el.at(
"CAS2");
62 if (identifiers[0] == CAS1 && identifiers[1] == CAS2) {
63 return std::make_tuple(el,
false);
65 if (identifiers[0] == CAS2 && identifiers[1] == CAS1) {
66 return std::make_tuple(el,
true);
71 if (flags.contains(
"estimate")) {
72 std::string scheme = flags.at(
"estimate");
73 if (scheme ==
"Lorentz-Berthelot") {
74 return std::make_tuple(nlohmann::json({
75 {
"betaT", 1.0}, {
"gammaT", 1.0}, {
"betaV", 1.0}, {
"gammaV", 1.0}, {
"F", 0.0}
79 throw std::invalid_argument(
"estimation scheme is not understood:" + scheme);
83 throw std::invalid_argument(
"Can't match the binary pair for: " + identifiers[0] +
"/" + identifiers[1]);
88 inline auto get_binary_interaction_double(
const nlohmann::json& collection,
const std::vector<std::string>& identifiers,
const nlohmann::json& flags,
const std::vector<double>& Tc,
const std::vector<double>& vc) {
89 auto [el, swap_needed] =
get_BIPdep(collection, identifiers, flags);
91 double betaT, gammaT, betaV, gammaV;
92 if (el.contains(
"betaT") && el.contains(
"gammaT") && el.contains(
"betaV") && el.contains(
"gammaV")) {
93 betaT = el[
"betaT"]; gammaT = el[
"gammaT"]; betaV = el[
"betaV"]; gammaV = el[
"gammaV"];
100 else if (el.contains(
"xi") && el.contains(
"zeta")) {
101 double xi = el[
"xi"], zeta = el[
"zeta"];
102 gammaT = 0.5 * (Tc[0] + Tc[1] + xi) / (2 * sqrt(Tc[0] * Tc[1]));
103 gammaV = 4.0 * (vc[0] + vc[1] + zeta) / (0.25 *
pow(1 /
pow(1 / vc[0], 1.0 / 3.0) + 1 /
pow(1 / vc[1], 1.0 / 3.0), 3));
108 throw std::invalid_argument(
"Could not understand what to do with this binary model specification: " + el.dump());
110 return std::make_tuple(betaT, gammaT, betaV, gammaV);
114 template <
typename Tcvec,
typename vcvec>
115 inline auto get_BIP_matrices(
const nlohmann::json& collection,
const std::vector<std::string>& components,
const nlohmann::json& flags,
const Tcvec& Tc,
const vcvec& vc) {
116 Eigen::MatrixXd betaT, gammaT, betaV, gammaV, YT, Yv;
117 auto N = components.size();
118 betaT.resize(N, N); betaT.setZero();
119 gammaT.resize(N, N); gammaT.setZero();
120 betaV.resize(N, N); betaV.setZero();
121 gammaV.resize(N, N); gammaV.setZero();
122 for (
auto i = 0U; i < N; ++i) {
123 for (
auto j = i + 1; j < N; ++j) {
124 auto [betaT_, gammaT_, betaV_, gammaV_] =
get_binary_interaction_double(collection, { components[i], components[j] }, flags, { Tc[i], Tc[j] }, { vc[i], vc[j] });
125 betaT(i, j) = betaT_; betaT(j, i) = 1.0 / betaT(i, j);
126 gammaT(i, j) = gammaT_; gammaT(j, i) = gammaT(i, j);
127 betaV(i, j) = betaV_; betaV(j, i) = 1.0 / betaV(i, j);
128 gammaV(i, j) = gammaV_; gammaV(j, i) = gammaV(i, j);
131 return std::make_tuple(betaT, gammaT, betaV, gammaV);
138 inline auto get_Tcvc(
const std::vector<nlohmann::json>& pureJSON) {
139 Eigen::ArrayXd Tc(pureJSON.size()), vc(pureJSON.size());
141 for (
auto& j : pureJSON) {
142 auto red = j[
"EOS"][0][
"STATES"][
"reducing"];
143 double Tc_ = red.at(
"T");
144 double rhoc_ = red.at(
"rhomolar");
149 return std::make_tuple(Tc, vc);
153 inline auto get_F_matrix(
const nlohmann::json& collection,
const std::vector<std::string>& identifiers,
const nlohmann::json& flags) {
154 auto N = identifiers.size();
155 Eigen::MatrixXd F(N, N);
156 for (
auto i = 0U; i < N; ++i) {
158 for (
auto j = i + 1; j < N; ++j) {
159 auto [el, swap_needed] =
get_BIPdep(collection, { identifiers[i], identifiers[j] }, flags);
176 Eigen::MatrixXd YT, Yv;
182 template<
typename ArrayLike>
184 const Eigen::MatrixXd&
betaT,
const Eigen::MatrixXd&
gammaT,
185 const Eigen::MatrixXd&
betaV,
const Eigen::MatrixXd&
gammaV,
186 const ArrayLike&
Tc,
const ArrayLike&
vc)
191 YT.resize(N, N); YT.setZero();
192 Yv.resize(N, N); Yv.setZero();
193 for (
auto i = 0; i < N; ++i) {
194 for (
auto j = i + 1; j < N; ++j) {
203 template <
typename MoleFractions>
204 auto Y(
const MoleFractions& z,
const Eigen::ArrayXd& Yc,
const Eigen::MatrixXd& beta,
const Eigen::MatrixXd& Yij)
const {
207 typename MoleFractions::value_type sum1 = 0.0;
208 for (
auto i = 0U; i < N; ++i) {
209 sum1 = sum1 +
pow2(z[i]) * Yc[i];
212 typename MoleFractions::value_type sum2 = 0.0;
213 for (
auto i = 0U; i < N - 1; ++i) {
214 for (
auto j = i + 1; j < N; ++j) {
215 sum2 = sum2 + 2.0 * z[i] * z[j] * (z[i] + z[j]) / (
pow2(beta(i, j)) * z[i] + z[j]) * Yij(i, j);
222 template<
typename MoleFractions>
auto get_Tr(
const MoleFractions& molefracs)
const {
return Y(molefracs,
Tc,
betaT, YT); }
223 template<
typename MoleFractions>
auto get_rhor(
const MoleFractions& molefracs)
const {
return 1.0 /
Y(molefracs,
vc,
betaV, Yv); }
225 const auto&
get_mat(
const std::string& key)
const {
226 if (key ==
"betaT"){
return betaT; }
227 if (key ==
"gammaT"){
return gammaT; }
228 if (key ==
"betaV"){
return betaV; }
229 if (key ==
"gammaV"){
return gammaV; }
230 throw std::invalid_argument(
"variable is not understood: " + key);
232 auto get_BIP(
const std::size_t& i,
const std::size_t& j,
const std::string& key)
const {
233 const auto& mat =
get_mat(key);
234 if (i <
static_cast<std::size_t
>(mat.rows()) && j <
static_cast<std::size_t
>(mat.cols())){
238 throw std::invalid_argument(
"Indices are out of bounds");
246 Eigen::MatrixXd YT, Yv;
252 template<
typename ArrayLike>
254 const Eigen::MatrixXd&
phiT,
const Eigen::MatrixXd&
lambdaT,
255 const Eigen::MatrixXd&
phiV,
const Eigen::MatrixXd&
lambdaV,
256 const ArrayLike&
Tc,
const ArrayLike&
vc)
261 YT.resize(N, N); YT.setZero();
262 Yv.resize(N, N); Yv.setZero();
263 for (
auto i = 0; i < N; ++i) {
264 for (
auto j = 0; j < N; ++j) {
265 YT(i, j) = sqrt(
Tc[i] *
Tc[j]);
266 YT(j, i) = sqrt(
Tc[i] *
Tc[j]);
267 Yv(i, j) = 1.0 / 8.0 *
pow3(cbrt(
vc[i]) + cbrt(
vc[j]));
268 Yv(j, i) = 1.0 / 8.0 *
pow3(cbrt(
vc[i]) + cbrt(
vc[j]));
272 template <
typename MoleFractions>
273 auto Y(
const MoleFractions& z,
const Eigen::MatrixXd& phi,
const Eigen::MatrixXd& lambda,
const Eigen::MatrixXd& Yij)
const {
275 typename MoleFractions::value_type sum = 0.0;
276 for (
auto i = 0U; i < N; ++i) {
277 for (
auto j = 0U; j < N; ++j) {
278 auto contrib = z[i] * z[j] * (phi(i, j) + z[j] * lambda(i, j)) * Yij(i, j);
284 template<
typename MoleFractions>
auto get_Tr(
const MoleFractions& molefracs)
const {
return Y(molefracs,
phiT,
lambdaT, YT); }
285 template<
typename MoleFractions>
auto get_rhor(
const MoleFractions& molefracs)
const {
return 1.0 /
Y(molefracs,
phiV,
lambdaV, Yv); }
287 const auto&
get_mat(
const std::string& key)
const {
288 if (key ==
"phiT"){
return phiT; }
289 if (key ==
"lambdaT"){
return lambdaT; }
290 if (key ==
"phiV"){
return phiV; }
291 if (key ==
"lambdaV"){
return lambdaV; }
292 throw std::invalid_argument(
"variable is not understood: " + key);
294 auto get_BIP(
const std::size_t& i,
const std::size_t& j,
const std::string& key)
const {
295 const auto& mat =
get_mat(key);
296 if (i <
static_cast<std::size_t
>(mat.rows()) && j <
static_cast<std::size_t
>(mat.cols())){
300 throw std::invalid_argument(
"Indices are out of bounds");
306 template<
typename... Args>
309 const std::variant<Args...> term;
310 auto get_Tc()
const {
return std::visit([](
const auto& t) {
return std::cref(t.Tc); }, term); }
311 auto get_vc()
const {
return std::visit([](
const auto& t) {
return std::cref(t.vc); }, term); }
315 template<
typename Instance>
318 template <
typename MoleFractions>
319 auto get_Tr(
const MoleFractions& molefracs)
const {
320 return std::visit([&](
auto& t) {
return t.get_Tr(molefracs); }, term);
323 template <
typename MoleFractions>
324 auto get_rhor(
const MoleFractions& molefracs)
const {
325 return std::visit([&](
auto& t) {
return t.get_rhor(molefracs); }, term);
328 auto get_BIP(
const std::size_t& i,
const std::size_t& j,
const std::string& key)
const {
329 return std::visit([&](
auto& t) {
return t.get_BIP(i, j, key); }, term);
const Eigen::MatrixXd lambdaT
MultiFluidInvariantReducingFunction(const Eigen::MatrixXd &phiT, const Eigen::MatrixXd &lambdaT, const Eigen::MatrixXd &phiV, const Eigen::MatrixXd &lambdaV, const ArrayLike &Tc, const ArrayLike &vc)
auto get_rhor(const MoleFractions &molefracs) const
auto get_Tr(const MoleFractions &molefracs) const
const Eigen::MatrixXd phiV
auto get_BIP(const std::size_t &i, const std::size_t &j, const std::string &key) const
const auto & get_mat(const std::string &key) const
const Eigen::MatrixXd lambdaV
auto Y(const MoleFractions &z, const Eigen::MatrixXd &phi, const Eigen::MatrixXd &lambda, const Eigen::MatrixXd &Yij) const
const Eigen::MatrixXd phiT
auto get_Tr(const MoleFractions &molefracs) const
const Eigen::MatrixXd gammaV
auto get_BIP(const std::size_t &i, const std::size_t &j, const std::string &key) const
const Eigen::MatrixXd betaT
auto Y(const MoleFractions &z, const Eigen::ArrayXd &Yc, const Eigen::MatrixXd &beta, const Eigen::MatrixXd &Yij) const
auto get_rhor(const MoleFractions &molefracs) const
MultiFluidReducingFunction(const Eigen::MatrixXd &betaT, const Eigen::MatrixXd &gammaT, const Eigen::MatrixXd &betaV, const Eigen::MatrixXd &gammaV, const ArrayLike &Tc, const ArrayLike &vc)
const Eigen::MatrixXd betaV
const auto & get_mat(const std::string &key) const
const Eigen::MatrixXd gammaT
auto get_rhor(const MoleFractions &molefracs) const
ReducingTermContainer(const Instance &instance)
auto get_Tr(const MoleFractions &molefracs) const
auto get_BIP(const std::size_t &i, const std::size_t &j, const std::string &key) const
auto get_binary_interaction_double(const nlohmann::json &collection, const std::vector< std::string > &identifiers, const nlohmann::json &flags, const std::vector< double > &Tc, const std::vector< double > &vc)
Get the binary interaction parameters for a given binary pair.
auto get_Tcvc(const std::vector< nlohmann::json > &pureJSON)
auto get_F_matrix(const nlohmann::json &collection, const std::vector< std::string > &identifiers, const nlohmann::json &flags)
Get the matrix F of Fij factors multiplying the departure functions.
auto get_BIPdep(const nlohmann::json &collection, const std::vector< std::string > &identifiers, const nlohmann::json &flags)
auto get_BIP_matrices(const nlohmann::json &collection, const std::vector< std::string > &components, const nlohmann::json &flags, const Tcvec &Tc, const vcvec &vc)
Build the matrices of betaT, gammaT, betaV, gammaV for the multi-fluid model.
auto pow(const double &x, const double &e)