From 1660167bd7dadec072a3d163ae2daff635bfae6c Mon Sep 17 00:00:00 2001 From: Ray Speth Date: Thu, 8 Sep 2016 22:24:38 -0400 Subject: [PATCH] Allow adding species / reactions from Python This enables the base functionality, but does not protect against any erroneous usage. --- interfaces/cython/cantera/kinetics.pyx | 4 +++ .../cython/cantera/test/test_kinetics.py | 28 +++++++++++++++++++ interfaces/cython/cantera/test/test_thermo.py | 18 ++++++++++++ interfaces/cython/cantera/thermo.pyx | 11 ++++++++ 4 files changed, 61 insertions(+) diff --git a/interfaces/cython/cantera/kinetics.pyx b/interfaces/cython/cantera/kinetics.pyx index a849ce1da5..a98eb29974 100644 --- a/interfaces/cython/cantera/kinetics.pyx +++ b/interfaces/cython/cantera/kinetics.pyx @@ -97,6 +97,10 @@ cdef class Kinetics(_SolutionBase): """ self.kinetics.modifyReaction(irxn, rxn._reaction) + def add_reaction(self, Reaction rxn): + """ Add a new reaction to this phase. """ + self.kinetics.addReaction(rxn._reaction) + def is_reversible(self, int i_reaction): """True if reaction `i_reaction` is reversible.""" self._check_reaction_index(i_reaction) diff --git a/interfaces/cython/cantera/test/test_kinetics.py b/interfaces/cython/cantera/test/test_kinetics.py index 81974570bc..99eb336351 100644 --- a/interfaces/cython/cantera/test/test_kinetics.py +++ b/interfaces/cython/cantera/test/test_kinetics.py @@ -223,6 +223,34 @@ def test_surface(self): k2 = surf2.kinetics_species_index(k) self.assertNear(rop1[k1], rop2[k2]) + def test_add_reaction(self): + gas1 = ct.Solution('h2o2.xml') + + S = ct.Species.listFromFile('h2o2.xml') + R = ct.Reaction.listFromFile('h2o2.xml') + gas2 = ct.Solution(thermo='IdealGas', kinetics='GasKinetics', + species=S, reactions=R[:5]) + + gas1.TPY = 800, 2*ct.one_atm, 'H2:0.3, O2:0.7, OH:2e-4, O:1e-3, H:5e-5' + gas2.TPY = gas1.TPY + + for r in R[5:]: + gas2.add_reaction(r) + + self.assertEqual(gas1.n_reactions, gas2.n_reactions) + + self.assertTrue((gas1.reactant_stoich_coeffs() == + gas2.reactant_stoich_coeffs()).all()) + self.assertTrue((gas1.product_stoich_coeffs() == + gas2.product_stoich_coeffs()).all()) + + self.assertArrayNear(gas1.delta_gibbs, + gas2.delta_gibbs) + self.assertArrayNear(gas1.reverse_rate_constants, + gas2.reverse_rate_constants) + self.assertArrayNear(gas1.net_production_rates, + gas2.net_production_rates) + class KineticsRepeatability(utilities.CanteraTest): """ diff --git a/interfaces/cython/cantera/test/test_thermo.py b/interfaces/cython/cantera/test/test_thermo.py index fdc57435ec..a9a5b04299 100644 --- a/interfaces/cython/cantera/test/test_thermo.py +++ b/interfaces/cython/cantera/test/test_thermo.py @@ -578,6 +578,24 @@ def test_uncopyable(self): with self.assertRaises(NotImplementedError): copy.copy(self.phase) + def test_add_species(self): + ref = ct.Solution('gri30.xml') + n_orig = self.phase.n_species + self.phase.add_species(ref.species('CO2')) + self.phase.add_species(ref.species('CO')) + + self.assertEqual(self.phase.n_species, n_orig + 2) + self.assertIn('CO2', self.phase.species_names) + self.assertIn('CO', self.phase.species_names) + + state = 400, 2e5, 'H2:0.7, CO2:0.2, CO:0.1' + ref.TPY = state + self.phase.TPY = state + self.assertNear(self.phase.enthalpy_mass, ref.enthalpy_mass) + self.assertNear(self.phase.entropy_mole, ref.entropy_mole) + self.assertArrayNear(ref[self.phase.species_names].partial_molar_cp, + self.phase.partial_molar_cp) + class TestThermo(utilities.CanteraTest): def setUp(self): diff --git a/interfaces/cython/cantera/thermo.pyx b/interfaces/cython/cantera/thermo.pyx index 179133b638..45eff40e1b 100644 --- a/interfaces/cython/cantera/thermo.pyx +++ b/interfaces/cython/cantera/thermo.pyx @@ -450,6 +450,17 @@ cdef class ThermoPhase(_SolutionBase): if self.kinetics: self.kinetics.invalidateCache() + def add_species(self, Species species): + """ + Add a new species to this phase. Missing elements will be added + automatically. + """ + self.thermo.addUndefinedElements() + self.thermo.addSpecies(species._species) + self.thermo.initThermo() + if self.kinetics: + self.kinetics.invalidateCache() + def n_atoms(self, species, element): """ Number of atoms of element *element* in species *species*. The element