1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
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::pressureInletVelocityFvPatchVectorField
93 const pressureInletVelocityFvPatchVectorField& pivpvf,
94 const DimensionedField<vector, volMesh>& iF
97 fixedValueFvPatchVectorField(pivpvf, iF),
98 phiName_(pivpvf.phiName_),
99 rhoName_(pivpvf.rhoName_)
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 void pressureInletVelocityFvPatchVectorField::updateCoeffs()
112 const surfaceScalarField& phi =
113 db().lookupObject<surfaceScalarField>(phiName_);
115 const fvsPatchField<scalar>& phip =
116 patch().patchField<surfaceScalarField, scalar>(phi);
118 vectorField n = patch().nf();
119 const Field<scalar>& magS = patch().magSf();
121 if (phi.dimensions() == dimVelocity*dimArea)
123 operator==(n*phip/magS);
125 else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
127 const fvPatchField<scalar>& rhop =
128 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
130 operator==(n*phip/(rhop*magS));
134 FatalErrorIn("pressureInletVelocityFvPatchVectorField::updateCoeffs()")
135 << "dimensions of phi are not correct"
136 << "\n on patch " << this->patch().name()
137 << " of field " << this->dimensionedInternalField().name()
138 << " in file " << this->dimensionedInternalField().objectPath()
142 fixedValueFvPatchVectorField::updateCoeffs();
146 void pressureInletVelocityFvPatchVectorField::write(Ostream& os) const
148 fvPatchVectorField::write(os);
149 if (phiName_ != "phi")
151 os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
153 if (rhoName_ != "rho")
155 os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
157 writeEntry("value", os);
161 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
163 void pressureInletVelocityFvPatchVectorField::operator=
165 const fvPatchField<vector>& pvf
168 fvPatchField<vector>::operator=(patch().nf()*(patch().nf() & pvf));
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 pressureInletVelocityFvPatchVectorField
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 } // End namespace Foam
184 // ************************************************************************* //