Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / postProcessing / velocityField / Lambda2 / Lambda2.C
blob8f04aa681c26963c4ed3a6a5c463a5679d4ea478
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  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));
60 //            (skew(gradU) & skew(gradU));
62         volScalarField Lambda2
63         (
64             IOobject
65             (
66                 "Lambda2",
67                 runTime.timeName(),
68                 mesh,
69                 IOobject::NO_READ,
70                 IOobject::NO_WRITE
71             ),
72             -eigenValues(SSplusWW)().component(vector::Y)
73         );
75         Info << "    Writing -Lambda2" << endl;
76         Lambda2.write();
77     }
78     else
79     {
80         Info<< "    No U" << endl;
81     }
83     Info<< "\nEnd\n" << endl;
87 // ************************************************************************* //