BUGFIX: Illegal use of uninitialised value (backport)
[foam-extend-3.2.git] / applications / utilities / postProcessing / velocityField / Lambda2 / Lambda2.C
blob154b5dde9c362e8d7e566840c0c5672b1aeba353
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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     Lambda2
28 Description
29     Calculates and writes the second largest eigenvalue of the sum of the
30     square of the symmetrical and anti-symmetrical parts of the velocity
31     gradient tensor.
33     The -noWrite option has no meaning.
35 \*---------------------------------------------------------------------------*/
37 #include "calc.H"
38 #include "fvc.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
44     IOobject Uheader
45     (
46         "U",
47         runTime.timeName(),
48         mesh,
49         IOobject::MUST_READ
50     );
52     if (Uheader.headerOk())
53     {
54         Info<< "    Reading U" << endl;
55         volVectorField U(Uheader, mesh);
57         const volTensorField gradU(fvc::grad(U));
59         volTensorField SSplusWW =
60             (symm(gradU) & symm(gradU));
61 //            (skew(gradU) & skew(gradU));
63         volScalarField Lambda2
64         (
65             IOobject
66             (
67                 "Lambda2",
68                 runTime.timeName(),
69                 mesh,
70                 IOobject::NO_READ,
71                 IOobject::NO_WRITE
72             ),
73             -eigenValues(SSplusWW)().component(vector::Y)
74         );
76         Info << "    Writing -Lambda2" << endl;
77         Lambda2.write();
78     }
79     else
80     {
81         Info<< "    No U" << endl;
82     }
84     Info<< "\nEnd\n" << endl;
88 // ************************************************************************* //