initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / velocityField / Q / Q.C
bloba5cb4ed5a54ef284bb23fe1dc798b0cfa3bfa2d4
1 /*---------------------------------------------------------------------------* \
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     Q
28 Description
29     Calculates and writes the second invariant of the velocity gradient tensor.
31     The -noWrite option just outputs the max/min values without writing
32     the field.
34 \*---------------------------------------------------------------------------*/
36 #include "calc.H"
37 #include "fvc.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
43     bool writeResults = !args.optionFound("noWrite");
45     IOobject Uheader
46     (
47         "U",
48         runTime.timeName(),
49         mesh,
50         IOobject::MUST_READ
51     );
53     if (Uheader.headerOk())
54     {
55         Info<< "    Reading U" << endl;
56         volVectorField U(Uheader, mesh);
57         volTensorField gradU = fvc::grad(U);
59         volScalarField Q
60         (
61             IOobject
62             (
63                 "Q",
64                 runTime.timeName(),
65                 mesh,
66                 IOobject::NO_READ,
67                 IOobject::NO_WRITE
68             ),
69             0.5*(sqr(tr(gradU)) - tr(((gradU)&(gradU))))
70         );
72         /*
73         // This is a second way of calculating Q, that delivers results
74         // very close, but not identical to the first approach.
76         volSymmTensorField S = symm(gradU);  // symmetric part of tensor
77         volTensorField W = skew(gradU);  // anti-symmetric part
79         volScalarField SS =  S&&S;
80         volScalarField WW =  W&&W;
82         volScalarField Q
83         (
84             IOobject
85             (
86                 "Q",
87                 runTime.timeName(),
88                 mesh,
89                 IOobject::NO_READ,
90                 IOobject::NO_WRITE
91             ),
92             0.5*(WW - SS)
93         );
94         */
96         Info<< "mag(Q) max/min : "
97             << max(Q).value() << " "
98             << min(Q).value() << endl;
100         if (writeResults)
101         {
102             Q.write();
103         }
104     }
105     else
106     {
107         Info<< "    No U" << endl;
108     }
110     Info<< "\nEnd\n" << endl;
114 // ************************************************************************* //