Skip to content

Commit

Permalink
Add get_names for PC-SAFT
Browse files Browse the repository at this point in the history
See #51
  • Loading branch information
ianhbell committed Sep 15, 2023
1 parent a5439f2 commit b93a8fb
Show file tree
Hide file tree
Showing 3 changed files with 17 additions and 1 deletion.
10 changes: 9 additions & 1 deletion include/teqp/models/pcsaft.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -327,6 +327,13 @@ class PCSAFTMixture {
}
return PCSAFTHardChainContribution(m, mminus1, sigma_Angstrom, epsilon_over_k, kmat);
}
auto extract_names(const std::vector<SAFTCoeffs> &coeffs){
std::vector<std::string> names;
for (const auto& c: coeffs){
names.push_back(c.name);
}
return names;
}
auto build_dipolar(const std::vector<SAFTCoeffs> &coeffs) -> std::optional<PCSAFTDipolarContribution>{
Eigen::ArrayXd mustar2(coeffs.size()), nmu(coeffs.size());
auto i = 0;
Expand Down Expand Up @@ -357,7 +364,7 @@ class PCSAFTMixture {
}
public:
PCSAFTMixture(const std::vector<std::string> &names, const Eigen::ArrayXXd& kmat = {}) : PCSAFTMixture(get_coeffs_from_names(names), kmat){};
PCSAFTMixture(const std::vector<SAFTCoeffs> &coeffs, const Eigen::ArrayXXd &kmat = {}) : kmat(kmat), hardchain(build_hardchain(coeffs)), dipolar(build_dipolar(coeffs)), quadrupolar(build_quadrupolar(coeffs)) {};
PCSAFTMixture(const std::vector<SAFTCoeffs> &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
Expand All @@ -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";
Expand Down
1 change: 1 addition & 0 deletions interface/pybind11_wrapper.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -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<PCSAFT_t>(o).get_m(); }), obj));
setattr("get_sigma_Angstrom", MethodType(py::cpp_function([](py::object& o){ return get_typed<PCSAFT_t>(o).get_sigma_Angstrom(); }), obj));
setattr("get_epsilon_over_k_K", MethodType(py::cpp_function([](py::object& o){ return get_typed<PCSAFT_t>(o).get_epsilon_over_k_K(); }), obj));
setattr("get_names", MethodType(py::cpp_function([](py::object& o){ return get_typed<PCSAFT_t>(o).get_names(); }), obj));
setattr("max_rhoN", MethodType(py::cpp_function([](py::object& o, double T, REArrayd& molefrac){ return get_typed<PCSAFT_t>(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<PCSAFT_t>(o).get_kmat(); }), obj));
}
Expand Down
7 changes: 7 additions & 0 deletions src/tests/catch_test_PCSAFT.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -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<std::string> names = { "Methane" };
auto model = PCSAFTMixture(names);
CHECK(model.get_names()[0] == "Methane");
}

TEST_CASE("Check 0n derivatives", "[PCSAFT]")
{
std::vector<std::string> names = { "Methane", "Ethane" };
Expand Down

0 comments on commit b93a8fb

Please sign in to comment.