initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / velocityField / Co / Co.C
blob2965d4727caafbfbe0a98ad4493f422af9e03c33
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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.
32     The -noWrite option just outputs the max values without writing the
33     field.
35 \*---------------------------------------------------------------------------*/
37 #include "calc.H"
38 #include "fvc.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 namespace Foam
44     tmp<volScalarField> Co(const surfaceScalarField& Cof)
45     {
46         const fvMesh& mesh = Cof.mesh();
48         tmp<volScalarField> tCo
49         (
50             new volScalarField
51             (
52                 IOobject
53                 (
54                     "Co",
55                     mesh.time().timeName(),
56                     mesh
57                 ),
58                 mesh,
59                 dimensionedScalar("0", Cof.dimensions(), 0)
60             )
61         );
63         volScalarField& Co = tCo();
65         // Set local references to mesh data
66         const unallocLabelList& owner = mesh.owner();
67         const unallocLabelList& neighbour = mesh.neighbour();
69         forAll(owner, facei)
70         {
71             label own = owner[facei];
72             label nei = neighbour[facei];
74             Co[own] = max(Co[own], Cof[facei]);
75             Co[nei] = max(Co[nei], Cof[facei]);
76         }
78         forAll(Co.boundaryField(), patchi)
79         {
80             Co.boundaryField()[patchi] = Cof.boundaryField()[patchi];
81         }
83         return tCo;
84     }
88 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
90     bool writeResults = !args.optionFound("noWrite");
92     IOobject phiHeader
93     (
94         "phi",
95         runTime.timeName(),
96         mesh,
97         IOobject::MUST_READ
98     );
100     if (phiHeader.headerOk())
101     {
102         autoPtr<surfaceScalarField> CoPtr;
104         Info<< "    Reading phi" << endl;
105         surfaceScalarField phi(phiHeader, mesh);
106         Info<< "    Calculating Co" << endl;
108         if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
109         {
110             // compressible
111             volScalarField rho
112             (
113                 IOobject
114                 (
115                     "rho",
116                     runTime.timeName(),
117                     mesh,
118                     IOobject::MUST_READ
119                 ),
120                 mesh
121             );
123             CoPtr.set
124             (
125                 new surfaceScalarField
126                 (
127                     IOobject
128                     (
129                         "Cof",
130                         runTime.timeName(),
131                         mesh,
132                         IOobject::NO_READ
133                     ),
134                     (
135                         mesh.surfaceInterpolation::deltaCoeffs()
136                       * (mag(phi)/(fvc::interpolate(rho)*mesh.magSf()))
137                       * runTime.deltaT()
138                     )
139                 )
140             );
141         }
142         else if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
143         {
144             // incompressible
145             CoPtr.set
146             (
147                 new surfaceScalarField
148                 (
149                     IOobject
150                     (
151                         "Cof",
152                         runTime.timeName(),
153                         mesh,
154                         IOobject::NO_READ
155                     ),
156                     (
157                         mesh.surfaceInterpolation::deltaCoeffs()
158                       * (mag(phi)/mesh.magSf())
159                       * runTime.deltaT()
160                     )
161                 )
162             );
163         }
164         else
165         {
166             FatalErrorIn(args.executable())
167                 << "Incorrect dimensions of phi: " << phi.dimensions()
168                     << abort(FatalError);
169         }
171         Info<< "Co max : " << max(CoPtr()).value() << endl;
173         if (writeResults)
174         {
175             CoPtr().write();
176             Co(CoPtr())().write();
177         }
178     }
179     else
180     {
181         Info<< "    No phi" << endl;
182     }
184     Info<< "\nEnd\n" << endl;
187 // ************************************************************************* //