1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::phaseModel
SourceFiles
phaseModel.C
\*---------------------------------------------------------------------------*/
#ifndef phaseModel_H
#define phaseModel_H
#include "dictionary.H"
#include "dictionaryEntry.H"
#include "dimensionedScalar.H"
#include "volFields.H"
#include "surfaceFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// Forward declarations
class diameterModel;
/*---------------------------------------------------------------------------*\
Class phaseModel Declaration
\*---------------------------------------------------------------------------*/
class phaseModel
:
public volScalarField
{
// Private data
//- Name of phase
word name_;
dictionary phaseDict_;
//- Kinematic viscosity
dimensionedScalar nu_;
//- Thermal conductivity
dimensionedScalar kappa_;
//- Heat capacity
dimensionedScalar Cp_;
//- Density
dimensionedScalar rho_;
//- Velocity
volVectorField U_;
//- Substantive derivative of the velocity
volVectorField DDtU_;
//- Volumetric flux of the phase
surfaceScalarField alphaPhi_;
//- Volumetric flux for the phase
autoPtr<surfaceScalarField> phiPtr_;
//- Diameter model
autoPtr<diameterModel> dPtr_;
public:
// Constructors
phaseModel
(
const word& phaseName,
const dictionary& phaseDict,
const fvMesh& mesh
);
//- Return clone
autoPtr<phaseModel> clone() const;
//- Return a pointer to a new phase created on freestore
// from Istream
class iNew
{
const fvMesh& mesh_;
public:
iNew
(
const fvMesh& mesh
)
:
mesh_(mesh)
{}
autoPtr<phaseModel> operator()(Istream& is) const
{
dictionaryEntry ent(dictionary::null, is);
return autoPtr<phaseModel>
(
new phaseModel(ent.keyword(), ent, mesh_)
);
}
};
//- Destructor
virtual ~phaseModel();
// Member Functions
const word& name() const
{
return name_;
}
const word& keyword() const
{
return name();
}
tmp<volScalarField> d() const;
const dimensionedScalar& nu() const
{
return nu_;
}
const dimensionedScalar& kappa() const
{
return kappa_;
}
const dimensionedScalar& Cp() const
{
return Cp_;
}
const dimensionedScalar& rho() const
{
return rho_;
}
const volVectorField& U() const
{
return U_;
}
volVectorField& U()
{
return U_;
}
const volVectorField& DDtU() const
{
return DDtU_;
}
volVectorField& DDtU()
{
return DDtU_;
}
const surfaceScalarField& phi() const
{
return phiPtr_();
}
surfaceScalarField& phi()
{
return phiPtr_();
}
const surfaceScalarField& alphaPhi() const
{
return alphaPhi_;
}
surfaceScalarField& alphaPhi()
{
return alphaPhi_;
}
//- Correct the phase properties
void correct();
//-Inherit read from volScalarField
using volScalarField::read;
//- Read base transportProperties dictionary
bool read(const dictionary& phaseDict);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
|