Fix tutorials: typo in tutorials/viscoelastic/viscoelasticFluidFoam/S-MDCPP/constant...
[OpenFOAM-1.6-ext.git] / src / thermophysicalModels / specie / transport / const / constTransportI.H
bloba1c79b67c779c43fb47f4d5f28f3ff4f36449dc6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 namespace Foam
32 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
34 // Construct from components
35 template<class thermo>
36 inline constTransport<thermo>::constTransport
38     const thermo& t,
39     const scalar mu,
40     const scalar Pr
43     thermo(t),
44     Mu(mu),
45     rPr(1.0/Pr)
49 // Construct as named copy
50 template<class thermo>
51 inline constTransport<thermo>::constTransport
53     const word& name,
54     const constTransport& ct
57     thermo(name, ct),
58     Mu(ct.Mu),
59     rPr(ct.rPr)
63 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
65 // Dynamic viscosity [kg/ms]
66 template<class thermo>
67 inline scalar constTransport<thermo>::mu(const scalar) const
69     return Mu;
73 // Thermal conductivity [W/mK]
74 template<class thermo>
75 inline scalar constTransport<thermo>::kappa(const scalar T) const
77     return this->Cp(T)*mu(T)*rPr;
81 // Thermal diffusivity for enthalpy [kg/ms]
82 template<class thermo>
83 inline scalar constTransport<thermo>::alpha(const scalar T) const
85     scalar Cp_ = this->Cp(T);
87     scalar deltaT = T - specie::Tstd;
88     scalar CpBar =
89         (deltaT*(this->H(T) - this->H(specie::Tstd)) + Cp_)/(sqr(deltaT) + 1);
91     return Cp_*mu(T)*rPr/CpBar;
95 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
97 template<class thermo>
98 inline constTransport<thermo>& constTransport<thermo>::operator=
100     const constTransport<thermo>& ct
103     thermo::operator=(ct);
105     Mu = ct.Mu;
106     rPr = ct.rPr;
108     return *this;
112 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
114 template<class thermo>
115 inline constTransport<thermo> operator+
117     const constTransport<thermo>& ct1,
118     const constTransport<thermo>& ct2
121     thermo t
122     (
123         static_cast<const thermo&>(ct1) + static_cast<const thermo&>(ct2)
124     );
126     scalar molr1 = ct1.nMoles()/t.nMoles();
127     scalar molr2 = ct2.nMoles()/t.nMoles();
129     return constTransport<thermo>
130     (
131         t,
132         molr1*ct1.Mu + molr2*ct2.Mu,
133         molr1*ct1.rPr + molr2*ct2.rPr
134     );
138 template<class thermo>
139 inline constTransport<thermo> operator-
141     const constTransport<thermo>& ct1,
142     const constTransport<thermo>& ct2
145     thermo t
146     (
147         static_cast<const thermo&>(ct1) - static_cast<const thermo&>(ct2)
148     );
150     scalar molr1 = ct1.nMoles()/t.nMoles();
151     scalar molr2 = ct2.nMoles()/t.nMoles();
153     return constTransport<thermo>
154     (
155         t,
156         molr1*ct1.Mu - molr2*ct2.Mu,
157         molr1*ct1.rPr - molr2*ct2.rPr
158     );
162 template<class thermo>
163 inline constTransport<thermo> operator*
165     const scalar s,
166     const constTransport<thermo>& ct
169     return constTransport<thermo>
170     (
171         s*static_cast<const thermo&>(ct),
172         ct.Mu,
173         ct.rPr
174     );
178 template<class thermo>
179 inline constTransport<thermo> operator==
181     const constTransport<thermo>& ct1,
182     const constTransport<thermo>& ct2
185     return ct2 - ct1;
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 } // End namespace Foam
193 // ************************************************************************* //