constIsoSolidTransportI.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-2020 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 {
57  (
59  );
60 }
61 
62 
63 template<class Thermo>
66 (
67  const dictionary& dict
68 )
69 {
71  (
73  );
74 }
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
79 template<class Thermo>
81 kappa(const scalar p, const scalar T) const
82 {
83  return kappa_;
84 }
85 
86 template<class Thermo>
88 Kappa(const scalar p, const scalar T) const
89 {
90  return vector(kappa_, kappa_, kappa_);
91 }
92 
93 
94 template<class Thermo>
96 mu(const scalar p, const scalar T) const
97 {
99  return scalar(0);
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
104 
105 template<class Thermo>
106 inline void Foam::constIsoSolidTransport<Thermo>::operator+=
107 (
109 )
110 {
111  scalar Y1 = this->Y();
112  Thermo::operator+=(ct);
113 
114  Y1 /= this->Y();
115  scalar Y2 = ct.Y()/this->Y();
116 
117  kappa_ = Y1*kappa_ + Y2*ct.kappa_;
118 }
119 
120 
121 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
122 
123 
124 template<class Thermo>
125 inline Foam::constIsoSolidTransport<Thermo> Foam::operator*
126 (
127  const scalar s,
129 )
130 {
132  (
133  s*static_cast<const Thermo&>(ct),
134  ct.kappa_
135  );
136 }
137 
138 
139 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
vector Kappa(const scalar p, const scalar T) const
Un-isotropic thermal conductivity [W/m/K].
static autoPtr< constIsoSolidTransport > New(const dictionary &dict)
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
autoPtr< constIsoSolidTransport > clone() const
Construct and return a clone.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
PtrList< volScalarField > & Y
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
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:370
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].