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
flowType
flowType.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
flowType
26
27
Description
28
Calculates and writes the flowType of velocity field U.
29
30
The -noWrite option has no meaning.
31
32
The flow type parameter is obtained according to the following equation:
33
@verbatim
34
|D| - |Omega|
35
lambda = -------------
36
|D| + |Omega|
37
38
-1 = rotational flow
39
0 = simple shear flow
40
1 = planar extensional flow
41
@endverbatim
42
43
Usage
44
45
- flowType [OPTIONS]
46
47
@param -dict <dictionary name>\n
48
Use named dictionary instead of system/controlDict.
49
50
@param -noZero \n
51
Ignore timestep 0.
52
53
@param -constant \n
54
Include the constant directory.
55
56
@param -time <time>\n
57
Apply only to specific time.
58
59
@param -latestTime \n
60
Only apply to latest time step.
61
62
@param -case <dir>\n
63
Case directory.
64
65
@param -parallel \n
66
Run in parallel.
67
68
@param -help \n
69
Display help message.
70
71
@param -doc \n
72
Display Doxygen API documentation page for this application.
73
74
@param -srcDoc \n
75
Display Doxygen source documentation page for this application.
76
77
\*---------------------------------------------------------------------------*/
78
79
#include <
postCalc/calc.H
>
80
#include <
finiteVolume/fvc.H
>
81
82
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83
84
void
Foam::calc
(
const
argList&
args
,
const
Time& runTime,
const
fvMesh&
mesh
)
85
{
86
IOobject Uheader
87
(
88
"U"
,
89
runTime.timeName(),
90
mesh
,
91
IOobject::MUST_READ
92
);
93
94
if
(Uheader.headerOk())
95
{
96
Info
<<
" Reading U"
<<
endl
;
97
volVectorField
U
(Uheader, mesh);
98
99
volTensorField
gradU =
fvc::grad
(
U
);
100
volScalarField
magD =
mag
(
symm
(gradU));
101
volScalarField
magOmega =
mag
(
skew
(gradU));
102
dimensionedScalar
smallMagD(
"smallMagD"
, magD.dimensions(), SMALL);
103
104
Info
<<
" Calculating flowType"
<<
endl
;
105
106
volScalarField
flowType
107
(
108
IOobject
109
(
110
"flowType"
,
111
runTime.timeName(),
112
mesh
,
113
IOobject::NO_READ
114
),
115
(magD - magOmega)/(magD + magOmega + smallMagD)
116
);
117
118
flowType.write();
119
}
120
else
121
{
122
Info
<<
" No U"
<<
endl
;
123
}
124
125
Info
<<
"\nEnd\n"
<<
endl
;
126
}
127
128
129
// ************************ vim: set sw=4 sts=4 et: ************************ //