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

dsmcFields.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) 2009-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 "dsmcFields.H"
00027 #include <finiteVolume/volFields.H>
00028 #include <OpenFOAM/dictionary.H>
00029 #include <dsmc/dsmcCloud.H>
00030 
00031 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035     defineTypeNameAndDebug(dsmcFields, 0);
00036 }
00037 
00038 
00039 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00040 
00041 Foam::dsmcFields::dsmcFields
00042 (
00043     const word& name,
00044     const objectRegistry& obr,
00045     const dictionary& dict,
00046     const bool loadFromFiles
00047 )
00048 :
00049     name_(name),
00050     obr_(obr),
00051     active_(true)
00052 {
00053     // Check if the available mesh is an fvMesh, otherwise deactivate
00054     if (!isA<fvMesh>(obr_))
00055     {
00056         active_ = false;
00057         WarningIn
00058         (
00059             "dsmcFields::dsmcFields"
00060             "(const objectRegistry&, const dictionary&)"
00061         )   << "No fvMesh available, deactivating." << nl
00062             << endl;
00063     }
00064 
00065     read(dict);
00066 }
00067 
00068 
00069 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00070 
00071 Foam::dsmcFields::~dsmcFields()
00072 {}
00073 
00074 
00075 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00076 
00077 void Foam::dsmcFields::read(const dictionary& dict)
00078 {
00079     if (active_)
00080     {
00081 
00082     }
00083 }
00084 
00085 
00086 void Foam::dsmcFields::execute()
00087 {
00088     // Do nothing - only valid on write
00089 }
00090 
00091 
00092 void Foam::dsmcFields::end()
00093 {
00094     // Do nothing - only valid on write
00095 }
00096 
00097 
00098 void Foam::dsmcFields::write()
00099 {
00100     if (active_)
00101     {
00102         word rhoNMeanName = "rhoNMean";
00103         word rhoMMeanName = "rhoMMean";
00104         word momentumMeanName = "momentumMean";
00105         word linearKEMeanName = "linearKEMean";
00106         word internalEMeanName = "internalEMean";
00107         word iDofMeanName = "iDofMean";
00108         word fDMeanName = "fDMean";
00109 
00110         const volScalarField& rhoNMean = obr_.lookupObject<volScalarField>
00111         (
00112             rhoNMeanName
00113         );
00114 
00115         const volScalarField& rhoMMean = obr_.lookupObject<volScalarField>
00116         (
00117             rhoMMeanName
00118         );
00119 
00120         const volVectorField& momentumMean = obr_.lookupObject<volVectorField>
00121         (
00122             momentumMeanName
00123         );
00124 
00125         const volScalarField& linearKEMean = obr_.lookupObject<volScalarField>
00126         (
00127             linearKEMeanName
00128         );
00129 
00130         const volScalarField& internalEMean = obr_.lookupObject<volScalarField>
00131         (
00132             internalEMeanName
00133         );
00134 
00135         const volScalarField& iDofMean = obr_.lookupObject<volScalarField>
00136         (
00137             iDofMeanName
00138         );
00139 
00140         volVectorField fDMean =  obr_.lookupObject<volVectorField>
00141         (
00142             fDMeanName
00143         );
00144 
00145         if (min(mag(rhoNMean)).value() > VSMALL)
00146         {
00147             Info<< "Calculating dsmcFields." << endl;
00148 
00149             Info<< "    Calculating UMean field." << endl;
00150             volVectorField UMean
00151             (
00152                 IOobject
00153                 (
00154                     "UMean",
00155                     obr_.time().timeName(),
00156                     obr_,
00157                     IOobject::NO_READ
00158                 ),
00159                 momentumMean/rhoMMean
00160             );
00161 
00162             Info<< "    Calculating translationalT field." << endl;
00163             volScalarField translationalT
00164             (
00165                 IOobject
00166                 (
00167                     "translationalT",
00168                     obr_.time().timeName(),
00169                     obr_,
00170                     IOobject::NO_READ
00171                 ),
00172                 2.0/(3.0*dsmcCloud::kb*rhoNMean)
00173                 *(linearKEMean - 0.5*rhoMMean*(UMean & UMean))
00174             );
00175 
00176             Info<< "    Calculating internalT field." << endl;
00177             volScalarField internalT
00178             (
00179                 IOobject
00180                 (
00181                     "internalT",
00182                     obr_.time().timeName(),
00183                     obr_,
00184                     IOobject::NO_READ
00185                 ),
00186                 (2.0/dsmcCloud::kb)*(internalEMean/iDofMean)
00187             );
00188 
00189             Info<< "    Calculating overallT field." << endl;
00190             volScalarField overallT
00191             (
00192                 IOobject
00193                 (
00194                     "overallT",
00195                     obr_.time().timeName(),
00196                     obr_,
00197                     IOobject::NO_READ
00198                 ),
00199                 2.0/(dsmcCloud::kb*(3.0*rhoNMean + iDofMean))
00200                *(linearKEMean - 0.5*rhoMMean*(UMean & UMean) + internalEMean)
00201             );
00202 
00203             Info<< "    Calculating pressure field." << endl;
00204             volScalarField p
00205             (
00206                 IOobject
00207                 (
00208                     "p",
00209                     obr_.time().timeName(),
00210                     obr_,
00211                     IOobject::NO_READ
00212                 ),
00213                 dsmcCloud::kb*rhoNMean*translationalT
00214             );
00215 
00216             const fvMesh& mesh = fDMean.mesh();
00217 
00218             forAll(mesh.boundaryMesh(), i)
00219             {
00220                 const polyPatch& patch = mesh.boundaryMesh()[i];
00221 
00222                 if (isA<wallPolyPatch>(patch))
00223                 {
00224                     p.boundaryField()[i] =
00225                         fDMean.boundaryField()[i]
00226                       & (patch.faceAreas()/mag(patch.faceAreas()));
00227                 }
00228             }
00229 
00230             Info<< "    mag(UMean) max/min : "
00231                 << max(mag(UMean)).value() << " "
00232                 << min(mag(UMean)).value() << endl;
00233 
00234             Info<< "    translationalT max/min : "
00235                 << max(translationalT).value() << " "
00236                 << min(translationalT).value() << endl;
00237 
00238             Info<< "    internalT max/min : "
00239                 << max(internalT).value() << " "
00240                 << min(internalT).value() << endl;
00241 
00242             Info<< "    overallT max/min : "
00243                 << max(overallT).value() << " "
00244                 << min(overallT).value() << endl;
00245 
00246             Info<< "    p max/min : "
00247                 << max(p).value() << " "
00248                 << min(p).value() << endl;
00249 
00250             UMean.write();
00251 
00252             translationalT.write();
00253 
00254             internalT.write();
00255 
00256             overallT.write();
00257 
00258             p.write();
00259 
00260             Info<< "dsmcFields written." << nl << endl;
00261         }
00262         else
00263         {
00264             Info<< "Small value (" << min(mag(rhoNMean))
00265                 << ") found in rhoNMean field. "
00266                 << "Not calculating dsmcFields to avoid division by zero."
00267                 << endl;
00268         }
00269     }
00270 }
00271 
00272 
00273 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines