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

flowType.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     flowType
00026 
00027 Description
00028     Calculates and writes the flowType of velocity field U.
00029 
00030     The -noWrite option has no meaning.
00031 
00032     The flow type parameter is obtained according to the following equation:
00033     @verbatim
00034                  |D| - |Omega|
00035         lambda = -------------
00036                  |D| + |Omega|
00037 
00038         -1 = rotational flow
00039          0 = simple shear flow
00040          1 = planar extensional flow
00041     @endverbatim
00042 
00043 Usage
00044 
00045     - flowType [OPTIONS]
00046 
00047     @param -dict <dictionary name>\n
00048     Use named dictionary instead of system/controlDict.
00049 
00050     @param -noZero \n
00051     Ignore timestep 0.
00052 
00053     @param -constant \n
00054     Include the constant directory.
00055 
00056     @param -time <time>\n
00057     Apply only to specific time.
00058 
00059     @param -latestTime \n
00060     Only apply to latest time step.
00061 
00062     @param -case <dir>\n
00063     Case directory.
00064 
00065     @param -parallel \n
00066     Run in parallel.
00067 
00068     @param -help \n
00069     Display help message.
00070 
00071     @param -doc \n
00072     Display Doxygen API documentation page for this application.
00073 
00074     @param -srcDoc \n
00075     Display Doxygen source documentation page for this application.
00076 
00077 \*---------------------------------------------------------------------------*/
00078 
00079 #include <postCalc/calc.H>
00080 #include <finiteVolume/fvc.H>
00081 
00082 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00083 
00084 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
00085 {
00086     IOobject Uheader
00087     (
00088         "U",
00089         runTime.timeName(),
00090         mesh,
00091         IOobject::MUST_READ
00092     );
00093 
00094     if (Uheader.headerOk())
00095     {
00096         Info<< "    Reading U" << endl;
00097         volVectorField U(Uheader, mesh);
00098 
00099         volTensorField gradU = fvc::grad(U);
00100         volScalarField magD = mag(symm(gradU));
00101         volScalarField magOmega = mag(skew(gradU));
00102         dimensionedScalar smallMagD("smallMagD", magD.dimensions(), SMALL);
00103 
00104         Info<< "    Calculating flowType" << endl;
00105 
00106         volScalarField flowType
00107         (
00108             IOobject
00109             (
00110                 "flowType",
00111                 runTime.timeName(),
00112                 mesh,
00113                 IOobject::NO_READ
00114             ),
00115             (magD - magOmega)/(magD + magOmega + smallMagD)
00116         );
00117 
00118         flowType.write();
00119     }
00120     else
00121     {
00122         Info<< "    No U" << endl;
00123     }
00124 
00125     Info<< "\nEnd\n" << endl;
00126 }
00127 
00128 
00129 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines