vorticity.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) 2011 OpenFOAM Foundation
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  vorticity
26 
27 Description
28  Calculates and writes the vorticity of velocity field U.
29 
30  The -noWrite option just outputs the max/min values without writing
31  the field.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "calc.H"
36 #include "fvc.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
41 {
42  bool writeResults = !args.optionFound("noWrite");
43 
44  IOobject Uheader
45  (
46  "U",
47  runTime.timeName(),
48  mesh,
49  IOobject::MUST_READ
50  );
51 
52  if (Uheader.headerOk())
53  {
54  Info<< " Reading U" << endl;
55  volVectorField U(Uheader, mesh);
56 
57  Info<< " Calculating vorticity" << endl;
58  volVectorField vorticity
59  (
60  IOobject
61  (
62  "vorticity",
63  runTime.timeName(),
64  mesh,
65  IOobject::NO_READ
66  ),
67  fvc::curl(U)
68  );
69 
70  volScalarField magVorticity
71  (
72  IOobject
73  (
74  "magVorticity",
75  runTime.timeName(),
76  mesh,
77  IOobject::NO_READ
78  ),
79  mag(vorticity)
80  );
81 
82  Info<< "vorticity max/min : "
83  << max(magVorticity).value() << " "
84  << min(magVorticity).value() << endl;
85 
86  if (writeResults)
87  {
88  vorticity.write();
89  magVorticity.write();
90  }
91  }
92  else
93  {
94  Info<< " No U" << endl;
95  }
96 
97  Info<< "\nEnd\n" << endl;
98 }
99 
100 
101 // ************************************************************************* //
Foam::fvc::curl
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcCurl.C:45
Foam::calc
void calc(const argList &args, const Time &runTime, const fvMesh &mesh)
Definition: Lambda2.C:38
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
U
U
Definition: pEqn.H:46
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::Info
messageStream Info
calc.H
Function prototype for all simple post-processing functions e.g. calcDivPhi, calcMagU etc.
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
fvc.H
args
Foam::argList args(argc, argv)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)