initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / postProcessing / velocityField / Co / Co.C
blobd69bc37a7987e04b351ef0b09a95111db561ac7d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     Co
28 Description
29     Calculates and writes the Co number as a surfaceScalarField obtained
30     from field phi.
31     The -noWrite option just outputs the max values without writing the
32     field.
34 \*---------------------------------------------------------------------------*/
36 #include "calc.H"
37 #include "fvc.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
43     bool writeResults = !args.options().found("noWrite");
45     IOobject phiHeader
46     (
47         "phi",
48         runTime.timeName(),
49         mesh,
50         IOobject::MUST_READ
51     );
53     if (phiHeader.headerOk())
54     {
55         autoPtr<surfaceScalarField> CoPtr;
57         Info<< "    Reading phi" << endl;
58         surfaceScalarField phi(phiHeader, mesh);
59         Info<< "    Calculating Co" << endl;
61         if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
62         {
63             // compressible
64             volScalarField rho
65             (
66                 IOobject
67                 (
68                     "rho",
69                     runTime.timeName(),
70                     mesh,
71                     IOobject::MUST_READ
72                 ),
73                 mesh
74             );
76             CoPtr.set
77             (
78                 new surfaceScalarField
79                 (
80                     IOobject
81                     (
82                         "Co",
83                         runTime.timeName(),
84                         mesh,
85                         IOobject::NO_READ
86                     ),
87                     (
88                         mesh.surfaceInterpolation::deltaCoeffs()
89                       * (mag(phi)/(fvc::interpolate(rho)*mesh.magSf()))
90                       * runTime.deltaT()
91                     )
92                 )
93             );
94         }
95         else if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
96         {
97             // incompressible
98             CoPtr.set
99             (
100                 new surfaceScalarField
101                 (
102                     IOobject
103                     (
104                         "Co",
105                         runTime.timeName(),
106                         mesh,
107                         IOobject::NO_READ
108                     ),
109                     (
110                         mesh.surfaceInterpolation::deltaCoeffs()
111                       * (mag(phi)/mesh.magSf())
112                       * runTime.deltaT()
113                     )
114                 )
115             );
116         }
117         else
118         {
119             FatalErrorIn(args.executable())
120                 << "Incorrect dimensions of phi: " << phi.dimensions()
121                     << abort(FatalError);
122         }
124         Info<< "Co max : " << max(CoPtr()).value() << endl;
126         if (writeResults)
127         {
128             CoPtr().write();
129         }
130     }
131     else
132     {
133         Info<< "    No phi" << endl;
134     }
137 // ************************************************************************* //