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

SRFModel.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 "SRFModel.H"
00027 #include <finiteVolume/SRFVelocityFvPatchVectorField.H>
00028 
00029 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033     namespace SRF
00034     {
00035         defineTypeNameAndDebug(SRFModel, 0);
00036         defineRunTimeSelectionTable(SRFModel, dictionary);
00037     }
00038 }
00039 
00040 
00041 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00042 
00043 Foam::SRF::SRFModel::SRFModel
00044 (
00045     const word& type,
00046     const volVectorField& Urel
00047 )
00048 :
00049     IOdictionary
00050     (
00051         IOobject
00052         (
00053             "SRFProperties",
00054             Urel.time().constant(),
00055             Urel.db(),
00056             IOobject::MUST_READ,
00057             IOobject::NO_WRITE
00058         )
00059     ),
00060     Urel_(Urel),
00061     mesh_(Urel_.mesh()),
00062     axis_(lookup("axis")),
00063     SRFModelCoeffs_(subDict(type + "Coeffs")),
00064     omega_(dimensionedVector("omega", dimless/dimTime, vector::zero))
00065 {
00066     // Normalise the axis
00067     axis_ /= mag(axis_);
00068 }
00069 
00070 
00071 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00072 
00073 Foam::SRF::SRFModel::~SRFModel()
00074 {}
00075 
00076 
00077 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00078 
00079 bool Foam::SRF::SRFModel::read()
00080 {
00081     if (regIOobject::read())
00082     {
00083         // Re-read axis
00084         SRFModelCoeffs_.lookup("axis") >> axis_;
00085         axis_ /= mag(axis_);
00086 
00087         // Re-read sub-model coeffs
00088         SRFModelCoeffs_ = subDict(type() + "Coeffs");
00089 
00090         return true;
00091     }
00092     else
00093     {
00094         return false;
00095     }
00096 }
00097 
00098 
00099 const Foam::vector& Foam::SRF::SRFModel::axis() const
00100 {
00101     return axis_;
00102 }
00103 
00104 
00105 const Foam::dimensionedVector& Foam::SRF::SRFModel::omega() const
00106 {
00107     return omega_;
00108 }
00109 
00110 
00111 Foam::tmp<Foam::DimensionedField<Foam::vector, Foam::volMesh> >
00112 Foam::SRF::SRFModel::Fcoriolis() const
00113 {
00114     return tmp<DimensionedField<vector, volMesh> >
00115     (
00116         new DimensionedField<vector, volMesh>
00117         (
00118             IOobject
00119             (
00120                 "Fcoriolis",
00121                 mesh_.time().timeName(),
00122                 mesh_,
00123                 IOobject::NO_READ,
00124                 IOobject::NO_WRITE
00125             ),
00126             2.0*omega_ ^ Urel_
00127         )
00128     );
00129 }
00130 
00131 
00132 Foam::tmp<Foam::DimensionedField<Foam::vector, Foam::volMesh> >
00133 Foam::SRF::SRFModel::Fcentrifugal() const
00134 {
00135     return tmp<DimensionedField<vector, volMesh> >
00136     (
00137         new DimensionedField<vector, volMesh>
00138         (
00139             IOobject
00140             (
00141                 "Fcentrifugal",
00142                 mesh_.time().timeName(),
00143                 mesh_,
00144                 IOobject::NO_READ,
00145                 IOobject::NO_WRITE
00146             ),
00147             omega_ ^ (omega_ ^ mesh_.C())
00148         )
00149     );
00150 }
00151 
00152 
00153 Foam::tmp<Foam::DimensionedField<Foam::vector, Foam::volMesh> >
00154 Foam::SRF::SRFModel::Su() const
00155 {
00156     return Fcoriolis() + Fcentrifugal();
00157 }
00158 
00159 
00160 Foam::vectorField Foam::SRF::SRFModel::velocity
00161 (
00162     const vectorField& positions
00163 ) const
00164 {
00165     return omega_.value() ^ (positions - axis_*(axis_ & positions));
00166 }
00167 
00168 
00169 Foam::tmp<Foam::volVectorField> Foam::SRF::SRFModel::U() const
00170 {
00171     return tmp<volVectorField>
00172     (
00173         new volVectorField
00174         (
00175             IOobject
00176             (
00177                 "Usrf",
00178                 mesh_.time().timeName(),
00179                 mesh_,
00180                 IOobject::NO_READ,
00181                 IOobject::NO_WRITE
00182             ),
00183             omega_ ^ (mesh_.C() - axis_*(axis_ & mesh_.C()))
00184         )
00185     );
00186 }
00187 
00188 
00189 Foam::tmp<Foam::volVectorField> Foam::SRF::SRFModel::Uabs() const
00190 {
00191     const volVectorField Usrf = U();
00192 
00193     tmp<volVectorField> tUabs
00194     (
00195         new volVectorField
00196         (
00197             IOobject
00198             (
00199                 "Uabs",
00200                 mesh_.time().timeName(),
00201                 mesh_,
00202                 IOobject::NO_READ,
00203                 IOobject::NO_WRITE,
00204                 false
00205             ),
00206             Usrf
00207         )
00208     );
00209 
00210     // Add SRF contribution to internal field
00211     tUabs().internalField() += Urel_.internalField();
00212 
00213     // Add Urel boundary contributions
00214     const volVectorField::GeometricBoundaryField& bvf = Urel_.boundaryField();
00215 
00216     forAll(bvf, i)
00217     {
00218         if (isA<SRFVelocityFvPatchVectorField>(bvf[i]))
00219         {
00220             // Only include relative contributions from
00221             // SRFVelocityFvPatchVectorField's
00222             const SRFVelocityFvPatchVectorField& UrelPatch =
00223                 refCast<const SRFVelocityFvPatchVectorField>(bvf[i]);
00224             if (UrelPatch.relative())
00225             {
00226                 tUabs().boundaryField()[i] += Urel_.boundaryField()[i];
00227             }
00228         }
00229         else
00230         {
00231             tUabs().boundaryField()[i] += Urel_.boundaryField()[i];
00232         }
00233     }
00234 
00235     return tUabs;
00236 }
00237 
00238 
00239 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines