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-2019 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 {
57  (
59  );
60 }
61 
62 
63 template<class Thermo>
66 (
67  const dictionary& dict
68 )
69 {
71  (
73  );
74 }
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
78 template<class Thermo>
80 kappa(const scalar p, const scalar T) const
81 {
82  return mag(kappa_);
83 }
84 
85 template<class Thermo>
87 Kappa(const scalar p, const scalar T) const
88 {
89  return kappa_;
90 }
91 
92 
93 template<class Thermo>
95 mu(const scalar p, const scalar T) const
96 {
98  return scalar(0);
99 }
100 
101 
102 template<class Thermo>
104 alphah(const scalar p, const scalar T) const
105 {
106  return kappa_/this->Cp(p, T);
107 }
108 
109 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
110 
111 template<class Thermo>
112 inline void Foam::constAnIsoSolidTransport<Thermo>::operator+=
113 (
115 )
116 {
117  scalar Y1 = this->Y();
118 
119  Y1 /= this->Y();
120  scalar Y2 = ct.Y()/this->Y();
121 
122  kappa_ = Y1*kappa_ + Y2*ct.kappa_;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
127 
128 
129 template<class Thermo>
130 inline Foam::constAnIsoSolidTransport<Thermo> Foam::operator*
131 (
132  const scalar s,
134 )
135 {
137  (
138  s*static_cast<const Thermo&>(ct),
139  ct.kappa_
140  );
141 }
142 
143 // ************************************************************************* //
dictionary dict
const dimensionedScalar & kappa
Coulomb constant: default SI units: [N.m2/C2].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
vector alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/m/s].
vector Kappa(const scalar p, const scalar T) const
Un-isotropic thermal conductivity [W/m/K].
static autoPtr< constAnIsoSolidTransport > New(const dictionary &dict)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
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/m/K].
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
autoPtr< constAnIsoSolidTransport > clone() const
Construct and return a clone.
PtrList< volScalarField > & Y
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
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...