Skip to content

Commit

Permalink
[Python] Prevent adding species to in-use Thermo objects
Browse files Browse the repository at this point in the history
  • Loading branch information
speth committed Sep 10, 2016
1 parent 178b452 commit e49f843
Show file tree
Hide file tree
Showing 6 changed files with 54 additions and 0 deletions.
4 changes: 4 additions & 0 deletions interfaces/cython/cantera/_cantera.pxd
Original file line number Diff line number Diff line change
Expand Up @@ -874,6 +874,7 @@ cdef class ThermoPhase(_SolutionBase):
cpdef int species_index(self, species) except *
cdef np.ndarray _getArray1(self, thermoMethod1d method)
cdef void _setArray1(self, thermoMethod1d method, values) except *
cdef public object _references

cdef class InterfacePhase(ThermoPhase):
cdef CxxSurfPhase* surf
Expand Down Expand Up @@ -906,6 +907,7 @@ cdef class DustyGasTransport(Transport):
cdef class Mixture:
cdef CxxMultiPhase* mix
cdef list _phases
cdef object _weakref_proxy
cpdef int element_index(self, element) except *

cdef class Func1:
Expand All @@ -919,6 +921,7 @@ cdef class ReactorBase:
cdef list _inlets
cdef list _outlets
cdef list _walls
cdef object _weakref_proxy

cdef class Reactor(ReactorBase):
cdef CxxReactor* reactor
Expand Down Expand Up @@ -982,6 +985,7 @@ cdef class ReactorNet:
cdef class Domain1D:
cdef CxxDomain1D* domain
cdef _SolutionBase gas
cdef object _weakref_proxy
cdef public pybool have_user_tolerances

cdef class Boundary1D(Domain1D):
Expand Down
6 changes: 6 additions & 0 deletions interfaces/cython/cantera/mixture.pyx
Original file line number Diff line number Diff line change
@@ -1,5 +1,9 @@
import warnings

# Need a pure-python class to store weakrefs to
class _WeakrefProxy(object):
pass

cdef class Mixture:
"""
Expand Down Expand Up @@ -32,13 +36,15 @@ cdef class Mixture:
def __cinit__(self, phases):
self.mix = new CxxMultiPhase()
self._phases = []
self._weakref_proxy = _WeakrefProxy()

cdef _SolutionBase phase
if isinstance(phases[0], _SolutionBase):
# Assign default composition to the list of phases
phases = [(p, 1 if i == 0 else 0) for i,p in enumerate(phases)]

for phase,moles in phases:
phase._references[self._weakref_proxy] = True
self.mix.addPhase(phase.thermo, moles)
self._phases.append(phase)

Expand Down
8 changes: 8 additions & 0 deletions interfaces/cython/cantera/onedim.pyx
Original file line number Diff line number Diff line change
@@ -1,18 +1,24 @@
import interrupts

# Need a pure-python class to store weakrefs to
class _WeakrefProxy(object):
pass

cdef class Domain1D:
def __cinit__(self, *args, **kwargs):
self.domain = NULL

# The signature of this function causes warnings for Sphinx documentation
def __init__(self, _SolutionBase phase, *args, name=None, **kwargs):
self._weakref_proxy = _WeakrefProxy()
if self.domain is NULL:
raise TypeError("Can't instantiate abstract class Domain1D.")

if name is not None:
self.name = name

self.gas = phase
self.gas._references[self._weakref_proxy] = True
self.have_user_tolerances = False

property index:
Expand Down Expand Up @@ -407,6 +413,8 @@ cdef class _FlowBase(Domain1D):
"""
Set the `Solution` object used for calculating transport properties.
"""
self._weakref_proxy = _WeakrefProxy()
self.gas._references[self._weakref_proxy] = True
self.gas = phase
self.flow.setTransport(deref(self.gas.transport))

Expand Down
6 changes: 6 additions & 0 deletions interfaces/cython/cantera/reactor.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,10 @@ import numbers as _numbers

_reactor_counts = _defaultdict(int)

# Need a pure-python class to store weakrefs to
class _WeakrefProxy(object):
pass

cdef class ReactorBase:
"""
Common base class for reactors and reservoirs.
Expand All @@ -13,6 +17,7 @@ cdef class ReactorBase:

# The signature of this function causes warnings for Sphinx documentation
def __init__(self, ThermoPhase contents=None, name=None, *, volume=None):
self._weakref_proxy = _WeakrefProxy()
self._inlets = []
self._outlets = []
self._walls = []
Expand All @@ -38,6 +43,7 @@ cdef class ReactorBase:
properties and kinetic rates for this reactor.
"""
self._thermo = solution
self._thermo._references[self._weakref_proxy] = True
self.rbase.setThermoMgr(deref(solution.thermo))

property name:
Expand Down
25 changes: 25 additions & 0 deletions interfaces/cython/cantera/test/test_thermo.py
Original file line number Diff line number Diff line change
@@ -1,5 +1,6 @@
from .utilities import unittest
import numpy as np
import gc

import cantera as ct
from . import utilities
Expand Down Expand Up @@ -596,6 +597,30 @@ def test_add_species(self):
self.assertArrayNear(ref[self.phase.species_names].partial_molar_cp,
self.phase.partial_molar_cp)

def test_add_species_disabled(self):
ref = ct.Solution('gri30.xml')

reactor = ct.IdealGasReactor(self.phase)
with self.assertRaises(RuntimeError):
self.phase.add_species(ref.species('CH4'))
del reactor
gc.collect()
self.phase.add_species(ref.species('CH4'))

flame = ct.FreeFlame(self.phase, width=0.1)
with self.assertRaises(RuntimeError):
self.phase.add_species(ref.species('CO'))
del flame
gc.collect()
self.phase.add_species(ref.species('CO'))

mix = ct.Mixture([(self.phase, 2.0)])
with self.assertRaises(RuntimeError):
self.phase.add_species(ref.species('CH2O'))
del mix
gc.collect()
self.phase.add_species(ref.species('CH2O'))


class TestThermo(utilities.CanteraTest):
def setUp(self):
Expand Down
5 changes: 5 additions & 0 deletions interfaces/cython/cantera/thermo.pyx
Original file line number Diff line number Diff line change
@@ -1,4 +1,5 @@
import warnings
import weakref

cdef enum ThermoBasis:
mass_basis = 0
Expand Down Expand Up @@ -223,6 +224,7 @@ cdef class ThermoPhase(_SolutionBase):
super().__init__(*args, **kwargs)
if 'source' not in kwargs:
self.thermo_basis = mass_basis
self._references = weakref.WeakKeyDictionary()

def report(self, show_thermo=True, float threshold=1e-14):
"""
Expand Down Expand Up @@ -455,6 +457,9 @@ cdef class ThermoPhase(_SolutionBase):
Add a new species to this phase. Missing elements will be added
automatically.
"""
if self._references:
raise RuntimeError('Cannot add species to ThermoPhase object if it'
' is linked to a Reactor, Domain1D (flame), or Mixture object.')
self.thermo.addUndefinedElements()
self.thermo.addSpecies(species._species)
self.thermo.initThermo()
Expand Down

0 comments on commit e49f843

Please sign in to comment.