Files
openfoam/applications/solvers/multiphase/twoPhaseEulerFoam/phaseModel/phaseModel/phaseModel.H
2009-02-05 15:28:32 +00:00

155 lines
3.6 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Class
Foam::phaseModel
SourceFiles
phaseModel.C
\*---------------------------------------------------------------------------*/
#ifndef phaseModel_H
#define phaseModel_H
#include "dictionary.H"
#include "dimensionedScalar.H"
#include "volFields.H"
#include "surfaceFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class phaseModel Declaration
\*---------------------------------------------------------------------------*/
class phaseModel
{
// Private data
dictionary dict_;
//- Name of phase
word name_;
//- Characteristic diameter of phase
dimensionedScalar d_;
//- kinematic viscosity
dimensionedScalar nu_;
//- density
dimensionedScalar rho_;
//- Velocity
volVectorField U_;
//- Fluxes
autoPtr<surfaceScalarField> phiPtr_;
public:
// Constructors
phaseModel
(
const fvMesh& mesh,
const dictionary& transportProperties,
const word& phaseName
);
// Selectors
//- Return a reference to the selected turbulence model
static autoPtr<phaseModel> New
(
const fvMesh& mesh,
const dictionary& transportProperties,
const word& phaseName
);
//- Destructor
virtual ~phaseModel();
// Member Functions
const word& name() const
{
return name_;
}
const dimensionedScalar& d() const
{
return d_;
}
const dimensionedScalar& nu() const
{
return nu_;
}
const dimensionedScalar& rho() const
{
return rho_;
}
const volVectorField& U() const
{
return U_;
}
volVectorField& U()
{
return U_;
}
const surfaceScalarField& phi() const
{
return phiPtr_();
}
surfaceScalarField& phi()
{
return phiPtr_();
}
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //