teqp 0.19.1
Loading...
Searching...
No Matches
multifluid_mutant.hpp
Go to the documentation of this file.
1#pragma once
2
3
4namespace teqp {
5
11 template<typename DepartureFunction, typename BaseClass>
13
14 private:
15 std::string meta = "";
16
17 public:
18 const BaseClass& base;
20 const DepartureFunction dep;
21
22 template<class VecType>
23 auto R(const VecType& molefrac) const { return base.R(molefrac); }
24
25 MultiFluidAdapter(const BaseClass& base, ReducingFunctions&& redfunc, DepartureFunction&& depfunc) : base(base), redfunc(redfunc), dep(depfunc) {};
26
28 void set_meta(const std::string& m) { meta = m; }
30 auto get_meta() const { return meta; }
32 const std::variant<double, std::string> get_BIP(const std::size_t &i, const std::size_t &j, const std::string& key) const{
33 if (key == "F" || key == "Fij"){
34 auto F = dep.get_F();
35 if (0 <= i && i < F.rows() && 0 <= j && j < F.cols()){
36 return F(i,j);
37 }
38 }
39 return redfunc.get_BIP(i, j, key);
40 }
41
42 template<typename TType, typename RhoType, typename MoleFracType>
43 auto alphar(const TType& T,
44 const RhoType& rho,
45 const MoleFracType& molefrac) const
46 {
47 auto Tred = forceeval(redfunc.get_Tr(molefrac));
48 auto rhored = forceeval(redfunc.get_rhor(molefrac));
49 auto delta = forceeval(rho / rhored);
50 auto tau = forceeval(Tred / T);
51 auto val = base.corr.alphar(tau, delta, molefrac) + dep.alphar(tau, delta, molefrac);
52 return forceeval(val);
53 }
54 };
55
56 template<class Model>
57 auto build_multifluid_mutant(const Model& model, const nlohmann::json& jj) {
58
59 auto N = model.redfunc.Tc.size();
60
61 // Allocate the matrices of default models and F factors
62 Eigen::MatrixXd F(N, N); F.setZero();
63 std::vector<std::vector<DepartureTerms>> funcs(N);
64 for (auto i = 0; i < N; ++i) { funcs[i].resize(N); }
65
66 // Build the F and departure function matrix
67 for (auto i = 0; i < N; ++i) {
68 for (auto j = i; j < N; ++j) {
69 if (i == j) {
70 funcs[i][i].add_term(NullEOSTerm());
71 }
72 else {
73 // Extract the given entry
74 auto entry = jj[std::to_string(i)][std::to_string(j)];
75
76 // Set the entry in the matrix of F and departure functions
77 auto dep = entry.at("departure");
78 auto BIP = entry.at("BIP");
79 F(i, j) = BIP.at("Fij");
80 F(j, i) = F(i, j);
81 funcs[i][j] = build_departure_function(dep);
82 funcs[j][i] = build_departure_function(dep);
83 }
84 }
85 }
86
87 // Determine what sort of reducing function is to be used
88 auto get_reducing = [&](const auto& deptype) {
89 auto red = model.redfunc;
90 auto Tc = red.Tc, vc = red.vc;
91 if (deptype == "invariant") {
92 using mat = std::decay_t<decltype(MultiFluidInvariantReducingFunction::phiT)>;
93 mat phiT = mat::Zero(N, N), lambdaT = mat::Zero(N, N), phiV = mat::Zero(N, N), lambdaV = mat::Zero(N, N);
94
95 for (auto i = 0; i < N; ++i) {
96 for (auto j = i+1; j < N; ++j) {
97 // Extract the given entry
98 auto entry = jj.at(std::to_string(i)).at(std::to_string(j));
99
100 auto BIP = entry.at("BIP");
101 // Set the reducing function parameters in the copy
102 phiT(i, j) = BIP.at("phiT");
103 phiT(j, i) = phiT(i, j);
104 lambdaT(i, j) = BIP.at("lambdaT");
105 lambdaT(j, i) = -lambdaT(i, j);
106
107 phiV(i, j) = BIP.at("phiV");
108 phiV(j, i) = phiV(i, j);
109 lambdaV(i, j) = BIP.at("lambdaV");
110 lambdaV(j, i) = -lambdaV(i, j);
111 }
112 }
113 return ReducingFunctions(MultiFluidInvariantReducingFunction(phiT, lambdaT, phiV, lambdaV, Tc, vc));
114 }
115 else {
116 using mat = std::decay_t<decltype(MultiFluidReducingFunction::betaT)>;
117 mat betaT = mat::Zero(N, N), gammaT = mat::Zero(N, N), betaV = mat::Zero(N, N), gammaV = mat::Zero(N, N);
118
119 for (auto i = 0; i < N; ++i) {
120 for (auto j = i+1; j < N; ++j) {
121 // Extract the given entry
122 auto entry = jj.at(std::to_string(i)).at(std::to_string(j));
123 auto BIP = entry.at("BIP");
124 // Set the reducing function parameters in the copy
125 betaT(i, j) = BIP.at("betaT");
126 betaT(j, i) = 1 / betaT(i, j);
127 betaV(i, j) = BIP.at("betaV");
128 betaV(j, i) = 1 / betaV(i, j);
129 gammaT(i, j) = BIP.at("gammaT"); gammaT(j, i) = gammaT(i, j);
130 gammaV(i, j) = BIP.at("gammaV"); gammaV(j, i) = gammaV(i, j);
131 }
132 }
133 return ReducingFunctions(MultiFluidReducingFunction(betaT, gammaT, betaV, gammaV, Tc, vc));
134 }
135 };
136 std::string deptype = (jj.at("0").at("1").at("BIP").contains("type")) ? jj.at("0").at("1").at("BIP")["type"] : "";
137 ReducingFunctions newred = get_reducing(deptype);
138
139 auto newdep = DepartureContribution(std::move(F), std::move(funcs));
140 auto mfa = MultiFluidAdapter(model, std::move(newred), std::move(newdep));
142 mfa.set_meta(jj.dump());
143 return mfa;
144 }
145
146}
MultiFluidAdapter(const BaseClass &base, ReducingFunctions &&redfunc, DepartureFunction &&depfunc)
auto R(const VecType &molefrac) const
auto alphar(const TType &T, const RhoType &rho, const MoleFracType &molefrac) const
auto get_meta() const
Get the metadata stored in string form.
const std::variant< double, std::string > get_BIP(const std::size_t &i, const std::size_t &j, const std::string &key) const
Return a binary interaction parameter.
const ReducingFunctions redfunc
const DepartureFunction dep
void set_meta(const std::string &m)
Store some sort of metadata in string form (perhaps a JSON representation of the model?...
auto get_rhor(const MoleFractions &molefracs) const
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 build_multifluid_mutant(const Model &model, const nlohmann::json &jj)
auto build_departure_function(const nlohmann::json &j)
ReducingTermContainer< MultiFluidReducingFunction, MultiFluidInvariantReducingFunction > ReducingFunctions
auto forceeval(T &&expr)
Definition types.hpp:46