1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software; you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "buoyantPressureFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "uniformDimensionedFields.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 buoyantPressureFvPatchScalarField::
41 buoyantPressureFvPatchScalarField
44 const DimensionedField<scalar, volMesh>& iF
47 fixedGradientFvPatchScalarField(p, iF),
52 buoyantPressureFvPatchScalarField::
53 buoyantPressureFvPatchScalarField
56 const DimensionedField<scalar, volMesh>& iF,
57 const dictionary& dict
60 fixedGradientFvPatchScalarField(p, iF),
61 rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
63 fvPatchField<scalar>::operator=(patchInternalField());
68 buoyantPressureFvPatchScalarField::
69 buoyantPressureFvPatchScalarField
71 const buoyantPressureFvPatchScalarField& ptf,
73 const DimensionedField<scalar, volMesh>& iF,
74 const fvPatchFieldMapper& mapper
77 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
78 rhoName_(ptf.rhoName_)
82 buoyantPressureFvPatchScalarField::
83 buoyantPressureFvPatchScalarField
85 const buoyantPressureFvPatchScalarField& ptf
88 fixedGradientFvPatchScalarField(ptf),
89 rhoName_(ptf.rhoName_)
93 buoyantPressureFvPatchScalarField::
94 buoyantPressureFvPatchScalarField
96 const buoyantPressureFvPatchScalarField& ptf,
97 const DimensionedField<scalar, volMesh>& iF
100 fixedGradientFvPatchScalarField(ptf, iF),
101 rhoName_(ptf.rhoName_)
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 void buoyantPressureFvPatchScalarField::updateCoeffs()
114 const uniformDimensionedVectorField& g =
115 db().lookupObject<uniformDimensionedVectorField>("g");
117 const fvPatchField<scalar>& rho =
118 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
120 // If the variable name is "p_rgh" or "pd" assume it is p - rho*g.h
121 // and set the gradient appropriately.
122 // Otherwise assume the variable is the static pressure.
125 dimensionedInternalField().name() == "p_rgh"
126 || dimensionedInternalField().name() == "pd"
129 gradient() = -rho.snGrad()*(g.value() & patch().Cf());
133 gradient() = rho*(g.value() & patch().nf());
136 fixedGradientFvPatchScalarField::updateCoeffs();
140 void buoyantPressureFvPatchScalarField::write(Ostream& os) const
142 fixedGradientFvPatchScalarField::write(os);
143 os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
144 writeEntry("value", os);
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 buoyantPressureFvPatchScalarField
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 } // End namespace Foam
160 // ************************************************************************* //