mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
428 lines
8.1 KiB
C
428 lines
8.1 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2012-2016 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 "MRFZoneList.H"
|
|
#include "volFields.H"
|
|
#include "fixedValueFvsPatchFields.H"
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::MRFZoneList::MRFZoneList
|
|
(
|
|
const fvMesh& mesh,
|
|
const dictionary& dict
|
|
)
|
|
:
|
|
PtrList<MRFZone>(),
|
|
mesh_(mesh)
|
|
{
|
|
reset(dict);
|
|
|
|
active(true);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
|
|
Foam::MRFZoneList::~MRFZoneList()
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
bool Foam::MRFZoneList::active(const bool warn) const
|
|
{
|
|
bool a = false;
|
|
forAll(*this, i)
|
|
{
|
|
a = a || this->operator[](i).active();
|
|
}
|
|
|
|
if (warn && this->size() && !a)
|
|
{
|
|
Info<< " No MRF zones active" << endl;
|
|
}
|
|
|
|
return a;
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::reset(const dictionary& dict)
|
|
{
|
|
label count = 0;
|
|
forAllConstIter(dictionary, dict, iter)
|
|
{
|
|
if (iter().isDict())
|
|
{
|
|
count++;
|
|
}
|
|
}
|
|
|
|
this->setSize(count);
|
|
label i = 0;
|
|
forAllConstIter(dictionary, dict, iter)
|
|
{
|
|
if (iter().isDict())
|
|
{
|
|
const word& name = iter().keyword();
|
|
const dictionary& modelDict = iter().dict();
|
|
|
|
Info<< " creating MRF zone: " << name << endl;
|
|
|
|
this->set
|
|
(
|
|
i++,
|
|
new MRFZone(name, mesh_, modelDict)
|
|
);
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
bool Foam::MRFZoneList::read(const dictionary& dict)
|
|
{
|
|
bool allOk = true;
|
|
forAll(*this, i)
|
|
{
|
|
MRFZone& pm = this->operator[](i);
|
|
bool ok = pm.read(dict.subDict(pm.name()));
|
|
allOk = (allOk && ok);
|
|
}
|
|
return allOk;
|
|
}
|
|
|
|
|
|
bool Foam::MRFZoneList::writeData(Ostream& os) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
os << nl;
|
|
this->operator[](i).writeData(os);
|
|
}
|
|
|
|
return os.good();
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::addAcceleration
|
|
(
|
|
const volVectorField& U,
|
|
volVectorField& ddtU
|
|
) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).addCoriolis(U, ddtU);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::addAcceleration(fvVectorMatrix& UEqn) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).addCoriolis(UEqn);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::addAcceleration
|
|
(
|
|
const volScalarField& rho,
|
|
fvVectorMatrix& UEqn
|
|
) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).addCoriolis(rho, UEqn);
|
|
}
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volVectorField> Foam::MRFZoneList::DDt
|
|
(
|
|
const volVectorField& U
|
|
) const
|
|
{
|
|
tmp<volVectorField> tacceleration
|
|
(
|
|
new volVectorField
|
|
(
|
|
IOobject
|
|
(
|
|
"MRFZoneList:acceleration",
|
|
U.mesh().time().timeName(),
|
|
U.mesh()
|
|
),
|
|
U.mesh(),
|
|
dimensionedVector("0", U.dimensions()/dimTime, vector::zero)
|
|
)
|
|
);
|
|
volVectorField& acceleration = tacceleration.ref();
|
|
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).addCoriolis(U, acceleration);
|
|
}
|
|
|
|
return tacceleration;
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volVectorField> Foam::MRFZoneList::DDt
|
|
(
|
|
const volScalarField& rho,
|
|
const volVectorField& U
|
|
) const
|
|
{
|
|
return rho*DDt(U);
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::makeRelative(volVectorField& U) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeRelative(U);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::makeRelative(surfaceScalarField& phi) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeRelative(phi);
|
|
}
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceScalarField> Foam::MRFZoneList::relative
|
|
(
|
|
const tmp<surfaceScalarField>& tphi
|
|
) const
|
|
{
|
|
if (size())
|
|
{
|
|
tmp<surfaceScalarField> rphi
|
|
(
|
|
New
|
|
(
|
|
tphi,
|
|
"relative(" + tphi().name() + ')',
|
|
tphi().dimensions(),
|
|
true
|
|
)
|
|
);
|
|
|
|
makeRelative(rphi.ref());
|
|
|
|
tphi.clear();
|
|
|
|
return rphi;
|
|
}
|
|
else
|
|
{
|
|
return tmp<surfaceScalarField>(tphi, true);
|
|
}
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::FieldField<Foam::fvsPatchField, Foam::scalar>>
|
|
Foam::MRFZoneList::relative
|
|
(
|
|
const tmp<FieldField<fvsPatchField, scalar>>& tphi
|
|
) const
|
|
{
|
|
if (size())
|
|
{
|
|
tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
|
|
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeRelative(rphi.ref());
|
|
}
|
|
|
|
tphi.clear();
|
|
|
|
return rphi;
|
|
}
|
|
else
|
|
{
|
|
return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
|
|
}
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::Field<Foam::scalar>>
|
|
Foam::MRFZoneList::relative
|
|
(
|
|
const tmp<Field<scalar>>& tphi,
|
|
const label patchi
|
|
) const
|
|
{
|
|
if (size())
|
|
{
|
|
tmp<Field<scalar>> rphi(New(tphi, true));
|
|
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeRelative(rphi.ref(), patchi);
|
|
}
|
|
|
|
tphi.clear();
|
|
|
|
return rphi;
|
|
}
|
|
else
|
|
{
|
|
return tmp<Field<scalar>>(tphi, true);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::makeRelative
|
|
(
|
|
const surfaceScalarField& rho,
|
|
surfaceScalarField& phi
|
|
) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeRelative(rho, phi);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::makeAbsolute(volVectorField& U) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeAbsolute(U);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::makeAbsolute(surfaceScalarField& phi) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeAbsolute(phi);
|
|
}
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceScalarField> Foam::MRFZoneList::absolute
|
|
(
|
|
const tmp<surfaceScalarField>& tphi
|
|
) const
|
|
{
|
|
if (size())
|
|
{
|
|
tmp<surfaceScalarField> rphi
|
|
(
|
|
New
|
|
(
|
|
tphi,
|
|
"absolute(" + tphi().name() + ')',
|
|
tphi().dimensions(),
|
|
true
|
|
)
|
|
);
|
|
|
|
makeAbsolute(rphi.ref());
|
|
|
|
tphi.clear();
|
|
|
|
return rphi;
|
|
}
|
|
else
|
|
{
|
|
return tmp<surfaceScalarField>(tphi, true);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::makeAbsolute
|
|
(
|
|
const surfaceScalarField& rho,
|
|
surfaceScalarField& phi
|
|
) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).makeAbsolute(rho, phi);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::correctBoundaryVelocity(volVectorField& U) const
|
|
{
|
|
forAll(*this, i)
|
|
{
|
|
operator[](i).correctBoundaryVelocity(U);
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::MRFZoneList::correctBoundaryFlux
|
|
(
|
|
const volVectorField& U,
|
|
surfaceScalarField& phi
|
|
) const
|
|
{
|
|
FieldField<fvsPatchField, scalar> phibf
|
|
(
|
|
relative(mesh_.Sf().boundaryField() & U.boundaryField())
|
|
);
|
|
|
|
forAll(mesh_.boundary(), patchi)
|
|
{
|
|
if
|
|
(
|
|
isA<fixedValueFvsPatchScalarField>(phi.boundaryField()[patchi])
|
|
)
|
|
{
|
|
phi.boundaryField()[patchi] == phibf[patchi];
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
|
|
|
|
Foam::Ostream& Foam::operator<<
|
|
(
|
|
Ostream& os,
|
|
const MRFZoneList& models
|
|
)
|
|
{
|
|
models.writeData(os);
|
|
return os;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|