constAnIsoSolidTransportI.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Thermo>
32 (
33  const Thermo& t,
34  const vector kappa
35 )
36 :
37  Thermo(t),
38  kappa_(kappa)
39 {}
40 
41 
42 template<class Thermo>
44 (
45  const word& name,
46  const constAnIsoSolidTransport& ct
47 )
48 :
49  Thermo(name, ct),
50  kappa_(ct.kappa_)
51 {}
52 
53 
54 template<class Thermo>
57 (
58  const dictionary& dict
59 )
60 {
62 }
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
66 template<class Thermo>
68 kappa(const scalar p, const scalar T) const
69 {
70  return mag(kappa_);
71 }
72 
73 template<class Thermo>
75 Kappa(const scalar p, const scalar T) const
76 {
77  return kappa_;
78 }
79 
80 
81 template<class Thermo>
83 mu(const scalar p, const scalar T) const
84 {
86  return 0;
87 }
88 
89 
90 template<class Thermo>
92 alphah(const scalar p, const scalar T) const
93 {
94  return kappa_/this->Cp(p, T);
95 }
96 
97 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
98 
99 template<class Thermo>
101 (
103 )
104 {
105  scalar Y1 = this->Y();
106 
107  Y1 /= this->Y();
108  scalar Y2 = ct.Y()/this->Y();
109 
110  kappa_ = Y1*kappa_ + Y2*ct.kappa_;
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
115 
116 
117 template<class Thermo>
118 inline Foam::constAnIsoSolidTransport<Thermo> Foam::operator*
119 (
120  const scalar s,
121  const constAnIsoSolidTransport<Thermo>& ct
122 )
123 {
124  return constAnIsoSolidTransport<Thermo>
125  (
126  s*static_cast<const Thermo&>(ct),
127  ct.kappa_
128  );
129 }
130 
131 // ************************************************************************* //
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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].
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
static autoPtr< constAnIsoSolidTransport > New(const dictionary &dict)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
Vector< scalar > vector
Definition: vector.H:57
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & Cp
Definition: EEqn.H:7
scalar kappa(const scalar p, const scalar T) const
Isotropic thermal conductivity [W/mK].
const volScalarField & T
scalar mu(const scalar p, const scalar T) const
Dynamic viscosity [kg/ms].
PtrList< volScalarField > & Y
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
volScalarField & p
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Constant properties Transport package. Templated into a given Thermodynamics package (needed for ther...