initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / intermediate / submodels / Kinematic / InjectionModel / ConeInjection / ConeInjection.C
blobb7dd300b45857b920c986b79e63726db8020f478
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 "ConeInjection.H"
28 #include "DataEntry.H"
30 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
32 template<class CloudType>
33 Foam::label Foam::ConeInjection<CloudType>::parcelsToInject
35     const scalar time0,
36     const scalar time1
37 ) const
39     if ((time0 >= 0.0) && (time0 < duration_))
40     {
41         return round((time1 - time0)*parcelsPerSecond_);
42     }
43     else
44     {
45         return 0;
46     }
50 template<class CloudType>
51 Foam::scalar Foam::ConeInjection<CloudType>::volumeToInject
53     const scalar time0,
54     const scalar time1
55 ) const
57     if ((time0 >= 0.0) && (time0 < duration_))
58     {
59         return volumeFlowRate_().integrate(time0, time1);
60     }
61     else
62     {
63         return 0.0;
64     }
68 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
70 template<class CloudType>
71 Foam::ConeInjection<CloudType>::ConeInjection
73     const dictionary& dict,
74     CloudType& owner
77     InjectionModel<CloudType>(dict, owner, typeName),
78     duration_(readScalar(this->coeffDict().lookup("duration"))),
79     position_(this->coeffDict().lookup("position")),
80     injectorCell_(-1),
81     direction_(this->coeffDict().lookup("direction")),
82     parcelsPerSecond_
83     (
84         readScalar(this->coeffDict().lookup("parcelsPerSecond"))
85     ),
86     volumeFlowRate_
87     (
88         DataEntry<scalar>::New
89         (
90             "volumeFlowRate",
91             this->coeffDict()
92         )
93     ),
94     Umag_
95     (
96         DataEntry<scalar>::New
97         (
98             "Umag",
99             this->coeffDict()
100         )
101     ),
102     thetaInner_
103     (
104         DataEntry<scalar>::New
105         (
106             "thetaInner",
107             this->coeffDict()
108         )
109     ),
110     thetaOuter_
111     (
112         DataEntry<scalar>::New
113         (
114             "thetaOuter",
115             this->coeffDict()
116         )
117     ),
118     parcelPDF_
119     (
120         pdf::New
121         (
122             this->coeffDict().subDict("parcelPDF"),
123             owner.rndGen()
124         )
125     ),
126     tanVec1_(vector::zero),
127     tanVec2_(vector::zero)
129     // Normalise direction vector
130     direction_ /= mag(direction_);
132     // Determine direction vectors tangential to direction
133     vector tangent = vector::zero;
134     scalar magTangent = 0.0;
136     while (magTangent < SMALL)
137     {
138         vector v = this->owner().rndGen().vector01();
140         tangent = v - (v & direction_)*direction_;
141         magTangent = mag(tangent);
142     }
144     tanVec1_ = tangent/magTangent;
145     tanVec2_ = direction_^tanVec1_;
147     // Set total volume to inject
148     this->volumeTotal_ = volumeFlowRate_().integrate(0.0, duration_);
150     // Set/cache the injector cell
151     this->findCellAtPosition(injectorCell_, position_);
155 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
157 template<class CloudType>
158 Foam::ConeInjection<CloudType>::~ConeInjection()
162 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
164 template<class CloudType>
165 bool Foam::ConeInjection<CloudType>::active() const
167     return true;
171 template<class CloudType>
172 Foam::scalar Foam::ConeInjection<CloudType>::timeEnd() const
174     return this->SOI_ + duration_;
178 template<class CloudType>
179 void Foam::ConeInjection<CloudType>::setPositionAndCell
181     const label,
182     const label,
183     const scalar,
184     vector& position,
185     label& cellOwner
188     position = position_;
189     cellOwner = injectorCell_;
193 template<class CloudType>
194 void Foam::ConeInjection<CloudType>::setProperties
196     const label parcelI,
197     const label,
198     const scalar time,
199     typename CloudType::parcelType& parcel
202     // set particle velocity
203     const scalar deg2Rad = mathematicalConstant::pi/180.0;
205     scalar t = time - this->SOI_;
206     scalar ti = thetaInner_().value(t);
207     scalar to = thetaOuter_().value(t);
208     scalar coneAngle = this->owner().rndGen().scalar01()*(to - ti) + ti;
210     coneAngle *= deg2Rad;
211     scalar alpha = sin(coneAngle);
212     scalar dcorr = cos(coneAngle);
213     scalar beta =
214         2.0*mathematicalConstant::pi*this->owner().rndGen().scalar01();
216     vector normal = alpha*(tanVec1_*cos(beta) + tanVec2_*sin(beta));
217     vector dirVec = dcorr*direction_;
218     dirVec += normal;
219     dirVec /= mag(dirVec);
221     parcel.U() = Umag_().value(t)*dirVec;
223     // set particle diameter
224     parcel.d() = parcelPDF_().sample();
228 template<class CloudType>
229 bool Foam::ConeInjection<CloudType>::fullyDescribed() const
231     return false;
235 template<class CloudType>
236 bool Foam::ConeInjection<CloudType>::validInjection(const label)
238     return true;
242 // ************************************************************************* //