Mach.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  Mach
26 
27 Description
28  Calculates and optionally writes the local Mach number from the velocity
29  field U at each time.
30 
31  The -nowrite option just outputs the max value without writing the field.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "calc.H"
36 #include "fluidThermo.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  IOobject Theader
53  (
54  "T",
55  runTime.timeName(),
56  mesh,
57  IOobject::MUST_READ
58  );
59 
60  // Check U and T exists
61  if (Uheader.headerOk() && Theader.headerOk())
62  {
63  autoPtr<volScalarField> MachPtr;
64 
65  volVectorField U(Uheader, mesh);
66 
67  if
68  (
69  IOobject
70  (
72  runTime.constant(),
73  mesh
74  ).headerOk()
75  )
76  {
77  // thermophysical Mach
78  autoPtr<fluidThermo> thermo
79  (
81  );
82 
83  volScalarField Cp(thermo->Cp());
84  volScalarField Cv(thermo->Cv());
85 
86  MachPtr.set
87  (
88  new volScalarField
89  (
90  IOobject
91  (
92  "Ma",
93  runTime.timeName(),
94  mesh
95  ),
96  mag(U)/(sqrt((Cp/Cv)*(Cp - Cv)*thermo->T()))
97  )
98  );
99  }
100  else
101  {
102  // thermodynamic Mach
103  IOdictionary thermoProps
104  (
105  IOobject
106  (
107  "thermodynamicProperties",
108  runTime.constant(),
109  mesh,
110  IOobject::MUST_READ_IF_MODIFIED,
111  IOobject::NO_WRITE
112  )
113  );
114 
115  dimensionedScalar R(thermoProps.lookup("R"));
116  dimensionedScalar Cv(thermoProps.lookup("Cv"));
117 
118  volScalarField T(Theader, mesh);
119 
120  MachPtr.set
121  (
122  new volScalarField
123  (
124  IOobject
125  (
126  "Ma",
127  runTime.timeName(),
128  mesh
129  ),
130  mag(U)/(sqrt(((Cv + R)/Cv)*R*T))
131  )
132  );
133  }
134 
135  Info<< "Mach max : " << max(MachPtr()).value() << endl;
136 
137  if (writeResults)
138  {
139  MachPtr().write();
140  }
141  }
142  else
143  {
144  Info<< " Missing U or T" << endl;
145  }
146 
147  Info<< nl << "End" << nl << endl;
148 }
149 
150 
151 // ************************************************************************* //
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
fluidThermo.H
Foam::calc
void calc(const argList &args, const Time &runTime, const fvMesh &mesh)
Definition: Lambda2.C:38
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
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 > &)
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
U
U
Definition: pEqn.H:46
dictName
const word dictName("particleTrackDict")
R
#define R(A, B, C, D, E, F, K, M)
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
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
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 > &)
T
const volScalarField & T
Definition: createFields.H:25
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
args
Foam::argList args(argc, argv)