192 lines
4.3 KiB
C++
192 lines
4.3 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2013 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 "dimensionedScalar.H"
|
|
#include "volFields.H"
|
|
#include "surfaceFields.H"
|
|
#include "transportModel.H"
|
|
#include "rhoThermo.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
|
|
// Forward declarations
|
|
class twoPhaseSystem;
|
|
class diameterModel;
|
|
|
|
|
|
/*---------------------------------------------------------------------------*\
|
|
Class phaseModel Declaration
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
class phaseModel
|
|
:
|
|
public volScalarField,
|
|
public transportModel
|
|
{
|
|
// Private data
|
|
|
|
//- Reference to the twoPhaseSystem to which this phase belongs
|
|
const twoPhaseSystem& fluid_;
|
|
|
|
//- Name of phase
|
|
word name_;
|
|
|
|
dictionary phaseDict_;
|
|
|
|
//- Thermophysical properties
|
|
autoPtr<rhoThermo> thermo_;
|
|
|
|
//- Velocity
|
|
volVectorField U_;
|
|
|
|
//- Fluxes
|
|
autoPtr<surfaceScalarField> phiPtr_;
|
|
|
|
//- Diameter model
|
|
autoPtr<diameterModel> dPtr_;
|
|
|
|
|
|
public:
|
|
|
|
// Constructors
|
|
|
|
phaseModel
|
|
(
|
|
const twoPhaseSystem& fluid,
|
|
const dictionary& phaseProperties,
|
|
const word& phaseName
|
|
);
|
|
|
|
|
|
//- Destructor
|
|
virtual ~phaseModel();
|
|
|
|
|
|
// Member Functions
|
|
|
|
//- Return the twoPhaseSystem to which this phase belongs
|
|
const twoPhaseSystem& fluid() const
|
|
{
|
|
return fluid_;
|
|
}
|
|
|
|
const word& name() const
|
|
{
|
|
return name_;
|
|
}
|
|
|
|
tmp<volScalarField> d() const;
|
|
|
|
tmp<volScalarField> nu() const
|
|
{
|
|
return thermo_->nu();
|
|
}
|
|
|
|
//- Return the laminar viscosity for patch
|
|
tmp<scalarField> nu(const label patchi) const
|
|
{
|
|
return thermo_->nu(patchi);
|
|
}
|
|
|
|
tmp<volScalarField> kappa() const
|
|
{
|
|
return thermo_->kappa();
|
|
}
|
|
|
|
tmp<volScalarField> Cp() const
|
|
{
|
|
return thermo_->Cp();
|
|
}
|
|
|
|
const volScalarField& rho() const
|
|
{
|
|
return thermo_->rho();
|
|
}
|
|
|
|
const rhoThermo& thermo() const
|
|
{
|
|
return thermo_();
|
|
}
|
|
|
|
rhoThermo& thermo()
|
|
{
|
|
return thermo_();
|
|
}
|
|
|
|
const volVectorField& U() const
|
|
{
|
|
return U_;
|
|
}
|
|
|
|
volVectorField& U()
|
|
{
|
|
return U_;
|
|
}
|
|
|
|
const surfaceScalarField& phi() const
|
|
{
|
|
return phiPtr_();
|
|
}
|
|
|
|
surfaceScalarField& phi()
|
|
{
|
|
return phiPtr_();
|
|
}
|
|
|
|
//- Dummy correct
|
|
void correct()
|
|
{}
|
|
|
|
//- Dummy read
|
|
bool read()
|
|
{
|
|
return true;
|
|
}
|
|
};
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
} // End namespace Foam
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
#endif
|
|
|
|
// ************************************************************************* //
|