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

Co.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 Application
00025     Co
00026 
00027 Description
00028     Calculates and writes the Co number as a surfaceScalarField obtained
00029     from field phi.
00030 
00031     The -noWrite option just outputs the max values without writing the
00032     field.
00033 
00034 Usage
00035 
00036     - Co [OPTIONS]
00037 
00038     @param -noWrite \n
00039     Suppress output to files.
00040 
00041     @param -dict <dictionary name>\n
00042     Use named dictionary instead of system/controlDict.
00043 
00044     @param -noZero \n
00045     Ignore timestep 0.
00046 
00047     @param -constant \n
00048     Include the constant directory.
00049 
00050     @param -time <time>\n
00051     Apply only to specific time.
00052 
00053     @param -latestTime \n
00054     Only apply to latest time step.
00055 
00056     @param -case <dir>\n
00057     Case directory.
00058 
00059     @param -parallel \n
00060     Run in parallel.
00061 
00062     @param -help \n
00063     Display help message.
00064 
00065     @param -doc \n
00066     Display Doxygen API documentation page for this application.
00067 
00068     @param -srcDoc \n
00069     Display Doxygen source documentation page for this application.
00070 
00071 \*---------------------------------------------------------------------------*/
00072 
00073 #include <postCalc/calc.H>
00074 #include <finiteVolume/fvc.H>
00075 
00076 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00077 
00078 namespace Foam
00079 {
00080     tmp<volScalarField> Co(const surfaceScalarField& Cof)
00081     {
00082         const fvMesh& mesh = Cof.mesh();
00083 
00084         tmp<volScalarField> tCo
00085         (
00086             new volScalarField
00087             (
00088                 IOobject
00089                 (
00090                     "Co",
00091                     mesh.time().timeName(),
00092                     mesh
00093                 ),
00094                 mesh,
00095                 dimensionedScalar("0", Cof.dimensions(), 0)
00096             )
00097         );
00098 
00099         volScalarField& Co = tCo();
00100 
00101         // Set local references to mesh data
00102         const unallocLabelList& owner = mesh.owner();
00103         const unallocLabelList& neighbour = mesh.neighbour();
00104 
00105         forAll(owner, facei)
00106         {
00107             label own = owner[facei];
00108             label nei = neighbour[facei];
00109 
00110             Co[own] = max(Co[own], Cof[facei]);
00111             Co[nei] = max(Co[nei], Cof[facei]);
00112         }
00113 
00114         forAll(Co.boundaryField(), patchi)
00115         {
00116             Co.boundaryField()[patchi] = Cof.boundaryField()[patchi];
00117         }
00118 
00119         return tCo;
00120     }
00121 }
00122 
00123 
00124 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
00125 {
00126     bool writeResults = !args.optionFound("noWrite");
00127 
00128     IOobject phiHeader
00129     (
00130         "phi",
00131         runTime.timeName(),
00132         mesh,
00133         IOobject::MUST_READ
00134     );
00135 
00136     if (phiHeader.headerOk())
00137     {
00138         autoPtr<surfaceScalarField> CoPtr;
00139 
00140         Info<< "    Reading phi" << endl;
00141         surfaceScalarField phi(phiHeader, mesh);
00142         Info<< "    Calculating Co" << endl;
00143 
00144         if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
00145         {
00146             // compressible
00147             volScalarField rho
00148             (
00149                 IOobject
00150                 (
00151                     "rho",
00152                     runTime.timeName(),
00153                     mesh,
00154                     IOobject::MUST_READ
00155                 ),
00156                 mesh
00157             );
00158 
00159             CoPtr.set
00160             (
00161                 new surfaceScalarField
00162                 (
00163                     IOobject
00164                     (
00165                         "Cof",
00166                         runTime.timeName(),
00167                         mesh,
00168                         IOobject::NO_READ
00169                     ),
00170                     (
00171                         mesh.surfaceInterpolation::deltaCoeffs()
00172                       * (mag(phi)/(fvc::interpolate(rho)*mesh.magSf()))
00173                       * runTime.deltaT()
00174                     )
00175                 )
00176             );
00177         }
00178         else if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
00179         {
00180             // incompressible
00181             CoPtr.set
00182             (
00183                 new surfaceScalarField
00184                 (
00185                     IOobject
00186                     (
00187                         "Cof",
00188                         runTime.timeName(),
00189                         mesh,
00190                         IOobject::NO_READ
00191                     ),
00192                     (
00193                         mesh.surfaceInterpolation::deltaCoeffs()
00194                       * (mag(phi)/mesh.magSf())
00195                       * runTime.deltaT()
00196                     )
00197                 )
00198             );
00199         }
00200         else
00201         {
00202             FatalErrorIn(args.executable())
00203                 << "Incorrect dimensions of phi: " << phi.dimensions()
00204                     << abort(FatalError);
00205         }
00206 
00207         Info<< "Co max : " << max(CoPtr()).value() << endl;
00208 
00209         if (writeResults)
00210         {
00211             CoPtr().write();
00212             Co(CoPtr())().write();
00213         }
00214     }
00215     else
00216     {
00217         Info<< "    No phi" << endl;
00218     }
00219 
00220     Info<< "\nEnd\n" << endl;
00221 }
00222 
00223 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines