constIsoSolidTransportI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
28 template<class thermo>
30 (
31  const thermo& t,
32  const scalar kappa
33 )
34 :
35  thermo(t),
36  kappa_(kappa)
37 {}
38 
39 
40 template<class thermo>
42 (
43  const word& name,
44  const constIsoSolidTransport& ct
45 )
46 :
47  thermo(name, ct),
48  kappa_(ct.kappa_)
49 {}
50 
51 
52 template<class Thermo>
55 (
56  const dictionary& dict
57 )
58 {
60  (
62  );
63 }
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 template<class thermo>
69 kappa(const scalar p, const scalar T) const
70 {
71  return kappa_;
72 }
73 
74 template<class thermo>
76 Kappa(const scalar p, const scalar T) const
77 {
78  return vector(kappa_, kappa_, kappa_);
79 }
80 
81 
82 template<class thermo>
84 mu(const scalar p, const scalar T) const
85 {
87  return scalar(0);
88 }
89 
90 
91 template<class thermo>
93 alphah(const scalar p, const scalar T) const
94 {
95  return kappa_/this->Cpv(p, T);
96 }
97 
98 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
99 
100 template<class thermo>
101 inline void Foam::constIsoSolidTransport<thermo>::operator=
102 (
104 )
105 {
106  thermo::operator=(ct);
107  kappa_ = ct.kappa_;
108 }
109 
110 
111 template<class thermo>
112 inline void Foam::constIsoSolidTransport<thermo>::operator+=
113 (
115 )
116 {
117  scalar molr1 = this->nMoles();
118  thermo::operator+=(ct);
119 
120  molr1 /= this->nMoles();
121  scalar molr2 = ct.nMoles()/this->nMoles();
122 
123  kappa_ = molr1*kappa_ + molr2*ct.kappa_;
124 }
125 
126 
127 template<class thermo>
128 inline void Foam::constIsoSolidTransport<thermo>::operator-=
129 (
131 )
132 {
133  scalar molr1 = this->nMoles();
134 
135  thermo::operator-=(ct);
136 
137  molr1 /= this->nMoles();
138  scalar molr2 = ct.nMoles()/this->nMoles();
139 
140  kappa_ = molr1*kappa_ - molr2*ct.kappa_;
141 }
142 
143 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
144 
145 
146 template<class thermo>
147 inline Foam::constIsoSolidTransport<thermo> Foam::operator*
148 (
149  const scalar s,
151 )
152 {
154  (
155  s*static_cast<const thermo&>(ct),
156  ct.kappa_
157  );
158 }
159 
160 
161 // ************************************************************************* //
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
scalar kappa(const scalar p, const scalar T) const
Isotropic thermal conductivity [W/mK].
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
psiReactionThermo & thermo
Definition: createFields.H:31
static autoPtr< constIsoSolidTransport > New(const dictionary &dict)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
A class for handling words, derived from string.
Definition: word.H:59
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
vector Kappa(const scalar p, const scalar T) const
Un-isotropic thermal conductivity [W/mK].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
volScalarField & p
Constant properties Transport package. Templated into a given thermodynamics package (needed for ther...
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366