Skip to content
Snippets Groups Projects
Commit cf85ee07 authored by Ian Bell's avatar Ian Bell
Browse files

Possible now also to modify Fij (though you can't fiddle with the departure function)

parent b213242f
No related branches found
No related tags found
No related merge requests found
...@@ -704,12 +704,71 @@ auto build_BIPmodified(Model& model, const nlohmann::json& j) { ...@@ -704,12 +704,71 @@ auto build_BIPmodified(Model& model, const nlohmann::json& j) {
return MultiFluidReducingFunctionAdapter(model, std::move(newred)); return MultiFluidReducingFunctionAdapter(model, std::move(newred));
} }
/**
This class holds a lightweight reference to the core parts of the model
The reducing and departure functions are moved into this class, while the donor class is used for the corresponding states portion
*/
template<typename ReducingFunction, typename DepartureFunction, typename BaseClass>
class MultiFluidAdapter {
public:
const BaseClass& base;
const ReducingFunction redfunc;
const DepartureFunction depfunc;
template<class VecType>
auto R(const VecType& molefrac) const { return base.R(molefrac); }
MultiFluidAdapter(const BaseClass& base, ReducingFunction&& redfunc, DepartureFunction &&depfunc) : base(base), redfunc(redfunc), depfunc(depfunc) {};
template<typename TType, typename RhoType, typename MoleFracType>
auto alphar(const TType& T,
const RhoType& rho,
const MoleFracType& molefrac) const
{
auto Tred = forceeval(redfunc.get_Tr(molefrac));
auto rhored = forceeval(redfunc.get_rhor(molefrac));
auto delta = forceeval(rho / rhored);
auto tau = forceeval(Tred / T);
auto val = base.corr.alphar(tau, delta, molefrac) + depfunc.alphar(tau, delta, molefrac);
return forceeval(val);
}
};
template<class Model>
auto build_multifluid_mutant(Model& model, const nlohmann::json& j) {
auto red = model.redfunc;
auto betaT = red.betaT;
betaT(0, 1) = j["betaT"];
betaT(1, 0) = 1 / betaT(0, 1);
auto betaV = red.betaV;
betaV(0, 1) = j["betaV"];
betaV(1, 0) = 1 / betaV(0, 1);
auto gammaT = red.gammaT, gammaV = red.gammaV;
gammaT(0, 1) = j["gammaT"]; gammaT(1, 0) = gammaT(0, 1);
gammaV(0, 1) = j["gammaV"]; gammaV(1, 0) = gammaV(0, 1);
auto Tc = red.Tc, vc = red.vc;
auto newred = MultiFluidReducingFunction(betaT, gammaT, betaV, gammaV, Tc, vc);
if (j.contains("Fij") && j["Fij"] != 0.0) {
throw std::invalid_argument("Don't support Fij != 0 for now");
}
auto N = 2;
// Allocate the matrix with default models
Eigen::MatrixXd F(2, 2); F.setZero();
std::vector<std::vector<MultiFluidDepartureFunction>> funcs(N);
for (auto i = 0; i < funcs.size(); ++i) {
funcs[i].resize(funcs.size());
for (auto j = 0; j < N; ++j) {
funcs[i][j].set_type("none");
}
}
auto newdep = DepartureContribution(std::move(F), std::move(funcs));
return MultiFluidAdapter(model, std::move(newred), std::move(newdep));
}
class DummyEOS { class DummyEOS {
......
...@@ -6,9 +6,19 @@ ...@@ -6,9 +6,19 @@
void add_multifluid_mutant(py::module& m) { void add_multifluid_mutant(py::module& m) {
using MultiFluid = decltype(build_multifluid_model(std::vector<std::string>{"", ""}, "", "")); using MultiFluid = decltype(build_multifluid_model(std::vector<std::string>{"", ""}, "", ""));
m.def("build_BIPmodified", &build_BIPmodified<MultiFluid>); {
using RedType = std::decay_t<decltype(MultiFluid::redfunc)>; m.def("build_BIPmodified", &build_BIPmodified<MultiFluid>);
using BIPmod = MultiFluidReducingFunctionAdapter<RedType, MultiFluid>; using RedType = std::decay_t<decltype(MultiFluid::redfunc)>;
auto wMFBIP = py::class_<BIPmod>(m, "MultiFluidBIP"); using BIPmod = MultiFluidReducingFunctionAdapter<RedType, MultiFluid>;
add_derivatives<BIPmod>(m, wMFBIP); auto wMFBIP = py::class_<BIPmod>(m, "MultiFluidBIP");
add_derivatives<BIPmod>(m, wMFBIP);
}
{
m.def("build_multifluid_mutant", &build_multifluid_mutant<MultiFluid>);
using RedType = std::decay_t<decltype(MultiFluid::redfunc)>;
using DepType = std::decay_t<decltype(MultiFluid::dep)>;
using BIPmod = MultiFluidAdapter<RedType, DepType, MultiFluid>;
auto wMFBIP = py::class_<BIPmod>(m, "MultiFluidMutant");
add_derivatives<BIPmod>(m, wMFBIP);
}
} }
\ No newline at end of file
...@@ -159,8 +159,11 @@ int main(){ ...@@ -159,8 +159,11 @@ int main(){
//time_calls(coolprop_root, BIPcollection); //time_calls(coolprop_root, BIPcollection);
{ {
nlohmann::json flags = { {"estimate", true} }; nlohmann::json flags = { {"estimate", true},{"another","key"} };
auto model = build_multifluid_model({ "Ethane", "R1234ze(E)" }, coolprop_root, BIPcollection, flags); auto model = build_multifluid_model({ "Ethane", "R1234ze(E)" }, coolprop_root, BIPcollection, flags);
nlohmann::json j = { {"betaT", 1.0},{"gammaT", 1.0},{"betaV", 1.0},{"gammaV", 1.0},{"Fij", 0.0} };
auto mutant = build_mutant(model, j);
} }
{ {
auto model = build_multifluid_model({ "methane", "ethane" }, coolprop_root, BIPcollection); auto model = build_multifluid_model({ "methane", "ethane" }, coolprop_root, BIPcollection);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment