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 | 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  engineSwirl
26 
27 Description
28  Generates a swirling flow for engine calulations.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "mathematicalConstants.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 int main(int argc, char *argv[])
38 {
39 
40  #include "setRootCase.H"
41  #include "createTime.H"
42  #include "createMesh.H"
43  #include "createFields.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47  scalar Vphi = (constant::mathematical::pi*swirlRPMRatio*rpm/30).value();
48  scalar b1 = j1(swirlProfile).value();
49  scalar b2 = 2.0*b1/swirlProfile.value() - j0(swirlProfile).value();
50 
51  scalar omega = 0.125*(Vphi*bore*swirlProfile/b2).value();
52 
53  scalar cylinderRadius = 0.5*bore.value();
54 
55  scalar Umax = 0.0;
56  forAll(mesh.C(), celli)
57  {
58  vector c = mesh.C()[celli] - swirlCenter;
59  scalar r = ::pow(sqr(c & xT) + sqr(c & yT), 0.5);
60 
61  if (r <= cylinderRadius)
62  {
63  scalar b = j1(swirlProfile*r/cylinderRadius).value();
64  scalar vEff = omega*b;
65  r = max(r, SMALL);
66  U[celli] = ((vEff/r)*(c & yT))*xT + (-(vEff/r)*(c & xT))*yT;
67  Umax = max(Umax, mag(U[celli]));
68  }
69  }
70 
71  Info<< "Umax = " << Umax << endl;
72 
73  U.write();
74 
75  Info<< "\n end\n";
76 
77  return 0;
78 }
79 
80 
81 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Type & value() const
Return const reference to value.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionedScalar j1(const dimensionedScalar &ds)
dynamicFvMesh & mesh
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar j0(const dimensionedScalar &ds)