ptot.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  ptot
26 
27 Description
28  For each time: calculate the total pressure.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 int main(int argc, char *argv[])
37 {
38  timeSelector::addOptions();
39  #include "addRegionOption.H"
40 
41  #include "setRootCase.H"
42  #include "createTime.H"
43 
44  instantList timeDirs = timeSelector::select0(runTime, args);
45 
46  #include "createNamedMesh.H"
47 
48  forAll(timeDirs, timeI)
49  {
50  runTime.setTime(timeDirs[timeI], timeI);
51 
52  Info<< "Time = " << runTime.timeName() << endl;
53 
54  IOobject pheader
55  (
56  "p",
57  runTime.timeName(),
58  mesh,
59  IOobject::MUST_READ
60  );
61 
62  IOobject Uheader
63  (
64  "U",
65  runTime.timeName(),
66  mesh,
67  IOobject::MUST_READ
68  );
69 
70 
71  // Check p and U exist
72  if (pheader.headerOk() && Uheader.headerOk())
73  {
74  mesh.readUpdate();
75 
76  Info<< " Reading p" << endl;
77  volScalarField p(pheader, mesh);
78 
79  Info<< " Reading U" << endl;
80  volVectorField U(Uheader, mesh);
81 
82  Info<< " Calculating ptot" << endl;
83  if (p.dimensions() == dimensionSet(0, 2, -2, 0, 0))
84  {
85  volScalarField ptot
86  (
87  IOobject
88  (
89  "ptot",
90  runTime.timeName(),
91  mesh,
92  IOobject::NO_READ
93  ),
94  p + 0.5*magSqr(U)
95  );
96  ptot.write();
97  }
98  else
99  {
100  IOobject rhoheader
101  (
102  "rho",
103  runTime.timeName(),
104  mesh,
105  IOobject::MUST_READ
106  );
107 
108  // Check rho exists
109  if (rhoheader.headerOk())
110  {
111  Info<< " Reading rho" << endl;
112  volScalarField rho(rhoheader, mesh);
113 
114  volScalarField ptot
115  (
116  IOobject
117  (
118  "ptot",
119  runTime.timeName(),
120  mesh,
121  IOobject::NO_READ
122  ),
123  p + 0.5*rho*magSqr(U)
124  );
125  ptot.write();
126  }
127  else
128  {
129  Info<< " No rho" << endl;
130  }
131  }
132  }
133  else
134  {
135  Info<< " No p or U" << endl;
136  }
137 
138  Info<< endl;
139  }
140 
141  Info<< "End" << nl << endl;
142 
143  return 0;
144 }
145 
146 
147 // ************************************************************************* //
p
p
Definition: pEqn.H:62
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:42
U
U
Definition: pEqn.H:46
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
addRegionOption.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
createTime.H
fvCFD.H
args
Foam::argList args(argc, argv)
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)