mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
137 lines
3.5 KiB
C
137 lines
3.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 "addToRunTimeSelectionTable.H"
|
|
#include "powerLaw.H"
|
|
#include "geometricOneField.H"
|
|
#include "fvMatrices.H"
|
|
|
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
namespace porosityModels
|
|
{
|
|
defineTypeNameAndDebug(powerLaw, 0);
|
|
addToRunTimeSelectionTable(porosityModel, powerLaw, mesh);
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::porosityModels::powerLaw::powerLaw
|
|
(
|
|
const word& name,
|
|
const word& modelType,
|
|
const fvMesh& mesh,
|
|
const dictionary& dict,
|
|
const word& cellZoneName
|
|
)
|
|
:
|
|
porosityModel(name, modelType, mesh, dict, cellZoneName),
|
|
C0_(readScalar(coeffs_.lookup("C0"))),
|
|
C1_(readScalar(coeffs_.lookup("C1"))),
|
|
rhoName_(coeffs_.lookupOrDefault<word>("rho", "rho"))
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
|
|
Foam::porosityModels::powerLaw::~powerLaw()
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
void Foam::porosityModels::powerLaw::correct
|
|
(
|
|
fvVectorMatrix& UEqn
|
|
) const
|
|
{
|
|
const vectorField& U = UEqn.psi();
|
|
const scalarField& V = mesh_.V();
|
|
scalarField& Udiag = UEqn.diag();
|
|
|
|
if (UEqn.dimensions() == dimForce)
|
|
{
|
|
const volScalarField& rho =
|
|
mesh_.lookupObject<volScalarField>(rhoName_);
|
|
|
|
apply(Udiag, V, rho, U);
|
|
}
|
|
else
|
|
{
|
|
apply(Udiag, V, geometricOneField(), U);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::porosityModels::powerLaw::correct
|
|
(
|
|
fvVectorMatrix& UEqn,
|
|
const volScalarField& rho,
|
|
const volScalarField& mu
|
|
) const
|
|
{
|
|
const vectorField& U = UEqn.psi();
|
|
const scalarField& V = mesh_.V();
|
|
scalarField& Udiag = UEqn.diag();
|
|
|
|
apply(Udiag, V, rho, U);
|
|
}
|
|
|
|
|
|
void Foam::porosityModels::powerLaw::correct
|
|
(
|
|
const fvVectorMatrix& UEqn,
|
|
volTensorField& AU
|
|
) const
|
|
{
|
|
const vectorField& U = UEqn.psi();
|
|
|
|
if (UEqn.dimensions() == dimForce)
|
|
{
|
|
const volScalarField& rho =
|
|
mesh_.lookupObject<volScalarField>(rhoName_);
|
|
|
|
apply(AU, rho, U);
|
|
}
|
|
else
|
|
{
|
|
apply(AU, geometricOneField(), U);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::porosityModels::powerLaw::writeData(Ostream& os) const
|
|
{
|
|
os << indent << name_ << endl;
|
|
dict_.write(os);
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|