Co.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-2015 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  Co
26 
27 Description
28  Calculates and writes the Co number as a volScalarField obtained
29  from field phi.
30 
31  The -noWrite option just outputs the max values without writing the
32  field.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "calc.H"
37 #include "fvc.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
42 {
43  bool writeResults = !args.optionFound("noWrite");
44 
45  IOobject phiHeader
46  (
47  "phi",
48  runTime.timeName(),
49  mesh,
50  IOobject::MUST_READ
51  );
52 
53  if (phiHeader.headerOk())
54  {
56  (
57  IOobject
58  (
59  "Co",
60  runTime.timeName(),
61  mesh,
62  IOobject::NO_READ
63  ),
64  mesh,
65  dimensionedScalar("0", dimless, 0),
66  zeroGradientFvPatchScalarField::typeName
67  );
68 
69  Info<< " Reading phi" << endl;
70  surfaceScalarField phi(phiHeader, mesh);
71 
72  if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
73  {
74  Info<< " Calculating compressible Co" << endl;
75 
76  Info<< " Reading rho" << endl;
78  (
79  IOobject
80  (
81  "rho",
82  runTime.timeName(),
83  mesh,
84  IOobject::MUST_READ
85  ),
86  mesh
87  );
88 
89  Co.dimensionedInternalField() =
90  (0.5*runTime.deltaT())
91  *fvc::surfaceSum(mag(phi))().dimensionedInternalField()
92  /(rho*mesh.V());
93  Co.correctBoundaryConditions();
94  }
95  else if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
96  {
97  Info<< " Calculating incompressible Co" << endl;
98 
99  Co.dimensionedInternalField() =
100  (0.5*runTime.deltaT())
101  *fvc::surfaceSum(mag(phi))().dimensionedInternalField()
102  /mesh.V();
103  Co.correctBoundaryConditions();
104  }
105  else
106  {
108  << "Incorrect dimensions of phi: " << phi.dimensions()
109  << abort(FatalError);
110  }
111 
112  Info<< "Co max : " << max(Co).value() << endl;
113 
114  if (writeResults)
115  {
116  Co.write();
117  }
118  }
119  else
120  {
121  Info<< " No phi" << endl;
122  }
123 
124  Info<< "\nEnd\n" << endl;
125 }
126 
127 
128 // ************************************************************************* //
Foam::fvc::surfaceSum
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcSurfaceIntegrate.C:138
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
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 > &)
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::Info
messageStream Info
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
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
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
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)