This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/compressible/rhopSonicFoam/BCs/rho/gradientRhoFvPatchScalarField.C
2016-06-21 15:04:12 +02:00

139 lines
3.8 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend 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.
foam-extend 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 foam-extend. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "gradientRhoFvPatchScalarField.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
gradientRhoFvPatchScalarField::gradientRhoFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedGradientFvPatchScalarField(p, iF)
{}
gradientRhoFvPatchScalarField::gradientRhoFvPatchScalarField
(
const gradientRhoFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedGradientFvPatchScalarField(ptf, p, iF, mapper)
{}
gradientRhoFvPatchScalarField::gradientRhoFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
fixedGradientFvPatchScalarField(p, iF)
{
if (dict.found("gradient"))
{
gradient() = scalarField("gradient", dict, p.size());
fixedGradientFvPatchScalarField::updateCoeffs();
fixedGradientFvPatchScalarField::evaluate();
}
else
{
fvPatchField<scalar>::operator=(patchInternalField());
gradient() = 0.0;
}
}
gradientRhoFvPatchScalarField::gradientRhoFvPatchScalarField
(
const gradientRhoFvPatchScalarField& wbppsf
)
:
fixedGradientFvPatchScalarField(wbppsf)
{}
gradientRhoFvPatchScalarField::gradientRhoFvPatchScalarField
(
const gradientRhoFvPatchScalarField& wbppsf,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedGradientFvPatchScalarField(wbppsf, iF)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void gradientRhoFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
const fvPatchField<scalar>& psip =
lookupPatchField<volScalarField, scalar>("psi");
const fvPatchField<scalar>& pp =
lookupPatchField<volScalarField, scalar>("p");
gradient() = psip*pp.snGrad() + psip.snGrad()*pp;
fixedGradientFvPatchScalarField::updateCoeffs();
}
void gradientRhoFvPatchScalarField::write(Ostream& os) const
{
fixedGradientFvPatchScalarField::write(os);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePatchTypeField(fvPatchScalarField, gradientRhoFvPatchScalarField);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// ************************************************************************* //