initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / fields / fvPatchFields / derived / freestreamPressure / freestreamPressureFvPatchScalarField.C
blobf6f94738c6792e08082f05d1bf5228786fbd896b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 "freestreamPressureFvPatchScalarField.H"
28 #include "freestreamFvPatchFields.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 namespace Foam
39 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
41 freestreamPressureFvPatchScalarField::freestreamPressureFvPatchScalarField
43     const fvPatch& p,
44     const DimensionedField<scalar, volMesh>& iF
47     zeroGradientFvPatchScalarField(p, iF)
51 freestreamPressureFvPatchScalarField::freestreamPressureFvPatchScalarField
53     const freestreamPressureFvPatchScalarField& ptf,
54     const fvPatch& p,
55     const DimensionedField<scalar, volMesh>& iF,
56     const fvPatchFieldMapper& mapper
59     zeroGradientFvPatchScalarField(ptf, p, iF, mapper)
63 freestreamPressureFvPatchScalarField::freestreamPressureFvPatchScalarField
65     const fvPatch& p,
66     const DimensionedField<scalar, volMesh>& iF,
67     const dictionary& dict
70     zeroGradientFvPatchScalarField(p, iF, dict)
74 freestreamPressureFvPatchScalarField::freestreamPressureFvPatchScalarField
76     const freestreamPressureFvPatchScalarField& wbppsf
79     zeroGradientFvPatchScalarField(wbppsf)
83 freestreamPressureFvPatchScalarField::freestreamPressureFvPatchScalarField
85     const freestreamPressureFvPatchScalarField& wbppsf,
86     const DimensionedField<scalar, volMesh>& iF
89     zeroGradientFvPatchScalarField(wbppsf, iF)
93 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
95 void freestreamPressureFvPatchScalarField::updateCoeffs()
97     if (updated())
98     {
99         return;
100     }
102     const freestreamFvPatchVectorField& Up = 
103         refCast<const freestreamFvPatchVectorField>
104         (
105             patch().lookupPatchField<volVectorField, vector>("U")
106         );
108     const surfaceScalarField& phi = 
109         db().lookupObject<surfaceScalarField>("phi");
111     fvsPatchField<scalar>& phip =
112         const_cast<fvsPatchField<scalar>&>
113         (
114             patch().patchField<surfaceScalarField, scalar>(phi)
115         );
117     if (phi.dimensions() == dimVelocity*dimArea)
118     {
119         phip = patch().Sf() & Up.freestreamValue();
120     }
121     else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
122     {
123         const fvPatchField<scalar>& rhop =
124             patch().lookupPatchField<volScalarField, scalar>("rho");
126         phip = rhop*(patch().Sf() & Up.freestreamValue());
127     }
128     else
129     {
130         FatalErrorIn("freestreamPressureFvPatchScalarField::updateCoeffs()")
131             << "dimensions of phi are not correct"
132             << "\n    on patch " << this->patch().name()
133             << " of field " << this->dimensionedInternalField().name()
134             << " in file " << this->dimensionedInternalField().objectPath()
135             << exit(FatalError);
136     }
138     zeroGradientFvPatchScalarField::updateCoeffs();
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 makePatchTypeField(fvPatchScalarField, freestreamPressureFvPatchScalarField);
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 } // End namespace Foam
150 // ************************************************************************* //