FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

MRFZoneTemplates.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "MRFZone.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <finiteVolume/surfaceFields.H>
00030 #include <OpenFOAM/geometricOneField.H>
00031 
00032 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00033 
00034 template<class RhoFieldType>
00035 void Foam::MRFZone::relativeRhoFlux
00036 (
00037     const RhoFieldType& rho,
00038     surfaceScalarField& phi
00039 ) const
00040 {
00041     const surfaceVectorField& Cf = mesh_.Cf();
00042     const surfaceVectorField& Sf = mesh_.Sf();
00043 
00044     const vector& origin = origin_.value();
00045     const vector& Omega = Omega_.value();
00046 
00047     // Internal faces
00048     forAll(internalFaces_, i)
00049     {
00050         label facei = internalFaces_[i];
00051         phi[facei] -= rho[facei]*(Omega ^ (Cf[facei] - origin)) & Sf[facei];
00052     }
00053 
00054     // Included patches
00055     forAll(includedFaces_, patchi)
00056     {
00057         forAll(includedFaces_[patchi], i)
00058         {
00059             label patchFacei = includedFaces_[patchi][i];
00060 
00061             phi.boundaryField()[patchi][patchFacei] = 0.0;
00062         }
00063     }
00064 
00065     // Excluded patches
00066     forAll(excludedFaces_, patchi)
00067     {
00068         forAll(excludedFaces_[patchi], i)
00069         {
00070             label patchFacei = excludedFaces_[patchi][i];
00071 
00072             phi.boundaryField()[patchi][patchFacei] -=
00073                 rho.boundaryField()[patchi][patchFacei]
00074               * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
00075               & Sf.boundaryField()[patchi][patchFacei];
00076         }
00077     }
00078 }
00079 
00080 
00081 template<class RhoFieldType>
00082 void Foam::MRFZone::absoluteRhoFlux
00083 (
00084     const RhoFieldType& rho,
00085     surfaceScalarField& phi
00086 ) const
00087 {
00088     const surfaceVectorField& Cf = mesh_.Cf();
00089     const surfaceVectorField& Sf = mesh_.Sf();
00090 
00091     const vector& origin = origin_.value();
00092     const vector& Omega = Omega_.value();
00093 
00094     // Internal faces
00095     forAll(internalFaces_, i)
00096     {
00097         label facei = internalFaces_[i];
00098         phi[facei] += rho[facei]*(Omega ^ (Cf[facei] - origin)) & Sf[facei];
00099     }
00100 
00101     // Included patches
00102     forAll(includedFaces_, patchi)
00103     {
00104         forAll(includedFaces_[patchi], i)
00105         {
00106             label patchFacei = includedFaces_[patchi][i];
00107 
00108             phi.boundaryField()[patchi][patchFacei] +=
00109                 rho.boundaryField()[patchi][patchFacei]
00110               * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
00111               & Sf.boundaryField()[patchi][patchFacei];
00112         }
00113     }
00114 
00115     // Excluded patches
00116     forAll(excludedFaces_, patchi)
00117     {
00118         forAll(excludedFaces_[patchi], i)
00119         {
00120             label patchFacei = excludedFaces_[patchi][i];
00121 
00122             phi.boundaryField()[patchi][patchFacei] +=
00123                 rho.boundaryField()[patchi][patchFacei]
00124               * (Omega ^ (Cf.boundaryField()[patchi][patchFacei] - origin))
00125               & Sf.boundaryField()[patchi][patchFacei];
00126         }
00127     }
00128 }
00129 
00130 
00131 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines