Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
applications
utilities
postProcessing
velocityField
Lambda2
Lambda2.C
Go to the documentation of this file.
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.
10
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.
15
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.
20
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/>.
23
24
Application
25
Lambda2
26
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.
31
32
The -noWrite option has no meaning.
33
34
Usage
35
36
- Lambda2 [OPTIONS]
37
38
@param -dict <dictionary name>\n
39
Use named dictionary instead of system/controlDict.
40
41
@param -noZero \n
42
Ignore timestep 0.
43
44
@param -constant \n
45
Include the constant directory.
46
47
@param -time <time>\n
48
Apply only to specific time.
49
50
@param -latestTime \n
51
Only apply to latest time step.
52
53
@param -case <dir>\n
54
Case directory.
55
56
@param -parallel \n
57
Run in parallel.
58
59
@param -help \n
60
Display help message.
61
62
@param -doc \n
63
Display Doxygen API documentation page for this application.
64
65
@param -srcDoc \n
66
Display Doxygen source documentation page for this application.
67
68
\*---------------------------------------------------------------------------*/
69
70
#include <
postCalc/calc.H
>
71
#include <
finiteVolume/fvc.H
>
72
73
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74
75
void
Foam::calc
(
const
argList&
args
,
const
Time& runTime,
const
fvMesh&
mesh
)
76
{
77
IOobject Uheader
78
(
79
"U"
,
80
runTime.timeName(),
81
mesh
,
82
IOobject::MUST_READ
83
);
84
85
if
(Uheader.headerOk())
86
{
87
Info
<<
" Reading U"
<<
endl
;
88
volVectorField
U
(Uheader, mesh);
89
90
const
volTensorField
gradU(
fvc::grad
(
U
));
91
92
volTensorField
SSplusWW =
93
(
symm
(gradU) &
symm
(gradU)) + (
skew
(gradU) &
skew
(gradU));
94
95
volScalarField
Lambda2
96
(
97
IOobject
98
(
99
"Lambda2"
,
100
runTime.timeName(),
101
mesh
,
102
IOobject::NO_READ,
103
IOobject::NO_WRITE
104
),
105
-
eigenValues
(SSplusWW)().
component
(
vector::Y
)
106
);
107
108
Info
<<
" Writing -Lambda2"
<<
endl
;
109
Lambda2.
write
();
110
}
111
else
112
{
113
Info
<<
" No U"
<<
endl
;
114
}
115
116
Info
<<
"\nEnd\n"
<<
endl
;
117
}
118
119
120
// ************************ vim: set sw=4 sts=4 et: ************************ //