fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / intermediate / submodels / Kinematic / InjectionModel / ConeInjectionMP / ConeInjectionMP.C
blob4f685c09b776617ebc742ef4c3b80a634a27ba13
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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"
29 #include "mathematicalConstants.H"
31 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
33 template<class CloudType>
34 Foam::label Foam::ConeInjectionMP<CloudType>::parcelsToInject
36     const scalar time0,
37     const scalar time1
38 ) const
40     if ((time0 >= 0.0) && (time0 < duration_))
41     {
42         const scalar targetVolume = volumeFlowRate_().integrate(0, time1);
44         const label targetParcels =
45             parcelsPerInjector_*targetVolume/this->volumeTotal_;
47         const label nToInject = targetParcels - nInjected_;
49         nInjected_ += nToInject;
51         return positions_.size()*nToInject;
52     }
53     else
54     {
55         return 0;
56     }
60 template<class CloudType>
61 Foam::scalar Foam::ConeInjectionMP<CloudType>::volumeToInject
63     const scalar time0,
64     const scalar time1
65 ) const
67     if ((time0 >= 0.0) && (time0 < duration_))
68     {
69         return volumeFlowRate_().integrate(time0, time1);
70     }
71     else
72     {
73         return 0.0;
74     }
78 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
80 template<class CloudType>
81 Foam::ConeInjectionMP<CloudType>::ConeInjectionMP
83     const dictionary& dict,
84     CloudType& owner
87     InjectionModel<CloudType>(dict, owner, typeName),
88     positionsFile_(this->coeffDict().lookup("positionsFile")),
89     positions_
90     (
91         IOobject
92         (
93             positionsFile_,
94             owner.db().time().constant(),
95             owner.mesh(),
96             IOobject::MUST_READ,
97             IOobject::NO_WRITE
98         )
99     ),
100     injectorCells_(positions_.size()),
101     axesFile_(this->coeffDict().lookup("axesFile")),
102     axes_
103     (
104         IOobject
105         (
106             axesFile_,
107             owner.db().time().constant(),
108             owner.mesh(),
109             IOobject::MUST_READ,
110             IOobject::NO_WRITE
111         )
112     ),
113     duration_(readScalar(this->coeffDict().lookup("duration"))),
114     parcelsPerInjector_
115     (
116         readScalar(this->coeffDict().lookup("parcelsPerInjector"))
117     ),
118     volumeFlowRate_
119     (
120         DataEntry<scalar>::New
121         (
122             "volumeFlowRate",
123             this->coeffDict()
124         )
125     ),
126     Umag_
127     (
128         DataEntry<scalar>::New
129         (
130             "Umag",
131             this->coeffDict()
132         )
133     ),
134     thetaInner_
135     (
136         DataEntry<scalar>::New
137         (
138             "thetaInner",
139             this->coeffDict()
140         )
141     ),
142     thetaOuter_
143     (
144         DataEntry<scalar>::New
145         (
146             "thetaOuter",
147             this->coeffDict()
148         )
149     ),
150     parcelPDF_
151     (
152         pdf::New
153         (
154             this->coeffDict().subDict("parcelPDF"),
155             owner.rndGen()
156         )
157     ),
158     nInjected_(this->parcelsAddedTotal()),
159     tanVec1_(positions_.size()),
160     tanVec2_(positions_.size())
162     // Normalise direction vector and determine direction vectors
163     // tangential to direction
164     forAll(axes_, i)
165     {
166         axes_[i] /= mag(axes_[i]);
168         vector tangent = vector::zero;
169         scalar magTangent = 0.0;
171         while (magTangent < SMALL)
172         {
173             vector v = this->owner().rndGen().vector01();
175             tangent = v - (v & axes_[i])*axes_[i];
176             magTangent = mag(tangent);
177         }
179         tanVec1_[i] = tangent/magTangent;
180         tanVec2_[i] = axes_[i]^tanVec1_[i];
181     }
183     // Set total volume to inject
184     this->volumeTotal_ = volumeFlowRate_().integrate(0.0, duration_);
186     // Set/cache the injector cells
187     forAll(positions_, i)
188     {
189         this->findCellAtPosition
190         (
191             injectorCells_[i],
192             positions_[i]
193         );
194     }
198 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
200 template<class CloudType>
201 Foam::ConeInjectionMP<CloudType>::~ConeInjectionMP()
205 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
207 template<class CloudType>
208 bool Foam::ConeInjectionMP<CloudType>::active() const
210     return true;
214 template<class CloudType>
215 Foam::scalar Foam::ConeInjectionMP<CloudType>::timeEnd() const
217     return this->SOI_ + duration_;
221 template<class CloudType>
222 void Foam::ConeInjectionMP<CloudType>::setPositionAndCell
224     const label parcelI,
225     const label,
226     const scalar,
227     vector& position,
228     label& cellOwner
231     const label i = parcelI%positions_.size();
233     position = positions_[i];
234     cellOwner = injectorCells_[i];
238 template<class CloudType>
239 void Foam::ConeInjectionMP<CloudType>::setProperties
241     const label parcelI,
242     const label,
243     const scalar time,
244     typename CloudType::parcelType& parcel
247     // set particle velocity
248     const label i = parcelI%positions_.size();
250     const scalar deg2Rad = mathematicalConstant::pi/180.0;
252     scalar t = time - this->SOI_;
253     scalar ti = thetaInner_().value(t);
254     scalar to = thetaOuter_().value(t);
255     scalar coneAngle = this->owner().rndGen().scalar01()*(to - ti) + ti;
257     coneAngle *= deg2Rad;
258     scalar alpha = sin(coneAngle);
259     scalar dcorr = cos(coneAngle);
260     scalar beta = mathematicalConstant::twoPi*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 // ************************************************************************* //