initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / dieselSpray / spray / sprayInject.C
blob394f361c68eec9f727b263b449c503391f95a81d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "spray.H"
28 #include "breakupModel.H"
29 #include "collisionModel.H"
30 #include "dispersionModel.H"
31 #include "injectorModel.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace Foam
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 void spray::inject()
42     scalar time = runTime_.value();
43     scalar time0 = time0_;
45     // Inject the parcels for each injector sequentially
46     forAll(injectors_, i)
47     {
48         autoPtr<injectorType>& it = injectors()[i].properties();
49         if (!it->pressureIndependentVelocity())
50         {
51             scalar referencePressure = p().average().value();
52             it->correctProfiles(fuels(), referencePressure);
53         }
55         const label nHoles = it->nHoles();
57         // parcels have the same mass during a timestep
58         scalar mass = it->mass(time0, time, twoD_, angleOfWedge_);
60         label Np = it->nParcelsToInject(time0, time);
62         if (mass > 0)
63         {
64             Np = max(1, Np);
65             scalar mp = mass/Np/nHoles;
67             // constT is only larger than zero for the first
68             // part of the injection
69             scalar constT = max(0.0, it->tsoi() - time0);
71             // deltaT is the duration of injection during this timestep
72             scalar deltaT = min
73             (
74                 runTime_.deltaT().value(),
75                 min
76                 (
77                     time - it->tsoi(),
78                     it->teoi() - time0
79                 )
80             );
82             for(label j=0; j<Np; j++)
83             {
84                 // calculate the time of injection for parcel 'j'
85                 scalar toi = time0 + constT + deltaT*j/scalar(Np);
87                 for(label n=0; n<nHoles; n++)
88                 {
90                     // calculate the velocity of the injected parcel
91                     vector injectionPosition = it->position
92                     (
93                         n,
94                         toi,
95                         twoD_,
96                         angleOfWedge_,
97                         axisOfSymmetry_,
98                         axisOfWedge_,
99                         axisOfWedgeNormal_,
100                         rndGen_
101                     );
103                     scalar diameter = injection().d0(i, toi);
104                     vector direction =
105                         injection().direction(i, n, toi, diameter);
106                     vector U = injection().velocity(i, toi)*direction;
108                     scalar symComponent = direction & axisOfSymmetry_;
109                     vector normal = direction - symComponent*axisOfSymmetry_;
110                     normal /= mag(normal);
112                     // should be set from dict or model
113                     scalar deviation = breakup().y0();
114                     scalar ddev = breakup().yDot0();
116                     label injectorCell = mesh_.findCell(injectionPosition);
118 #                   include "findInjectorCell.H"
120                     if (injectorCell >= 0)
121                     {
122                         scalar liquidCore = 1.0;
124                         // construct the parcel that is to be injected
126                         parcel* pPtr = new parcel
127                         (
128                             *this,
129                             injectionPosition,
130                             injectorCell,
131                             normal,
132                             diameter,
133                             it->T(toi),
134                             mp,
135                             deviation,
136                             ddev,
137                             0.0,
138                             0.0,
139                             0.0,
140                             liquidCore,
141                             scalar(i),
142                             U,
143                             vector::zero,
144                             it->X(),
145                             fuels_->components()
146                         );
148                         injectedLiquidKE_ += 0.5*pPtr->m()*magSqr(U);
150                         scalar dt = time - toi;
152                         pPtr->stepFraction() =
153                             (runTime_.deltaT().value() - dt)
154                            /runTime_.deltaT().value();
156                         bool keepParcel = pPtr->move(*this);
158                         if (keepParcel)
159                         {
160                             addParticle(pPtr);
161                         }
162                         else
163                         {
164                             delete pPtr;
165                         }
166                     } // if (injectorCell....
167                 } // for(label n=0...
168             } // for(label j=0....
169         } // if (mass>0)...
170     } // forAll(injectors)...
172     time0_ = time;
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 } // End namespace Foam
180 // ************************************************************************* //