Pe.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  Pe
26 
27 Description
28  Calculates the Peclet number Pe from the flux phi and writes the maximum
29  value, the surfaceScalarField Pef and volScalarField Pe.
30 
31  With the -noWrite option just outputs the max value without writing
32  the fields.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "calc.H"
37 #include "surfaceInterpolate.H"
38 #include "fvcAverage.H"
39 
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
47 {
48  bool writeResults = !args.optionFound("noWrite");
49 
50  IOobject phiHeader
51  (
52  "phi",
53  runTime.timeName(),
54  mesh,
55  IOobject::MUST_READ
56  );
57 
58  if (phiHeader.headerOk())
59  {
60  autoPtr<surfaceScalarField> PePtr;
61 
62  Info<< " Reading phi" << endl;
63  surfaceScalarField phi(phiHeader, mesh);
64 
66  (
67  IOobject
68  (
69  "U",
70  runTime.timeName(),
71  mesh,
72  IOobject::MUST_READ
73  ),
74  mesh
75  );
76 
77  IOobject turbulencePropertiesHeader
78  (
79  "turbulenceProperties",
80  runTime.constant(),
81  mesh,
82  IOobject::MUST_READ_IF_MODIFIED,
83  IOobject::NO_WRITE
84  );
85 
86  Info<< " Calculating Pe" << endl;
87 
88  if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
89  {
90  if (turbulencePropertiesHeader.headerOk())
91  {
92  singlePhaseTransportModel laminarTransport(U, phi);
93 
94  autoPtr<incompressible::turbulenceModel> turbulenceModel
95  (
97  (
98  U,
99  phi,
101  )
102  );
103 
104  PePtr.set
105  (
107  (
108  IOobject
109  (
110  "Pef",
111  runTime.timeName(),
112  mesh
113  ),
114  mag(phi)
115  /(
116  mesh.magSf()
117  * mesh.surfaceInterpolation::deltaCoeffs()
119  )
120  )
121  );
122  }
123  else
124  {
125  IOdictionary transportProperties
126  (
127  IOobject
128  (
129  "transportProperties",
130  runTime.constant(),
131  mesh,
132  IOobject::MUST_READ_IF_MODIFIED,
133  IOobject::NO_WRITE
134  )
135  );
136 
137  dimensionedScalar nu(transportProperties.lookup("nu"));
138 
139  PePtr.set
140  (
142  (
143  IOobject
144  (
145  "Pef",
146  runTime.timeName(),
147  mesh
148  ),
149  mag(phi)
150  /(
151  mesh.magSf()
152  * mesh.surfaceInterpolation::deltaCoeffs()
153  * nu
154  )
155  )
156  );
157  }
158  }
159  else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
160  {
161  if (turbulencePropertiesHeader.headerOk())
162  {
163  autoPtr<fluidThermo> thermo(fluidThermo::New(mesh));
164 
166  (
167  IOobject
168  (
169  "rho",
170  runTime.timeName(),
171  mesh
172  ),
173  thermo->rho()
174  );
175 
176  autoPtr<compressible::turbulenceModel> turbulenceModel
177  (
179  (
180  rho,
181  U,
182  phi,
183  thermo()
184  )
185  );
186 
187  PePtr.set
188  (
190  (
191  IOobject
192  (
193  "Pef",
194  runTime.timeName(),
195  mesh
196  ),
197  mag(phi)
198  /(
199  mesh.magSf()
200  * mesh.surfaceInterpolation::deltaCoeffs()
202  )
203  )
204  );
205  }
206  else
207  {
208  IOdictionary transportProperties
209  (
210  IOobject
211  (
212  "transportProperties",
213  runTime.constant(),
214  mesh,
215  IOobject::MUST_READ_IF_MODIFIED,
216  IOobject::NO_WRITE
217  )
218  );
219 
220  dimensionedScalar mu(transportProperties.lookup("mu"));
221 
222  PePtr.set
223  (
225  (
226  IOobject
227  (
228  "Pef",
229  runTime.timeName(),
230  mesh
231  ),
232  mag(phi)
233  /(
234  mesh.magSf()
235  * mesh.surfaceInterpolation::deltaCoeffs()
236  * mu
237  )
238  )
239  );
240  }
241  }
242  else
243  {
245  << "Incorrect dimensions of phi: " << phi.dimensions()
246  << abort(FatalError);
247  }
248 
249  Info<< " Pe max : " << max(PePtr()).value() << endl;
250 
251  if (writeResults)
252  {
253  Info<< " Writing surfaceScalarField : "
254  << PePtr().name() << endl;
255  PePtr().write();
256 
257  volScalarField Pe
258  (
259  IOobject
260  (
261  "Pe",
262  runTime.timeName(),
263  mesh
264  ),
265  fvc::average(PePtr())
266  );
267 
268  Info<< " Writing volScalarField : "
269  << Pe.name() << endl;
270  Pe.write();
271  }
272  }
273  else
274  {
275  Info<< " No phi" << endl;
276  }
277 
278  Info<< "End" << nl << endl;
279 }
280 
281 
282 // ************************************************************************* //
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
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
singlePhaseTransportModel.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::calc
void calc(const argList &args, const Time &runTime, const fvMesh &mesh)
Definition: Lambda2.C:38
turbulentTransportModel.H
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
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 > &)
Foam::fvc::average
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
U
U
Definition: pEqn.H:46
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:60
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
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
surfaceInterpolate.H
Surface Interpolation.
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
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
fvcAverage.H
Area-weighted average a surfaceField creating a volField.
args
Foam::argList args(argc, argv)
turbulentFluidThermoModel.H