New solver: rhoPorousMRFPimpleFoam
[OpenFOAM-1.6.x.git] / src / finiteVolume / cfdTools / general / MRF / MRFZoneTemplates.C
blobb195f111a4457b78991afe1a2144871e4007b04c
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 "MRFZone.H"
28 #include "fvMesh.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "geometricOneField.H"
33 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
35 template<class RhoFieldType>
36 void Foam::MRFZone::relativeRhoFlux
38     const RhoFieldType& rho,
39     surfaceScalarField& phi
40 ) const
42     const surfaceVectorField& Cf = mesh_.Cf();
43     const surfaceVectorField& Sf = mesh_.Sf();
45     const vector& origin = origin_.value();
46     const vector& Omega = Omega_.value();
48     // Internal faces
49     forAll(internalFaces_, i)
50     {
51         label facei = internalFaces_[i];
52         phi[facei] -= rho[facei]*(Omega ^ (Cf[facei] - origin)) & Sf[facei];
53     }
55     // Included patches
56     forAll(includedFaces_, patchi)
57     {
58         forAll(includedFaces_[patchi], i)
59         {
60             label patchFacei = includedFaces_[patchi][i];
62             phi.boundaryField()[patchi][patchFacei] = 0.0;
63         }
64     }
66     // Excluded patches
67     forAll(excludedFaces_, patchi)
68     {
69         forAll(excludedFaces_[patchi], i)
70         {
71             label patchFacei = excludedFaces_[patchi][i];
73             phi.boundaryField()[patchi][patchFacei] -=
74                 rho.boundaryField()[patchi][patchFacei]
75                *(Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
76               & Sf.boundaryField()[patchi][patchFacei];
77         }
78     }
82 template<class RhoFieldType>
83 void Foam::MRFZone::absoluteRhoFlux
85     const RhoFieldType& rho,
86     surfaceScalarField& phi
87 ) const
89     const surfaceVectorField& Cf = mesh_.Cf();
90     const surfaceVectorField& Sf = mesh_.Sf();
92     const vector& origin = origin_.value();
93     const vector& Omega = Omega_.value();
95     // Internal faces
96     forAll(internalFaces_, i)
97     {
98         label facei = internalFaces_[i];
99         phi[facei] += (Omega ^ (Cf[facei] - origin)) & Sf[facei];
100     }
102     // Included patches
103     forAll(includedFaces_, patchi)
104     {
105         forAll(includedFaces_[patchi], i)
106         {
107             label patchFacei = includedFaces_[patchi][i];
109             phi.boundaryField()[patchi][patchFacei] +=
110                 (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
111               & Sf.boundaryField()[patchi][patchFacei];
112         }
113     }
115     // Excluded patches
116     forAll(excludedFaces_, patchi)
117     {
118         forAll(excludedFaces_[patchi], i)
119         {
120             label patchFacei = excludedFaces_[patchi][i];
122             phi.boundaryField()[patchi][patchFacei] +=
123                 (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
124               & Sf.boundaryField()[patchi][patchFacei];
125         }
126     }
130 // ************************************************************************* //