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

MRFZones.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 "MRFZones.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <finiteVolume/fvMesh.H>
00029 
00030 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034     defineTemplateTypeNameAndDebug(IOPtrList<MRFZone>, 0);
00035 }
00036 
00037 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00038 
00039 Foam::MRFZones::MRFZones(const fvMesh& mesh)
00040 :
00041     IOPtrList<MRFZone>
00042     (
00043         IOobject
00044         (
00045             "MRFZones",
00046             mesh.time().constant(),
00047             mesh,
00048             IOobject::MUST_READ,
00049             IOobject::NO_WRITE
00050         ),
00051         MRFZone::iNew(mesh)
00052     )
00053 {}
00054 
00055 
00056 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00057 
00058 void Foam::MRFZones::addCoriolis(fvVectorMatrix& UEqn) const
00059 {
00060     forAll(*this, i)
00061     {
00062         operator[](i).addCoriolis(UEqn);
00063     }
00064 }
00065 
00066 
00067 void Foam::MRFZones::addCoriolis
00068 (
00069     const volScalarField& rho,
00070     fvVectorMatrix& UEqn
00071 ) const
00072 {
00073     forAll(*this, i)
00074     {
00075         operator[](i).addCoriolis(rho, UEqn);
00076     }
00077 }
00078 
00079 
00080 void Foam::MRFZones::relativeVelocity(volVectorField& U) const
00081 {
00082     forAll(*this, i)
00083     {
00084         operator[](i).relativeVelocity(U);
00085     }
00086 }
00087 
00088 
00089 void Foam::MRFZones::absoluteVelocity(volVectorField& U) const
00090 {
00091     forAll(*this, i)
00092     {
00093         operator[](i).absoluteVelocity(U);
00094     }
00095 }
00096 
00097 
00098 void Foam::MRFZones::relativeFlux(surfaceScalarField& phi) const
00099 {
00100     forAll(*this, i)
00101     {
00102         operator[](i).relativeFlux(phi);
00103     }
00104 }
00105 
00106 
00107 void Foam::MRFZones::relativeFlux
00108 (
00109     const surfaceScalarField& rho,
00110     surfaceScalarField& phi
00111 ) const
00112 {
00113     forAll(*this, i)
00114     {
00115         operator[](i).relativeFlux(rho, phi);
00116     }
00117 }
00118 
00119 
00120 void Foam::MRFZones::absoluteFlux(surfaceScalarField& phi) const
00121 {
00122     forAll(*this, i)
00123     {
00124         operator[](i).absoluteFlux(phi);
00125     }
00126 }
00127 
00128 
00129 void Foam::MRFZones::absoluteFlux
00130 (
00131     const surfaceScalarField& rho,
00132     surfaceScalarField& phi
00133 ) const
00134 {
00135     forAll(*this, i)
00136     {
00137         operator[](i).absoluteFlux(rho, phi);
00138     }
00139 }
00140 
00141 
00142 void Foam::MRFZones::correctBoundaryVelocity(volVectorField& U) const
00143 {
00144     forAll(*this, i)
00145     {
00146         operator[](i).correctBoundaryVelocity(U);
00147     }
00148 }
00149 
00150 
00151 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines