/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2015-2017 OpenFOAM Foundation
Copyright (C) 2018-2020 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 .
\*---------------------------------------------------------------------------*/
#include "prghTotalPressureFvPatchScalarField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "gravityMeshObject.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::prghTotalPressureFvPatchScalarField::prghTotalPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField& iF
)
:
fixedValueFvPatchScalarField(p, iF),
UName_("U"),
phiName_("phi"),
rhoName_("rho"),
p0_(p.size(), Zero)
{}
Foam::prghTotalPressureFvPatchScalarField::prghTotalPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField& iF,
const dictionary& dict
)
:
fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
UName_(dict.getOrDefault("U", "U")),
phiName_(dict.getOrDefault("phi", "phi")),
rhoName_(dict.getOrDefault("rho", "rho")),
p0_("p0", dict, p.size())
{
if (!this->readValueEntry(dict))
{
fvPatchField::operator=(p0_);
}
}
Foam::prghTotalPressureFvPatchScalarField::prghTotalPressureFvPatchScalarField
(
const prghTotalPressureFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchScalarField(ptf, p, iF, mapper),
UName_(ptf.UName_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
p0_(ptf.p0_, mapper)
{}
Foam::prghTotalPressureFvPatchScalarField::prghTotalPressureFvPatchScalarField
(
const prghTotalPressureFvPatchScalarField& ptf
)
:
fixedValueFvPatchScalarField(ptf),
UName_(ptf.UName_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
p0_(ptf.p0_)
{}
Foam::prghTotalPressureFvPatchScalarField::prghTotalPressureFvPatchScalarField
(
const prghTotalPressureFvPatchScalarField& ptf,
const DimensionedField& iF
)
:
fixedValueFvPatchScalarField(ptf, iF),
UName_(ptf.UName_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
p0_(ptf.p0_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::prghTotalPressureFvPatchScalarField::autoMap
(
const fvPatchFieldMapper& m
)
{
fixedValueFvPatchScalarField::autoMap(m);
p0_.autoMap(m);
}
void Foam::prghTotalPressureFvPatchScalarField::rmap
(
const fvPatchScalarField& ptf,
const labelList& addr
)
{
fixedValueFvPatchScalarField::rmap(ptf, addr);
const prghTotalPressureFvPatchScalarField& tiptf =
refCast(ptf);
p0_.rmap(tiptf.p0_, addr);
}
void Foam::prghTotalPressureFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
const scalarField& rhop =
patch().lookupPatchField(rhoName_);
const scalarField& phip =
patch().lookupPatchField(phiName_);
const vectorField& Up =
patch().lookupPatchField(UName_);
const uniformDimensionedVectorField& g =
meshObjects::gravity::New(db().time());
const uniformDimensionedScalarField& hRef =
db().lookupObject("hRef");
dimensionedScalar ghRef
(
mag(g.value()) > SMALL
? g & (cmptMag(g.value())/mag(g.value()))*hRef
: dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
);
operator==
(
p0_
- 0.5*rhop*(neg(phip))*magSqr(Up)
- rhop*((g.value() & patch().Cf()) - ghRef.value())
);
fixedValueFvPatchScalarField::updateCoeffs();
}
void Foam::prghTotalPressureFvPatchScalarField::write(Ostream& os) const
{
fvPatchField::write(os);
os.writeEntryIfDifferent("U", "U", UName_);
os.writeEntryIfDifferent("phi", "phi", phiName_);
os.writeEntryIfDifferent("rho", "rho", rhoName_);
p0_.writeEntry("p0", os);
fvPatchField::writeValueEntry(os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
prghTotalPressureFvPatchScalarField
);
}
// ************************************************************************* //