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
Q
Q.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
Q
26
27
Description
28
Calculates and writes the second invariant of the velocity gradient tensor.
29
30
The -noWrite option just outputs the max/min values without writing
31
the field.
32
33
Usage
34
35
- Q [OPTIONS]
36
37
@param -noWrite \n
38
Suppress output to files.
39
40
@param -dict <dictionary name>\n
41
Use named dictionary instead of system/controlDict.
42
43
@param -noZero \n
44
Ignore timestep 0.
45
46
@param -constant \n
47
Include the constant directory.
48
49
@param -time <time>\n
50
Apply only to specific time.
51
52
@param -latestTime \n
53
Only apply to latest time step.
54
55
@param -case <dir>\n
56
Case directory.
57
58
@param -parallel \n
59
Run in parallel.
60
61
@param -help \n
62
Display help message.
63
64
@param -doc \n
65
Display Doxygen API documentation page for this application.
66
67
@param -srcDoc \n
68
Display Doxygen source documentation page for this application.
69
70
\*---------------------------------------------------------------------------*/
71
72
#include <
postCalc/calc.H
>
73
#include <
finiteVolume/fvc.H
>
74
75
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76
77
void
Foam::calc
(
const
argList&
args
,
const
Time& runTime,
const
fvMesh&
mesh
)
78
{
79
bool
writeResults = !args.optionFound(
"noWrite"
);
80
81
IOobject Uheader
82
(
83
"U"
,
84
runTime.timeName(),
85
mesh
,
86
IOobject::MUST_READ
87
);
88
89
if
(Uheader.headerOk())
90
{
91
Info
<<
" Reading U"
<<
endl
;
92
volVectorField
U
(Uheader, mesh);
93
volTensorField
gradU =
fvc::grad
(
U
);
94
95
volScalarField
Q
96
(
97
IOobject
98
(
99
"Q"
,
100
runTime.timeName(),
101
mesh
,
102
IOobject::NO_READ,
103
IOobject::NO_WRITE
104
),
105
0.5*(
sqr
(
tr
(gradU)) -
tr
(((gradU)&(gradU))))
106
);
107
108
/*
109
// This is a second way of calculating Q, that delivers results
110
// very close, but not identical to the first approach.
111
112
volSymmTensorField S = symm(gradU); // symmetric part of tensor
113
volTensorField W = skew(gradU); // anti-symmetric part
114
115
volScalarField SS = S&&S;
116
volScalarField WW = W&&W;
117
118
volScalarField Q
119
(
120
IOobject
121
(
122
"Q",
123
runTime.timeName(),
124
mesh,
125
IOobject::NO_READ,
126
IOobject::NO_WRITE
127
),
128
0.5*(WW - SS)
129
);
130
*/
131
132
Info
<<
"mag(Q) max/min : "
133
<<
max
(Q).value() <<
" "
134
<<
min
(Q).value() <<
endl
;
135
136
if
(writeResults)
137
{
138
Q.
write
();
139
}
140
}
141
else
142
{
143
Info
<<
" No U"
<<
endl
;
144
}
145
146
Info
<<
"\nEnd\n"
<<
endl
;
147
}
148
149
150
// ************************ vim: set sw=4 sts=4 et: ************************ //