Files
openfoam/src/thermophysicalModels/basic/basicThermo/basicThermo.C

237 lines
5.5 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 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/>.
\*---------------------------------------------------------------------------*/
#include "basicThermo.H"
#include "fvMesh.H"
#include "HashTable.H"
#include "zeroGradientFvPatchFields.H"
#include "fixedEnergyFvPatchScalarField.H"
#include "gradientEnergyFvPatchScalarField.H"
#include "mixedEnergyFvPatchScalarField.H"
#include "temperatureJumpFvPatchScalarField.H"
#include "energyJumpFvPatchScalarField.H"
/* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
namespace Foam
{
defineTypeNameAndDebug(basicThermo, 0);
defineRunTimeSelectionTable(basicThermo, fvMesh);
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
Foam::wordList Foam::basicThermo::heBoundaryTypes()
{
const volScalarField::GeometricBoundaryField& tbf = T_.boundaryField();
wordList hbt = tbf.types();
forAll(tbf, patchi)
{
if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
{
hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
}
else if
(
isA<zeroGradientFvPatchScalarField>(tbf[patchi])
|| isA<fixedGradientFvPatchScalarField>(tbf[patchi])
)
{
hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
}
else if(isA<mixedFvPatchScalarField>(tbf[patchi]))
{
hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
}
else if (isA<temperatureJumpFvPatchScalarField>(tbf[patchi]))
{
hbt[patchi] = energyJumpFvPatchScalarField::typeName;
}
}
return hbt;
}
void Foam::basicThermo::heBoundaryCorrection(volScalarField& h)
{
volScalarField::GeometricBoundaryField& hbf = h.boundaryField();
forAll(hbf, patchi)
{
if (isA<gradientEnergyFvPatchScalarField>(hbf[patchi]))
{
refCast<gradientEnergyFvPatchScalarField>(hbf[patchi]).gradient()
= hbf[patchi].fvPatchField::snGrad();
}
else if (isA<mixedEnergyFvPatchScalarField>(hbf[patchi]))
{
refCast<mixedEnergyFvPatchScalarField>(hbf[patchi]).refGrad()
= hbf[patchi].fvPatchField::snGrad();
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::basicThermo::basicThermo(const fvMesh& mesh)
:
IOdictionary
(
IOobject
(
"thermophysicalProperties",
mesh.time().constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::NO_WRITE
)
),
p_
(
IOobject
(
"p",
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
),
psi_
(
IOobject
(
"psi",
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionSet(0, -2, 2, 0, 0)
),
T_
(
IOobject
(
"T",
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
),
mu_
(
IOobject
(
"mu",
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionSet(1, -1, -1, 0, 0)
),
alpha_
(
IOobject
(
"alpha",
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionSet(1, -1, -1, 0, 0)
)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::basicThermo::~basicThermo()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::volScalarField& Foam::basicThermo::p()
{
return p_;
}
const Foam::volScalarField& Foam::basicThermo::p() const
{
return p_;
}
const Foam::volScalarField& Foam::basicThermo::psi() const
{
return psi_;
}
const Foam::volScalarField& Foam::basicThermo::T() const
{
return T_;
}
const Foam::volScalarField& Foam::basicThermo::mu() const
{
return mu_;
}
const Foam::volScalarField& Foam::basicThermo::alpha() const
{
return alpha_;
}
bool Foam::basicThermo::read()
{
return regIOobject::read();
}
// ************************************************************************* //