Formatting
[foam-extend-3.2.git] / src / lagrangian / dieselSpray / spray / sprayInject.C
blobc1163787bdece96bfff0f145fc155040276e94bb
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "spray.H"
27 #include "breakupModel.H"
28 #include "collisionModel.H"
29 #include "dispersionModel.H"
30 #include "injectorModel.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 void spray::inject()
41     scalar time = runTime_.value();
42     scalar time0 = time0_;
44     // Inject the parcels for each injector sequentially
45     forAll(injectors_, i)
46     {
47         autoPtr<injectorType>& it = injectors()[i].properties();
48         if (!it->pressureIndependentVelocity())
49         {
50             scalar referencePressure = p().average().value();
51             it->correctProfiles(fuels(), referencePressure);
52         }
54         const label nHoles = it->nHoles();
56         // parcels have the same mass during a timestep
57         scalar mass = it->mass(time0, time, twoD_, angleOfWedge_);
59         label Np = it->nParcelsToInject(time0, time);
61         if (mass > 0)
62         {
63             Np = max(1, Np);
64             scalar mp = mass/Np/nHoles;
66             // constT is only larger than zero for the first
67             // part of the injection
68             scalar constT = max(0.0, it->tsoi() - time0);
70             // deltaT is the duration of injection during this timestep
71             scalar deltaT = min
72             (
73                 runTime_.deltaT().value(),
74                 min
75                 (
76                     time - it->tsoi(),
77                     it->teoi() - time0
78                 )
79             );
81             for(label j=0; j<Np; j++)
82             {
83                 // calculate the time of injection for parcel 'j'
84                 scalar toi = time0 + constT + deltaT*j/scalar(Np);
86                 for(label n=0; n<nHoles; n++)
87                 {
89                     // calculate the velocity of the injected parcel
90                     vector injectionPosition = it->position
91                     (
92                         n,
93                         toi,
94                         twoD_,
95                         angleOfWedge_,
96                         axisOfSymmetry_,
97                         axisOfWedge_,
98                         axisOfWedgeNormal_,
99                         rndGen_
100                     );
102                     scalar diameter = injection().d0(i, toi);
103                     vector direction = injection().direction(i, n, toi, diameter);
104                     vector U = injection().velocity(i, toi)*direction;
106                     scalar symComponent = direction & axisOfSymmetry_;
107                     vector normal = direction - symComponent*axisOfSymmetry_;
108                     normal /= mag(normal);
110                     // should be set from dict or model
111                     scalar deviation = breakup().y0();
112                     scalar ddev = breakup().yDot0();
114                     label injectorCell = mesh_.findCell(injectionPosition);
116 #                   include "findInjectorCell.H"
118                     if (injectorCell >= 0)
119                     {
120                         scalar liquidCore = 1.0;
122                         // construct the parcel that is to be injected
124                         parcel* pPtr = new parcel
125                         (
126                             *this,
127                             injectionPosition,
128                             injectorCell,
129                             normal,
130                             diameter,
131                             it->T(toi),
132                             mp,
133                             deviation,
134                             ddev,
135                             0.0,
136                             0.0,
137                             0.0,
138                             liquidCore,
139                             scalar(i),
140                             U,
141                             vector::zero,
142                             it->X(),
143                             fuels_->components()
144                         );
146                         injectedLiquidKE_ += 0.5*pPtr->m()*magSqr(U);
148                         scalar dt = time - toi;
150                         pPtr->stepFraction() =
151                             (runTime_.deltaT().value() - dt)
152                            /runTime_.deltaT().value();
154                         bool keepParcel = pPtr->move(*this);
156                         if (keepParcel)
157                         {
158                             addParticle(pPtr);
159                         }
160                         else
161                         {
162                             delete pPtr;
163                         }
164                     } // if (injectorCell....
165                 } // for(label n=0...
166             } // for(label j=0....
167         } // if (mass>0)...
168     } // forAll(injectors)...
170     time0_ = time;
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 } // End namespace Foam
178 // ************************************************************************* //