GRINS-0.6.0
cantera_transport.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
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 "cantera/transport.h"
44 #include "libmesh/restore_warnings.h"
45 
46 namespace GRINS
47 {
48 
50  : _cantera_gas( mixture.get_chemistry() ),
51  _cantera_transport( mixture.get_transport() )
52  {
53  return;
54  }
55 
57  {
58  return;
59  }
60 
61  libMesh::Real CanteraTransport::mu( const CachedValues& cache, unsigned int qp ) 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  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
68 
69  libMesh::Real mu = 0.0;
70 
71  {
72  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
73 
76  try
77  {
78  _cantera_gas.setState_TPY(T, P, &Y[0]);
79  mu = _cantera_transport.viscosity();
80  }
81  catch(Cantera::CanteraError)
82  {
83  Cantera::showErrors(std::cerr);
84  libmesh_error();
85  }
86 
87  }
88 
89  return mu;
90  }
91 
92  libMesh::Real CanteraTransport::k( const CachedValues& cache, unsigned int qp ) const
93  {
94  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
95  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
96  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
97 
98  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
99 
100  libMesh::Real k = 0.0;
101 
102  {
103  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
104 
107  try
108  {
109  _cantera_gas.setState_TPY(T, P, &Y[0]);
110  k = _cantera_transport.thermalConductivity();
111  }
112  catch(Cantera::CanteraError)
113  {
114  Cantera::showErrors(std::cerr);
115  libmesh_error();
116  }
117 
118  }
119 
120  return k;
121  }
122 
123  void CanteraTransport::D( const CachedValues& cache, unsigned int qp,
124  std::vector<libMesh::Real>& D ) const
125  {
126  const libMesh::Real T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
127  const libMesh::Real P = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];
128  const std::vector<libMesh::Real>& Y = cache.get_cached_vector_values(Cache::MASS_FRACTIONS)[qp];
129 
130  libmesh_assert_equal_to( Y.size(), D.size() );
131  libmesh_assert_equal_to( Y.size(), _cantera_gas.nSpecies() );
132 
133  {
134  libMesh::Threads::spin_mutex::scoped_lock lock(cantera_mutex);
135 
138  try
139  {
140  _cantera_gas.setState_TPY(T, P, &Y[0]);
141  _cantera_transport.getMixDiffCoeffsMass(&D[0]);
142  }
143  catch(Cantera::CanteraError)
144  {
145  Cantera::showErrors(std::cerr);
146  libmesh_error();
147  }
148 
149  }
150 
151  return;
152  }
153 
154 } // namespace GRINS
155 
156 #endif //GRINS_HAVE_CANTERA
void D(const CachedValues &cache, unsigned int qp, std::vector< libMesh::Real > &D) const
Cantera::IdealGasMix & _cantera_gas
Cantera::Transport & _cantera_transport
GRINS namespace.
const std::vector< std::vector< libMesh::Number > > & get_cached_vector_values(unsigned int quantity) const
libMesh::Real mu(const CachedValues &cache, unsigned int qp) const
libMesh::Real k(const CachedValues &cache, unsigned int qp) const
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