ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / finiteVolume / cfdTools / general / MRF / MRFZoneTemplates.C
blob55821dbfde6ab56fa5ff4370bb808722bd0b7597
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "MRFZone.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "geometricOneField.H"
32 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
34 template<class RhoFieldType>
35 void Foam::MRFZone::relativeRhoFlux
37     const RhoFieldType& rho,
38     surfaceScalarField& phi
39 ) const
41     const surfaceVectorField& Cf = mesh_.Cf();
42     const surfaceVectorField& Sf = mesh_.Sf();
44     const vector& origin = origin_.value();
45     const vector& Omega = Omega_.value();
47     // Internal faces
48     forAll(internalFaces_, i)
49     {
50         label facei = internalFaces_[i];
51         phi[facei] -= rho[facei]*(Omega ^ (Cf[facei] - origin)) & Sf[facei];
52     }
54     // Included patches
55     forAll(includedFaces_, patchi)
56     {
57         forAll(includedFaces_[patchi], i)
58         {
59             label patchFacei = includedFaces_[patchi][i];
61             phi.boundaryField()[patchi][patchFacei] = 0.0;
62         }
63     }
65     // Excluded patches
66     forAll(excludedFaces_, patchi)
67     {
68         forAll(excludedFaces_[patchi], i)
69         {
70             label patchFacei = excludedFaces_[patchi][i];
72             phi.boundaryField()[patchi][patchFacei] -=
73                 rho.boundaryField()[patchi][patchFacei]
74               * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
75               & Sf.boundaryField()[patchi][patchFacei];
76         }
77     }
81 template<class RhoFieldType>
82 void Foam::MRFZone::absoluteRhoFlux
84     const RhoFieldType& rho,
85     surfaceScalarField& phi
86 ) const
88     const surfaceVectorField& Cf = mesh_.Cf();
89     const surfaceVectorField& Sf = mesh_.Sf();
91     const vector& origin = origin_.value();
92     const vector& Omega = Omega_.value();
94     // Internal faces
95     forAll(internalFaces_, i)
96     {
97         label facei = internalFaces_[i];
98         phi[facei] += rho[facei]*(Omega ^ (Cf[facei] - origin)) & Sf[facei];
99     }
101     // Included patches
102     forAll(includedFaces_, patchi)
103     {
104         forAll(includedFaces_[patchi], i)
105         {
106             label patchFacei = includedFaces_[patchi][i];
108             phi.boundaryField()[patchi][patchFacei] +=
109                 rho.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                 rho.boundaryField()[patchi][patchFacei]
124               * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
125               & Sf.boundaryField()[patchi][patchFacei];
126         }
127     }
131 // ************************************************************************* //