initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / intermediate / submodels / Kinematic / InjectionModel / ConeInjectionMP / ConeInjectionMP.C
blob8d1b57f2317001594da9a904742f28524bb93c6e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-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 "ConeInjectionMP.H"
28 #include "DataEntry.H"
30 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
32 template<class CloudType>
33 Foam::label Foam::ConeInjectionMP<CloudType>::parcelsToInject
35     const scalar time0,
36     const scalar time1
37 ) const
39     if ((time0 >= 0.0) && (time0 < duration_))
40     {
41         const scalar targetVolume = volumeFlowRate_().integrate(0, time1);
43         const label targetParcels =
44             parcelsPerInjector_*targetVolume/this->volumeTotal_;
46         const label nToInject = targetParcels - nInjected_;
48         nInjected_ += nToInject;
50         return positions_.size()*nToInject;
51     }
52     else
53     {
54         return 0;
55     }
59 template<class CloudType>
60 Foam::scalar Foam::ConeInjectionMP<CloudType>::volumeToInject
62     const scalar time0,
63     const scalar time1
64 ) const
66     if ((time0 >= 0.0) && (time0 < duration_))
67     {
68         return volumeFlowRate_().integrate(time0, time1);
69     }
70     else
71     {
72         return 0.0;
73     }
77 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
79 template<class CloudType>
80 Foam::ConeInjectionMP<CloudType>::ConeInjectionMP
82     const dictionary& dict,
83     CloudType& owner
86     InjectionModel<CloudType>(dict, owner, typeName),
87     positionsFile_(this->coeffDict().lookup("positionsFile")),
88     positions_
89     (
90         IOobject
91         (
92             positionsFile_,
93             owner.db().time().constant(),
94             owner.mesh(),
95             IOobject::MUST_READ,
96             IOobject::NO_WRITE
97         )
98     ),
99     injectorCells_(positions_.size()),
100     axesFile_(this->coeffDict().lookup("axesFile")),
101     axes_
102     (
103         IOobject
104         (
105             axesFile_,
106             owner.db().time().constant(),
107             owner.mesh(),
108             IOobject::MUST_READ,
109             IOobject::NO_WRITE
110         )
111     ),
112     duration_(readScalar(this->coeffDict().lookup("duration"))),
113     parcelsPerInjector_
114     (
115         readScalar(this->coeffDict().lookup("parcelsPerInjector"))
116     ),
117     volumeFlowRate_
118     (
119         DataEntry<scalar>::New
120         (
121             "volumeFlowRate",
122             this->coeffDict()
123         )
124     ),
125     Umag_
126     (
127         DataEntry<scalar>::New
128         (
129             "Umag",
130             this->coeffDict()
131         )
132     ),
133     thetaInner_
134     (
135         DataEntry<scalar>::New
136         (
137             "thetaInner",
138             this->coeffDict()
139         )
140     ),
141     thetaOuter_
142     (
143         DataEntry<scalar>::New
144         (
145             "thetaOuter",
146             this->coeffDict()
147         )
148     ),
149     parcelPDF_
150     (
151         pdf::New
152         (
153             this->coeffDict().subDict("parcelPDF"),
154             owner.rndGen()
155         )
156     ),
157     nInjected_(this->parcelsAddedTotal()),
158     tanVec1_(positions_.size()),
159     tanVec2_(positions_.size())
161     // Normalise direction vector and determine direction vectors
162     // tangential to direction
163     forAll(axes_, i)
164     {
165         axes_[i] /= mag(axes_[i]);
167         vector tangent = vector::zero;
168         scalar magTangent = 0.0;
170         while (magTangent < SMALL)
171         {
172             vector v = this->owner().rndGen().vector01();
174             tangent = v - (v & axes_[i])*axes_[i];
175             magTangent = mag(tangent);
176         }
178         tanVec1_[i] = tangent/magTangent;
179         tanVec2_[i] = axes_[i]^tanVec1_[i];
180     }
182     // Set total volume to inject
183     this->volumeTotal_ = volumeFlowRate_().integrate(0.0, duration_);
185     // Set/cache the injector cells
186     forAll(positions_, i)
187     {
188         this->findCellAtPosition
189         (
190             injectorCells_[i],
191             positions_[i]
192         );
193     }
197 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
199 template<class CloudType>
200 Foam::ConeInjectionMP<CloudType>::~ConeInjectionMP()
204 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
206 template<class CloudType>
207 bool Foam::ConeInjectionMP<CloudType>::active() const
209     return true;
213 template<class CloudType>
214 Foam::scalar Foam::ConeInjectionMP<CloudType>::timeEnd() const
216     return this->SOI_ + duration_;
220 template<class CloudType>
221 void Foam::ConeInjectionMP<CloudType>::setPositionAndCell
223     const label parcelI,
224     const label,
225     const scalar,
226     vector& position,
227     label& cellOwner
230     const label i = parcelI%positions_.size();
232     position = positions_[i];
233     cellOwner = injectorCells_[i];
237 template<class CloudType>
238 void Foam::ConeInjectionMP<CloudType>::setProperties
240     const label parcelI,
241     const label,
242     const scalar time,
243     typename CloudType::parcelType& parcel
246     // set particle velocity
247     const label i = parcelI%positions_.size();
249     const scalar deg2Rad = mathematicalConstant::pi/180.0;
251     scalar t = time - this->SOI_;
252     scalar ti = thetaInner_().value(t);
253     scalar to = thetaOuter_().value(t);
254     scalar coneAngle = this->owner().rndGen().scalar01()*(to - ti) + ti;
256     coneAngle *= deg2Rad;
257     scalar alpha = sin(coneAngle);
258     scalar dcorr = cos(coneAngle);
259     scalar beta =
260         2.0*mathematicalConstant::pi*this->owner().rndGen().scalar01();
262     vector normal = alpha*(tanVec1_[i]*cos(beta) + tanVec2_[i]*sin(beta));
263     vector dirVec = dcorr*axes_[i];
264     dirVec += normal;
266     dirVec /= mag(dirVec);
268     parcel.U() = Umag_().value(t)*dirVec;
270     // set particle diameter
271     parcel.d() = parcelPDF_().sample();
275 template<class CloudType>
276 bool Foam::ConeInjectionMP<CloudType>::fullyDescribed() const
278     return false;
282 template<class CloudType>
283 bool Foam::ConeInjectionMP<CloudType>::validInjection(const label)
285     return true;
289 // ************************************************************************* //