From b93a8fbaf211f1dd96d33205a84247dd677e4c81 Mon Sep 17 00:00:00 2001 From: Ian Bell Date: Fri, 15 Sep 2023 17:32:10 -0400 Subject: [PATCH] Add get_names for PC-SAFT See #51 --- include/teqp/models/pcsaft.hpp | 10 +++++++++- interface/pybind11_wrapper.cpp | 1 + src/tests/catch_test_PCSAFT.cxx | 7 +++++++ 3 files changed, 17 insertions(+), 1 deletion(-) diff --git a/include/teqp/models/pcsaft.hpp b/include/teqp/models/pcsaft.hpp index 9251efff..0fb5d468 100644 --- a/include/teqp/models/pcsaft.hpp +++ b/include/teqp/models/pcsaft.hpp @@ -327,6 +327,13 @@ class PCSAFTMixture { } return PCSAFTHardChainContribution(m, mminus1, sigma_Angstrom, epsilon_over_k, kmat); } + auto extract_names(const std::vector &coeffs){ + std::vector names; + for (const auto& c: coeffs){ + names.push_back(c.name); + } + return names; + } auto build_dipolar(const std::vector &coeffs) -> std::optional{ Eigen::ArrayXd mustar2(coeffs.size()), nmu(coeffs.size()); auto i = 0; @@ -357,7 +364,7 @@ class PCSAFTMixture { } public: PCSAFTMixture(const std::vector &names, const Eigen::ArrayXXd& kmat = {}) : PCSAFTMixture(get_coeffs_from_names(names), kmat){}; - PCSAFTMixture(const std::vector &coeffs, const Eigen::ArrayXXd &kmat = {}) : kmat(kmat), hardchain(build_hardchain(coeffs)), dipolar(build_dipolar(coeffs)), quadrupolar(build_quadrupolar(coeffs)) {}; + PCSAFTMixture(const std::vector &coeffs, const Eigen::ArrayXXd &kmat = {}) : names(extract_names(coeffs)), kmat(kmat), hardchain(build_hardchain(coeffs)), dipolar(build_dipolar(coeffs)), quadrupolar(build_quadrupolar(coeffs)) {}; // PCSAFTMixture( const PCSAFTMixture& ) = delete; // non construction-copyable PCSAFTMixture& operator=( const PCSAFTMixture& ) = delete; // non copyable @@ -366,6 +373,7 @@ class PCSAFTMixture { auto get_sigma_Angstrom() const { return sigma_Angstrom; } auto get_epsilon_over_k_K() const { return epsilon_over_k; } auto get_kmat() const { return kmat; } + auto get_names() const { return names;} auto print_info() { std::string s = std::string("i m sigma / A e/kB / K \n ++++++++++++++") + "\n"; diff --git a/interface/pybind11_wrapper.cpp b/interface/pybind11_wrapper.cpp index d282689e..1115e1ca 100644 --- a/interface/pybind11_wrapper.cpp +++ b/interface/pybind11_wrapper.cpp @@ -137,6 +137,7 @@ void attach_model_specific_methods(py::object& obj){ setattr("get_m", MethodType(py::cpp_function([](py::object& o){ return get_typed(o).get_m(); }), obj)); setattr("get_sigma_Angstrom", MethodType(py::cpp_function([](py::object& o){ return get_typed(o).get_sigma_Angstrom(); }), obj)); setattr("get_epsilon_over_k_K", MethodType(py::cpp_function([](py::object& o){ return get_typed(o).get_epsilon_over_k_K(); }), obj)); + setattr("get_names", MethodType(py::cpp_function([](py::object& o){ return get_typed(o).get_names(); }), obj)); setattr("max_rhoN", MethodType(py::cpp_function([](py::object& o, double T, REArrayd& molefrac){ return get_typed(o).max_rhoN(T, molefrac); }, "self"_a, "T"_a, "molefrac"_a), obj)); setattr("get_kmat", MethodType(py::cpp_function([](py::object& o){ return get_typed(o).get_kmat(); }), obj)); } diff --git a/src/tests/catch_test_PCSAFT.cxx b/src/tests/catch_test_PCSAFT.cxx index 65020485..e5d9f16a 100644 --- a/src/tests/catch_test_PCSAFT.cxx +++ b/src/tests/catch_test_PCSAFT.cxx @@ -26,6 +26,13 @@ TEST_CASE("Single alphar check value", "[PCSAFT]") CHECK(tdx::get_Ar00(model, T, Dmolar, z) == Approx(-0.032400020930842724)); } +TEST_CASE("Check get_names", "[PCSAFT]") +{ + std::vector names = { "Methane" }; + auto model = PCSAFTMixture(names); + CHECK(model.get_names()[0] == "Methane"); +} + TEST_CASE("Check 0n derivatives", "[PCSAFT]") { std::vector names = { "Methane", "Ethane" };