Fix tutorials: coupled/conjugateHeatFoam/conjugateCavity: fix Allrun file
[OpenFOAM-1.6-ext.git] / src / finiteVolume / fields / fvPatchFields / derived / fan / fanFvPatchFields.C
blobab258f1055b6f9fde7073a4ba2f0d1bba306e1c3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 namespace Foam
36     makePatchTypeField(fvPatchScalarField, fanFvPatchScalarField);
39 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
41 //- Specialisation of the jump-condition for the pressure
42 template<>
43 void Foam::fanFvPatchField<Foam::scalar>::updateCoeffs()
45     if (updated())
46     {
47         return;
48     }
50     jump_ = f_[0];
52     if (f_.size() > 1)
53     {
54         const surfaceScalarField& phi =
55             db().lookupObject<surfaceScalarField>("phi");
57         const fvsPatchField<scalar>& phip =
58             patch().patchField<surfaceScalarField, scalar>(phi);
60         scalarField Un = max
61         (
62             scalarField::subField(phip, size()/2)
63            /scalarField::subField(patch().magSf(), size()/2),
64             scalar(0)
65         );
67         if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
68         {
69             Un /=
70                 scalarField::subField
71                 (
72                     patch().lookupPatchField<volScalarField, scalar>("rho"),
73                     size()/2
74                 );
75         }
77         for(label i=1; i<f_.size(); i++)
78         {
79             jump_ += f_[i]*pow(Un, i);
80         }
82         jump_ = max(jump_, scalar(0));
83     }
85     jumpCyclicFvPatchField<scalar>::updateCoeffs();
89 // ************************************************************************* //