GRINS-0.6.0
cantera_thermo.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2015 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 #include "grins_config.h"
27 
28 #ifdef GRINS_HAVE_CANTERA
29 
30 // This class
31 #include "grins/cantera_thermo.h"
32 
33 // GRINS
34 #include "grins/cantera_mixture.h"
35 #include "grins/cached_values.h"
36 
37 // libMesh
38 #include "libmesh/getpot.h"
39 
40 namespace GRINS
41 {
42 
44  : _cantera_mixture(mixture),
45  _cantera_gas(mixture.get_chemistry())
46  {
47  return;
48  }
49 
51  {
52  return;
53  }
54 
55  libMesh::Real CanteraThermodynamics::cp( const CachedValues& cache, unsigned int qp ) const
56  {
57  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
58  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
59  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
60 
61  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
62 
63  libMesh::Real cp = 0.0;
64 
65  {
68  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
69 
70  try
71  {
72  _cantera_gas.setState_TPY( T, P, &Y[0] );
73 
74  cp = _cantera_gas.cp_mass();
75  }
76  catch(Cantera::CanteraError)
77  {
78  Cantera::showErrors(std::cerr);
79  libmesh_error();
80  }
81 
82  }
83 
84  return cp;
85  }
86 
87  libMesh::Real CanteraThermodynamics::cv( const CachedValues& cache, unsigned int qp ) const
88  {
89  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
90  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
91  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
92 
93  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
94 
95  libMesh::Real cv = 0.0;
96 
97  {
100  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
101 
102  try
103  {
104  _cantera_gas.setState_TPY( T, P, &Y[0] );
105 
106  cv = _cantera_gas.cv_mass();
107  }
108  catch(Cantera::CanteraError)
109  {
110  Cantera::showErrors(std::cerr);
111  libmesh_error();
112  }
113 
114  }
115 
116  return cv;
117  }
118 
119  libMesh::Real CanteraThermodynamics::h( const CachedValues& cache, unsigned int qp,
120  unsigned int species ) const
121  {
122  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
123  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
124  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
125 
126  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
127 
128  std::vector<libMesh::Real> h_RT( Y.size(), 0.0 );
129 
130  {
133  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
134 
135  try
136  {
137  _cantera_gas.setState_TPY( T, P, &Y[0] );
138 
139  _cantera_gas.getEnthalpy_RT( &h_RT[0] );
140  }
141  catch(Cantera::CanteraError)
142  {
143  Cantera::showErrors(std::cerr);
144  libmesh_error();
145  }
146 
147  }
148 
149  return h_RT[species]*_cantera_mixture.R(species)*T;
150  }
151 
152  void CanteraThermodynamics::h( const CachedValues& cache, unsigned int qp,
153  std::vector<libMesh::Real>& h) const
154  {
155  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
156  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
157  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
158 
159  libmesh_assert_equal_to( Y.size(), h.size() );
160  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
161 
162  {
165  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
166 
167  try
168  {
169  _cantera_gas.setState_TPY( T, P, &Y[0] );
170 
171  _cantera_gas.getEnthalpy_RT( &h[0] );
172  }
173  catch(Cantera::CanteraError)
174  {
175  Cantera::showErrors(std::cerr);
176  libmesh_error();
177  }
178 
179  for( unsigned int s = 0; s < h.size(); s++ )
180  {
181  h[s] *= _cantera_mixture.R(s)*T;
182  }
183 
184  }
185 
186  return;
187  }
188 
189  libMesh::Real CanteraThermodynamics::h( const libMesh::Real& T, unsigned int species ) const
190  {
191  std::vector<libMesh::Real> h_RT( _cantera_gas.nSpecies(), 0.0 );
192 
193  try
194  {
195  _cantera_gas.setTemperature( T );
196 
197  _cantera_gas.getEnthalpy_RT( &h_RT[0] );
198  }
199  catch(Cantera::CanteraError)
200  {
201  Cantera::showErrors(std::cerr);
202  libmesh_error();
203  }
204 
205  return h_RT[species]*_cantera_mixture.R(species)*T;
206  }
207 
208 } // namespace GRINS
209 
210 #endif //GRINS_HAVE_CANTERA
libMesh::Real R(unsigned int species) const
libMesh::Real cp(const CachedValues &cache, unsigned int qp) const
GRINS namespace.
const std::vector< std::vector< libMesh::Number > > & get_cached_vector_values(unsigned int quantity) const
Cantera::IdealGasMix & _cantera_gas
libMesh::Real cv(const CachedValues &cache, unsigned int qp) const
Wrapper class for storing state for computing thermochemistry and transport properties using Cantera...
libMesh::Real h(const CachedValues &cache, unsigned int qp, unsigned int species) const
const std::vector< libMesh::Number > & get_cached_values(unsigned int quantity) const
Definition: cached_values.C:99
CanteraMixture & _cantera_mixture

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1