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

Revamped the entire C++ interface from recommendations of Tobias Loew

Compilation of user code is very fast now
parent 041db1c3
Branches
No related tags found
No related merge requests found
...@@ -6,39 +6,40 @@ ...@@ -6,39 +6,40 @@
using namespace teqp::cppinterface; using namespace teqp::cppinterface;
#include "teqp/derivs.hpp" #include "teqp/derivs.hpp"
#include "teqp/models/multifluid.hpp"
using namespace teqp; using namespace teqp;
TEST_CASE("Test C++ interface", "[C++]") TEST_CASE("Test C++ interface", "[C++]")
{ {
auto modelnovar = teqp::build_multifluid_model({ "Methane","Ethane" }, "../mycp"); auto modelnovar = teqp::build_multifluid_model({ "Methane","Ethane" }, "../mycp");
teqp::AllowedModels model = modelnovar; auto model = teqp::cppinterface::make_multifluid_model({ "Methane","Ethane" }, "../mycp");
auto z = (Eigen::ArrayXd(2) << 0.5, 0.5).finished(); auto z = (Eigen::ArrayXd(2) << 0.5, 0.5).finished();
SECTION("Ar01") { SECTION("Ar01") {
double Ar01 = get_Arxy(model, 0, 1, 300, 3, z); double Ar01 = model->get_Arxy(0, 1, 300, 3, z);
} }
SECTION("critical trace") { SECTION("critical trace") {
double Tc1 = modelnovar.redfunc.Tc(0); double Tc1 = modelnovar.redfunc.Tc(0);
auto rhovec0 = (Eigen::ArrayXd(2) << 1/modelnovar.redfunc.vc(0), 0).finished(); auto rhovec0 = (Eigen::ArrayXd(2) << 1/modelnovar.redfunc.vc(0), 0).finished();
auto cr = trace_critical_arclength_binary(model, Tc1, rhovec0); auto cr = model->trace_critical_arclength_binary(Tc1, rhovec0);
std::cout << cr.dump(1) << std::endl; std::cout << cr.dump(1) << std::endl;
} }
} }
TEST_CASE("Benchmark C++ interface", "[C++]") TEST_CASE("Benchmark C++ interface", "[C++]")
{ {
teqp::AllowedModels model = teqp::build_multifluid_model({ "Methane", "Ethane"}, "../mycp"); auto model = teqp::cppinterface::make_multifluid_model({ "Methane", "Ethane"}, "../mycp");
auto z = (Eigen::ArrayXd(2) << 0.5, 0.5).finished(); auto z = (Eigen::ArrayXd(2) << 0.5, 0.5).finished();
auto model1novar = teqp::build_multifluid_model({ "Methane" }, "../mycp"); auto model1novar = teqp::build_multifluid_model({ "Methane" }, "../mycp");
teqp::AllowedModels model1 = model1novar; auto model1 = teqp::cppinterface::make_multifluid_model({ "Methane" }, "../mycp");
auto z1 = (Eigen::ArrayXd(1) << 1).finished(); auto z1 = (Eigen::ArrayXd(1) << 1).finished();
BENCHMARK("Ar01 two components") { BENCHMARK("Ar01 two components") {
return get_Arxy(model, 0, 1, 300, 3, z); return model->get_Arxy(0, 1, 300, 3, z);
}; };
BENCHMARK("Ar01 one component w/ Arxy (runtime lookup)") { BENCHMARK("Ar01 one component w/ Arxy (C++ interface, runtime lookup)") {
return get_Arxy(model1, 0, 1, 300, 3, z1); return model->get_Arxy(0, 1, 300, 3, z1);
}; };
BENCHMARK("Ar01 one component w/ Ar01 directly") { BENCHMARK("Ar01 one component w/ Ar01 directly") {
return TDXDerivatives<decltype(model1novar)>::get_Arxy<0,1,ADBackends::autodiff>(model1novar, 300, 3, z1); return TDXDerivatives<decltype(model1novar)>::get_Arxy<0,1,ADBackends::autodiff>(model1novar, 300, 3, z1);
......
#include "teqp/models/fwd.hpp"
#include "teqpcpp.hpp" #include "teqpcpp.hpp"
#include "teqp/derivs.hpp" #include "teqp/derivs.hpp"
#include "teqp/json_builder.hpp"
#include "teqp/algorithms/critical_tracing.hpp" #include "teqp/algorithms/critical_tracing.hpp"
double teqp::cppinterface::get_Arxy(const teqp::AllowedModels& model, const int NT, const int ND, const double T, const double rho, const Eigen::ArrayXd &molefracs){ using namespace teqp;
return std::visit([&](const auto& model) { using namespace teqp::cppinterface;
using tdx = teqp::TDXDerivatives<decltype(model), double, std::decay_t<decltype(molefracs)>>;
return tdx::template get_Ar(NT, ND, model, T, rho, molefracs); namespace teqp {
}, model); namespace cppinterface {
}
class ModelImplementer : public AbstractModel {
protected:
const AllowedModels m_model;
public:
ModelImplementer(AllowedModels&& model) : m_model(model) {};
double get_Arxy(const int NT, const int ND, const double T, const double rho, const Eigen::ArrayXd& molefracs) const override {
return std::visit([&](const auto& model) {
using tdx = teqp::TDXDerivatives<decltype(model), double, Eigen::ArrayXd>;
return tdx::template get_Ar(NT, ND, model, T, rho, molefracs);
}, m_model);
}
nlohmann::json trace_critical_arclength_binary(const double T0, const Eigen::ArrayXd& rhovec0) const override {
return std::visit([&](const auto& model) {
using crit = teqp::CriticalTracing<decltype(model), double, std::decay_t<decltype(rhovec0)>>;
return crit::trace_critical_arclength_binary(model, T0, rhovec0, "");
}, m_model);
}
};
std::unique_ptr<AbstractModel> make_model(const nlohmann::json& j) {
return std::make_unique<ModelImplementer>(build_model(j));
}
nlohmann::json teqp::cppinterface::trace_critical_arclength_binary(const teqp::AllowedModels& model, const double T0, const Eigen::ArrayXd &rhovec0) { std::unique_ptr<AbstractModel> make_multifluid_model(const std::vector<std::string>& components, const std::string& coolprop_root, const std::string& BIPcollectionpath, const nlohmann::json& flags, const std::string& departurepath) {
return std::visit([&](const auto& model) { return std::make_unique<ModelImplementer>(build_multifluid_model(components, coolprop_root, BIPcollectionpath, flags, departurepath));
using crit = teqp::CriticalTracing<decltype(model), double, std::decay_t<decltype(rhovec0)>>; }
return crit::trace_critical_arclength_binary(model, T0, rhovec0, ""); }
}, model);
} }
\ No newline at end of file
#pragma once #pragma once
#include <memory>
#include <Eigen/Dense> #include <Eigen/Dense>
#include "nlohmann/json.hpp" #include "nlohmann/json.hpp"
// Pulls in the AllowedModels variant
#include "teqp/models/fwd.hpp"
namespace teqp { namespace teqp {
namespace cppinterface { namespace cppinterface {
// Wrapper functions class AbstractModel {
double get_Arxy(const teqp::AllowedModels& model, const int NT, const int ND, const double T, const double rho, const Eigen::ArrayXd& molefracs); public:
nlohmann::json trace_critical_arclength_binary(const teqp::AllowedModels& model, const double T0, const Eigen::ArrayXd& rhovec0); virtual double get_Arxy(const int, const int, const double, const double, const Eigen::ArrayXd&) const = 0;
virtual nlohmann::json trace_critical_arclength_binary(const double T0, const Eigen::ArrayXd& rhovec0) const = 0;
};
// Generic JSON-based interface where the model description is encoded as JSON
std::unique_ptr<AbstractModel> make_model(const nlohmann::json &);
// Expose factory functions for different models
// ....
std::unique_ptr<AbstractModel> make_multifluid_model(
const std::vector<std::string>& components,
const std::string& coolprop_root,
const std::string& BIPcollectionpath = {},
const nlohmann::json& flags = {},
const std::string& departurepath = {}
);
} }
} }
\ No newline at end of file
...@@ -8,11 +8,22 @@ ...@@ -8,11 +8,22 @@
* all the templates every time this file is compiled * all the templates every time this file is compiled
*/ */
#include "teqpcpp.hpp" #include "teqpcpp.hpp"
#include <iostream>
int main() { int main() {
teqp::AllowedModels model = teqp::build_multifluid_model({ "Methane", "Ethane" }, "../mycp"); nlohmann::json j = {
{"kind", "multifluid"},
{"model", {
{"components", {"../mycp/dev/fluids/Methane.json","../mycp/dev/fluids/Ethane.json"}},
{"BIP", "../mycp/dev/mixtures/mixture_binary_pairs.json"},
{"departure", "../mycp/dev/mixtures/mixture_departure_functions.json"}
}
}};
//std::cout << j.dump(2);
auto am = teqp::cppinterface::make_model(j);
auto z = (Eigen::ArrayXd(2) << 0.5, 0.5).finished(); auto z = (Eigen::ArrayXd(2) << 0.5, 0.5).finished();
double Ar01 = teqp::cppinterface::get_Arxy(model, 0, 1, 300, 3, z); double Ar01 = am->get_Arxy(0, 1, 300, 3, z);
std::cout << Ar01 << std::endl;
} }
\ No newline at end of file
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment