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 "pressureInletVelocityFvPatchVectorField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 pressureInletVelocityFvPatchVectorField::pressureInletVelocityFvPatchVectorField
42 const DimensionedField<vector, volMesh>& iF
45 fixedValueFvPatchVectorField(p, iF),
51 pressureInletVelocityFvPatchVectorField::pressureInletVelocityFvPatchVectorField
53 const pressureInletVelocityFvPatchVectorField& ptf,
55 const DimensionedField<vector, volMesh>& iF,
56 const fvPatchFieldMapper& mapper
59 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
60 phiName_(ptf.phiName_),
61 rhoName_(ptf.rhoName_)
65 pressureInletVelocityFvPatchVectorField::pressureInletVelocityFvPatchVectorField
68 const DimensionedField<vector, volMesh>& iF,
69 const dictionary& dict
72 fixedValueFvPatchVectorField(p, iF),
73 phiName_(dict.lookupOrDefault<word>("phi", "phi")),
74 rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
76 fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
80 pressureInletVelocityFvPatchVectorField::pressureInletVelocityFvPatchVectorField
82 const pressureInletVelocityFvPatchVectorField& pivpvf
85 fixedValueFvPatchVectorField(pivpvf),
86 phiName_(pivpvf.phiName_),
87 rhoName_(pivpvf.rhoName_)
91 pressureInletVelocityFvPatchVectorField::
92 pressureInletVelocityFvPatchVectorField
94 const pressureInletVelocityFvPatchVectorField& pivpvf,
95 const DimensionedField<vector, volMesh>& iF
98 fixedValueFvPatchVectorField(pivpvf, iF),
99 phiName_(pivpvf.phiName_),
100 rhoName_(pivpvf.rhoName_)
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 void pressureInletVelocityFvPatchVectorField::updateCoeffs()
113 const surfaceScalarField& phi =
114 db().lookupObject<surfaceScalarField>(phiName_);
116 const fvsPatchField<scalar>& phip =
117 patch().patchField<surfaceScalarField, scalar>(phi);
119 vectorField n = patch().nf();
120 const Field<scalar>& magS = patch().magSf();
122 if (phi.dimensions() == dimVelocity*dimArea)
124 operator==(n*phip/magS);
126 else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
128 const fvPatchField<scalar>& rhop =
129 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
131 operator==(n*phip/(rhop*magS));
135 FatalErrorIn("pressureInletVelocityFvPatchVectorField::updateCoeffs()")
136 << "dimensions of phi are not correct"
137 << "\n on patch " << this->patch().name()
138 << " of field " << this->dimensionedInternalField().name()
139 << " in file " << this->dimensionedInternalField().objectPath()
143 fixedValueFvPatchVectorField::updateCoeffs();
147 void pressureInletVelocityFvPatchVectorField::write(Ostream& os) const
149 fvPatchVectorField::write(os);
150 if (phiName_ != "phi")
152 os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
154 if (rhoName_ != "rho")
156 os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
158 writeEntry("value", os);
162 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
164 void pressureInletVelocityFvPatchVectorField::operator=
166 const fvPatchField<vector>& pvf
169 fvPatchField<vector>::operator=(patch().nf()*(patch().nf() & pvf));
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 pressureInletVelocityFvPatchVectorField
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 } // End namespace Foam
185 // ************************************************************************* //