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 "fixedFluxPressureFvPatchScalarField.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
38 const DimensionedField<scalar, volMesh>& iF
41 fixedGradientFvPatchScalarField(p, iF),
49 Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
51 const fixedFluxPressureFvPatchScalarField& ptf,
53 const DimensionedField<scalar, volMesh>& iF,
54 const fvPatchFieldMapper& mapper
57 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
59 phiName_(ptf.phiName_),
60 rhoName_(ptf.rhoName_),
61 adjoint_(ptf.adjoint_)
65 Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
68 const DimensionedField<scalar, volMesh>& iF,
69 const dictionary& dict
72 fixedGradientFvPatchScalarField(p, iF),
73 UName_(dict.lookupOrDefault<word>("U", "U")),
74 phiName_(dict.lookupOrDefault<word>("phi", "phi")),
75 rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
76 adjoint_(dict.lookup("adjoint"))
78 if (dict.found("gradient"))
80 gradient() = scalarField("gradient", dict, p.size());
81 fixedGradientFvPatchScalarField::updateCoeffs();
82 fixedGradientFvPatchScalarField::evaluate();
86 fvPatchField<scalar>::operator=(patchInternalField());
92 Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
94 const fixedFluxPressureFvPatchScalarField& wbppsf
97 fixedGradientFvPatchScalarField(wbppsf),
98 UName_(wbppsf.UName_),
99 phiName_(wbppsf.phiName_),
100 rhoName_(wbppsf.rhoName_),
101 adjoint_(wbppsf.adjoint_)
105 Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
107 const fixedFluxPressureFvPatchScalarField& wbppsf,
108 const DimensionedField<scalar, volMesh>& iF
111 fixedGradientFvPatchScalarField(wbppsf, iF),
112 UName_(wbppsf.UName_),
113 phiName_(wbppsf.phiName_),
114 rhoName_(wbppsf.rhoName_),
115 adjoint_(wbppsf.adjoint_)
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 void Foam::fixedFluxPressureFvPatchScalarField::updateCoeffs()
128 const fvPatchField<vector>& Up =
129 patch().lookupPatchField<volVectorField, vector>(UName_);
131 const surfaceScalarField& phi =
132 db().lookupObject<surfaceScalarField>(phiName_);
134 fvsPatchField<scalar> phip =
135 patch().patchField<surfaceScalarField, scalar>(phi);
137 if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
139 const fvPatchField<scalar>& rhop =
140 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
145 const fvPatchField<scalar>& rAp =
146 patch().lookupPatchField<volScalarField, scalar>("(1|A("+UName_+"))");
150 gradient() = ((patch().Sf() & Up) - phip)/patch().magSf()/rAp;
154 gradient() = (phip - (patch().Sf() & Up))/patch().magSf()/rAp;
157 fixedGradientFvPatchScalarField::updateCoeffs();
161 void Foam::fixedFluxPressureFvPatchScalarField::write(Ostream& os) const
163 fvPatchScalarField::write(os);
166 os.writeKeyword("U") << UName_ << token::END_STATEMENT << nl;
168 if (phiName_ != "phi")
170 os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
172 if (rhoName_ != "rho")
174 os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
176 os.writeKeyword("adjoint") << adjoint_ << token::END_STATEMENT << nl;
177 gradient().writeEntry("gradient", os);
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 fixedFluxPressureFvPatchScalarField
192 // ************************************************************************* //