ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / applications / utilities / postProcessing / velocityField / Lambda2 / Lambda2.C
blob1c2b6a2c4a09d3be4eb35a032e8fc323691864ac
1 /*---------------------------------------------------------------------------* \
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Application
25     Lambda2
27 Description
28     Calculates and writes the second largest eigenvalue of the sum of the
29     square of the symmetrical and anti-symmetrical parts of the velocity
30     gradient tensor.
32     The -noWrite option has no meaning.
34 \*---------------------------------------------------------------------------*/
36 #include "calc.H"
37 #include "fvc.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
43     IOobject Uheader
44     (
45         "U",
46         runTime.timeName(),
47         mesh,
48         IOobject::MUST_READ
49     );
51     if (Uheader.headerOk())
52     {
53         Info<< "    Reading U" << endl;
54         volVectorField U(Uheader, mesh);
56         const volTensorField gradU(fvc::grad(U));
58         volTensorField SSplusWW =
59             (symm(gradU) & symm(gradU)) + (skew(gradU) & skew(gradU));
61         volScalarField Lambda2
62         (
63             IOobject
64             (
65                 "Lambda2",
66                 runTime.timeName(),
67                 mesh,
68                 IOobject::NO_READ,
69                 IOobject::NO_WRITE
70             ),
71             -eigenValues(SSplusWW)().component(vector::Y)
72         );
74         Info << "    Writing -Lambda2" << endl;
75         Lambda2.write();
76     }
77     else
78     {
79         Info<< "    No U" << endl;
80     }
82     Info<< "\nEnd\n" << endl;
86 // ************************************************************************* //