initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / fields / fvPatchFields / derived / fan / fanFvPatchFields.C
blob100fe188a20b3043123e7aab6f9e36912465b4f2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 "fanFvPatchFields.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 makePatchTypeField(fvPatchScalarField, fanFvPatchScalarField);
42 //- Specialisation of the jump-condition for the pressure
43 template<>
44 void fanFvPatchField<scalar>::updateCoeffs()
46     if (updated())
47     {
48         return;
49     }
51     jump_ = f_[0];
53     if (f_.size() > 1)
54     {
55         const surfaceScalarField& phi =
56             db().lookupObject<surfaceScalarField>("phi");
58         const fvsPatchField<scalar>& phip =
59             patch().patchField<surfaceScalarField, scalar>(phi);
61         scalarField Un =
62             scalarField::subField(phip, size()/2)
63            /scalarField::subField(patch().magSf(), size()/2);
65         if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
66         {
67             Un /= patch().lookupPatchField<volScalarField, scalar>("rho");
68         }
70         for(label i=1; i<f_.size(); i++)
71         {
72             jump_ += f_[i]*pow(Un, i);
73         }
74     }
76     jumpCyclicFvPatchField<scalar>::updateCoeffs();
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 } // End namespace Foam
84 // ************************************************************************* //