constAnIsoSolidTransportI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 vector kappa
33 )
34 :
35  Thermo(t),
36  kappa_(kappa)
37 {}
38 
39 
40 template<class Thermo>
42 (
43  const word& name,
44  const constAnIsoSolidTransport& 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 mag(kappa_);
72 }
73 
74 template<class Thermo>
76 Kappa(const scalar p, const scalar T) const
77 {
78  return 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->Cp(p, T);
96 }
97 
98 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
99 
100 template<class Thermo>
101 inline void Foam::constAnIsoSolidTransport<Thermo>::operator=
102 (
104 )
105 {
106  kappa_ = ct.kappa_;
107 }
108 
109 
110 template<class Thermo>
111 inline void Foam::constAnIsoSolidTransport<Thermo>::operator+=
112 (
114 )
115 {
116  scalar Y1 = this->Y();
117 
118  Y1 /= this->Y();
119  scalar Y2 = ct.Y()/this->Y();
120 
121  kappa_ = Y1*kappa_ + Y2*ct.kappa_;
122 }
123 
124 
125 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
126 
127 
128 template<class Thermo>
129 inline Foam::constAnIsoSolidTransport<Thermo> Foam::operator*
130 (
131  const scalar s,
133 )
134 {
136  (
137  s*static_cast<const Thermo&>(ct),
138  ct.kappa_
139  );
140 }
141 
142 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
vector alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
vector Kappa(const scalar p, const scalar T) const
Un-isotropic thermal conductivity [W/mK].
static autoPtr< constAnIsoSolidTransport > New(const dictionary &dict)
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].
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))
A class for handling words, derived from string.
Definition: word.H:59
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar kappa(const scalar p, const scalar T) const
Isotropic thermal conductivity [W/mK].
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField & p
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
Constant properties Transport package. Templated into a given Thermodynamics package (needed for ther...