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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
30template<class thermo>
32(
33 const thermo& t,
34 const scalar kappa
35)
36:
37 thermo(t),
38 kappa_(kappa)
39{}
40
41
42template<class thermo>
44(
45 const word& name,
47)
48:
49 thermo(name, ct),
50 kappa_(ct.kappa_)
51{}
52
53
54template<class Thermo>
57(
58 const dictionary& dict
59)
60{
62}
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
66template<class thermo>
68kappa(const scalar p, const scalar T) const
69{
70 return kappa_;
71}
72
73template<class thermo>
75Kappa(const scalar p, const scalar T) const
76{
77 return vector(kappa_, kappa_, kappa_);
78}
79
80
81template<class thermo>
83mu(const scalar p, const scalar T) const
84{
86 return 0;
87}
88
89
90template<class thermo>
92alphah(const scalar p, const scalar T) const
93{
94 return kappa_/this->Cp(p, T);
95}
96
97// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
98
99template<class thermo>
101(
103)
104{
105 scalar Y1 = this->Y();
106 thermo::operator+=(ct);
107
108 Y1 /= this->Y();
109 scalar Y2 = ct.Y()/this->Y();
110
111 kappa_ = Y1*kappa_ + Y2*ct.kappa_;
112}
113
114
115// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
116
117template<class thermo>
118inline Foam::constIsoSolidTransport<thermo> Foam::operator*
119(
120 const scalar s,
122)
123{
125 (
126 s*static_cast<const thermo&>(ct),
127 ct.kappa_
128 );
129}
130
131
132// ************************************************************************* //
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Constant properties Transport package. Templated into a given thermodynamics package (needed for ther...
scalar alphah(const scalar p, const scalar T) const
Thermal diffusivity of enthalpy [kg/ms].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & T
const volScalarField & Cp
Definition: EEqn.H:7
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
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))
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar mu
Atomic mass unit.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Vector< scalar > vector
Definition: vector.H:61
dictionary dict