1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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 \*---------------------------------------------------------------------------*/
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
41 const surfaceVectorField& Cf = mesh_.Cf();
42 const surfaceVectorField& Sf = mesh_.Sf();
44 const vector& origin = origin_.value();
45 const vector& Omega = Omega_.value();
48 forAll(internalFaces_, i)
50 label facei = internalFaces_[i];
51 phi[facei] -= rho[facei]*(Omega ^ (Cf[facei] - origin)) & Sf[facei];
55 forAll(includedFaces_, patchi)
57 forAll(includedFaces_[patchi], i)
59 label patchFacei = includedFaces_[patchi][i];
61 phi.boundaryField()[patchi][patchFacei] = 0.0;
66 forAll(excludedFaces_, patchi)
68 forAll(excludedFaces_[patchi], i)
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];
81 template<class RhoFieldType>
82 void Foam::MRFZone::absoluteRhoFlux
84 const RhoFieldType& rho,
85 surfaceScalarField& phi
88 const surfaceVectorField& Cf = mesh_.Cf();
89 const surfaceVectorField& Sf = mesh_.Sf();
91 const vector& origin = origin_.value();
92 const vector& Omega = Omega_.value();
95 forAll(internalFaces_, i)
97 label facei = internalFaces_[i];
98 phi[facei] += rho[facei]*(Omega ^ (Cf[facei] - origin)) & Sf[facei];
102 forAll(includedFaces_, patchi)
104 forAll(includedFaces_[patchi], i)
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];
116 forAll(excludedFaces_, patchi)
118 forAll(excludedFaces_[patchi], i)
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];
131 // ************************************************************************* //