FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

sprayInject.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "spray.H"
00027 #include <dieselSpray/breakupModel.H>
00028 #include <dieselSpray/collisionModel.H>
00029 #include <dieselSpray/dispersionModel.H>
00030 #include <dieselSpray/injectorModel.H>
00031 
00032 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00033 
00034 namespace Foam
00035 {
00036 
00037 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00038 
00039 void spray::inject()
00040 {
00041     scalar time = runTime_.value();
00042     scalar time0 = time0_;
00043 
00044     // Inject the parcels for each injector sequentially
00045     forAll(injectors_, i)
00046     {
00047         autoPtr<injectorType>& it = injectors()[i].properties();
00048         if (!it->pressureIndependentVelocity())
00049         {
00050             scalar referencePressure = p().average().value();
00051             it->correctProfiles(fuels(), referencePressure);
00052         }
00053 
00054         const label nHoles = it->nHoles();
00055 
00056         // parcels have the same mass during a timestep
00057         scalar mass = it->mass(time0, time, twoD_, angleOfWedge_);
00058 
00059         label Np = it->nParcelsToInject(time0, time);
00060 
00061         if (mass > 0)
00062         {
00063             Np = max(1, Np);
00064             scalar mp = mass/Np/nHoles;
00065 
00066             // constT is only larger than zero for the first
00067             // part of the injection
00068             scalar constT = max(0.0, it->tsoi() - time0);
00069 
00070             // deltaT is the duration of injection during this timestep
00071             scalar deltaT = min
00072             (
00073                 runTime_.deltaT().value(),
00074                 min
00075                 (
00076                     time - it->tsoi(),
00077                     it->teoi() - time0
00078                 )
00079             );
00080 
00081             for(label j=0; j<Np; j++)
00082             {
00083                 // calculate the time of injection for parcel 'j'
00084                 scalar toi = time0 + constT + deltaT*j/scalar(Np);
00085 
00086                 for(label n=0; n<nHoles; n++)
00087                 {
00088 
00089                 // calculate the velocity of the injected parcel
00090                     vector injectionPosition = it->position
00091                 (
00092                         n,
00093                         toi,
00094                         twoD_,
00095                         angleOfWedge_,
00096                         axisOfSymmetry_,
00097                         axisOfWedge_,
00098                         axisOfWedgeNormal_,
00099                         rndGen_
00100                     );
00101 
00102                     scalar diameter = injection().d0(i, toi);
00103                     vector direction =
00104                         injection().direction(i, n, toi, diameter);
00105                     vector U = injection().velocity(i, toi)*direction;
00106 
00107                 scalar symComponent = direction & axisOfSymmetry_;
00108                 vector normal = direction - symComponent*axisOfSymmetry_;
00109                 normal /= mag(normal);
00110 
00111                 // should be set from dict or model
00112                 scalar deviation = breakup().y0();
00113                 scalar ddev = breakup().yDot0();
00114 
00115                     label injectorCell = mesh_.findCell(injectionPosition);
00116 
00117 #                   include "findInjectorCell.H"
00118 
00119                     if (injectorCell >= 0)
00120                     {
00121                         scalar liquidCore = 1.0;
00122 
00123                         // construct the parcel that is to be injected
00124 
00125                     parcel* pPtr = new parcel
00126                     (
00127                         *this,
00128                         injectionPosition,
00129                         injectorCell,
00130                         normal,
00131                         diameter,
00132                             it->T(toi),
00133                         mp,
00134                         deviation,
00135                         ddev,
00136                         0.0,
00137                         0.0,
00138                         0.0,
00139                         liquidCore,
00140                         scalar(i),
00141                         U,
00142                         vector::zero,
00143                             it->X(),
00144                         fuels_->components()
00145                     );
00146 
00147                         injectedLiquidKE_ += 0.5*pPtr->m()*magSqr(U);
00148 
00149                         scalar dt = time - toi;
00150 
00151                     pPtr->stepFraction() =
00152                         (runTime_.deltaT().value() - dt)
00153                        /runTime_.deltaT().value();
00154 
00155                         bool keepParcel = pPtr->move(*this);
00156 
00157                         if (keepParcel)
00158                         {
00159                             addParticle(pPtr);
00160                         }
00161                         else
00162                         {
00163                             delete pPtr;
00164                         }
00165                     } // if (injectorCell....
00166                 } // for(label n=0...
00167             } // for(label j=0....
00168         } // if (mass>0)...
00169     } // forAll(injectors)...
00170 
00171     time0_ = time;
00172 }
00173 
00174 
00175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00176 
00177 } // End namespace Foam
00178 
00179 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines