GRINS-0.6.0
cantera_kinetics.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_kinetics.h"
32 
33 // GRINS
34 #include "grins/cached_values.h"
35 #include "grins/cantera_mixture.h"
36 
37 // libMesh
38 #include "libmesh/getpot.h"
39 
40 // Cantera (with compiler warnings disabled)
41 #include "libmesh/ignore_warnings.h"
42 #include "cantera/IdealGasMix.h"
43 #include "libmesh/restore_warnings.h"
44 
45 namespace GRINS
46 {
47 
49  : _cantera_gas( mixture.get_chemistry() )
50  {
51  return;
52  }
53 
55  {
56  return;
57  }
58 
60  unsigned int qp,
61  std::vector<libMesh::Real>& omega_dot ) const
62  {
63  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
64  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
65  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
66 
67  this->omega_dot_TPY( T, P, Y, omega_dot );
68 
69  return;
70  }
71 
72  void CanteraKinetics::omega_dot_TPY( const libMesh::Real T, const libMesh::Real P,
73  const std::vector<libMesh::Real>& mass_fractions,
74  std::vector<libMesh::Real>& omega_dot ) const
75  {
76  libmesh_assert_equal_to( mass_fractions.size(), omega_dot.size() );
77  libmesh_assert_equal_to( mass_fractions.size(), _cantera_gas.nSpecies() );
78  libmesh_assert_greater(T,0.0);
79  libmesh_assert_greater(P,0.0);
80 
81  {
84  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
85 
86  try
87  {
88  _cantera_gas.setState_TPY(T, P, &mass_fractions[0]);
89  _cantera_gas.getNetProductionRates(&omega_dot[0]);
90  }
91  catch(Cantera::CanteraError)
92  {
93  Cantera::showErrors(std::cerr);
94  libmesh_error();
95  }
96  }
97 
98 #ifdef DEBUG
99  for( unsigned int s = 0; s < omega_dot.size(); s++ )
100  {
101  if( libMesh::libmesh_isnan(omega_dot[s]) )
102  {
103  std::cout << "T = " << T << std::endl
104  << "P = " << P << std::endl;
105  for( unsigned int s = 0; s < omega_dot.size(); s++ )
106  {
107  std::cout << "Y[" << s << "] = " << mass_fractions[s] << std::endl;
108  }
109  for( unsigned int s = 0; s < omega_dot.size(); s++ )
110  {
111  std::cout << "omega_dot[" << s << "] = " << omega_dot[s] << std::endl;
112  }
113 
114  libmesh_error();
115  }
116  }
117 #endif
118 
119  for( unsigned int s = 0; s < omega_dot.size(); s++ )
120  {
121  // convert [kmol/m^3-s] to [kg/m^3-s]
122  omega_dot[s] *= this->_cantera_gas.molecularWeight(s);
123  }
124 
125  return;
126  }
127 
128  void CanteraKinetics::omega_dot_TRY( const libMesh::Real& T, const libMesh::Real rho,
129  const std::vector<libMesh::Real>& mass_fractions,
130  std::vector<libMesh::Real>& omega_dot ) const
131  {
132  libmesh_assert_equal_to( mass_fractions.size(), omega_dot.size() );
133  libmesh_assert_equal_to( mass_fractions.size(), _cantera_gas.nSpecies() );
134  libmesh_assert_greater(T,0.0);
135  libmesh_assert_greater(rho,0.0);
136 
137  {
140  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
141 
142  try
143  {
144  _cantera_gas.setState_TRY(T, rho, &mass_fractions[0]);
145  _cantera_gas.getNetProductionRates(&omega_dot[0]);
146  }
147  catch(Cantera::CanteraError)
148  {
149  Cantera::showErrors(std::cerr);
150  libmesh_error();
151  }
152  }
153 
154  for( unsigned int s = 0; s < omega_dot.size(); s++ )
155  {
156  // convert [kmol/m^3-s] to [kg/m^3-s]
157  omega_dot[s] *= this->_cantera_gas.molecularWeight(s);
158  }
159 
160  return;
161  }
162 
163 } // namespace GRINS
164 
165 #endif //GRINS_HAVE_CANTERA
void omega_dot_TRY(const libMesh::Real &T, const libMesh::Real rho, const std::vector< libMesh::Real > &mass_fractions, std::vector< libMesh::Real > &omega_dot) const
void omega_dot(const CachedValues &cache, unsigned int qp, std::vector< libMesh::Real > &omega_dot) const
void omega_dot_TPY(const libMesh::Real T, const libMesh::Real P, const std::vector< libMesh::Real > &mass_fractions, std::vector< libMesh::Real > &omega_dot) const
GRINS namespace.
const std::vector< std::vector< libMesh::Number > > & get_cached_vector_values(unsigned int quantity) const
Cantera::IdealGasMix & _cantera_gas
Wrapper class for storing state for computing thermochemistry and transport properties using Cantera...
const std::vector< libMesh::Number > & get_cached_values(unsigned int quantity) const
Definition: cached_values.C:99

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