engineSwirl.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  engineSwirl
28 
29 Group
30  grpPreProcessingUtilities
31 
32 Description
33  Generate a swirl flow for engine calculations.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
38 #include "unitConversion.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 int main(int argc, char *argv[])
43 {
44  argList::addNote
45  (
46  "Generate a swirl flow for engine calculations"
47  );
48 
49  #include "setRootCase.H"
50  #include "createTime.H"
51  #include "createNamedMesh.H"
52  #include "createFields.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56  scalar Vphi = (swirlRPMRatio * rpm * rpmToRads()).value();
57  scalar b1 = j1(swirlProfile).value();
58  scalar b2 = 2.0*b1/swirlProfile.value() - j0(swirlProfile).value();
59 
60  scalar omega = 0.125*(Vphi*bore*swirlProfile/b2).value();
61 
62  scalar cylinderRadius = 0.5*bore.value();
63 
64  scalar Umax = 0.0;
65  forAll(mesh.C(), celli)
66  {
67  vector c = mesh.C()[celli] - swirlCenter;
68  scalar r = ::pow(sqr(c & xT) + sqr(c & yT), 0.5);
69 
70  if (r <= cylinderRadius)
71  {
72  scalar b = j1(swirlProfile*r/cylinderRadius).value();
73  scalar vEff = omega*b;
74  r = max(r, SMALL);
75  U[celli] = ((vEff/r)*(c & yT))*xT + (-(vEff/r)*(c & xT))*yT;
76  Umax = max(Umax, mag(U[celli]));
77  }
78  }
79 
80  Info<< "Umax = " << Umax << endl;
81 
82  U.write();
83 
84  Info<< nl << "End" << nl << endl;
85 
86  return 0;
87 }
88 
89 
90 // ************************************************************************* //
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::rpmToRads
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
Definition: unitConversion.H:73
U
U
Definition: pEqn.H:72
setRootCase.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::j0
dimensionedScalar j0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:279
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
createTime.H
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
fvCFD.H
Foam::j1
dimensionedScalar j1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:280