initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / lagrangian / intermediate / submodels / Kinematic / InjectionModel / ManualInjection / ManualInjection.C
blob2af9bc20923394535b48061b2f2bc38dbf4be562
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 "ManualInjection.H"
29 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
31 template<class CloudType>
32 Foam::label Foam::ManualInjection<CloudType>::nParcelsToInject
34     const scalar time0,
35     const scalar time1
36 ) const
38     if ((0.0 >= time0) && (0.0 < time1))
39     {
40         return positions_.size();
41     }
42     else
43     {
44         return 0;
45     }
49 template<class CloudType>
50 Foam::scalar Foam::ManualInjection<CloudType>::volumeToInject
52     const scalar time0,
53     const scalar time1
54 ) const
56     // All parcels introduced at SOI
57     if ((0.0 >= time0) && (0.0 < time1))
58     {
59         return this->volumeTotal_;
60     }
61     else
62     {
63         return 0.0;
64     }
68 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
70 template<class CloudType>
71 Foam::ManualInjection<CloudType>::ManualInjection
73     const dictionary& dict,
74     CloudType& owner
77     InjectionModel<CloudType>(dict, owner, typeName),
78     positionsFile_(this->coeffDict().lookup("positionsFile")),
79     positions_
80     (
81         IOobject
82         (
83             positionsFile_,
84             owner.db().time().constant(),
85             owner.mesh(),
86             IOobject::MUST_READ,
87             IOobject::NO_WRITE
88         )
89     ),
90     diameters_(positions_.size()),
91     U0_(this->coeffDict().lookup("U0")),
92     parcelPDF_
93     (
94         pdf::New
95         (
96             this->coeffDict().subDict("parcelPDF"),
97             owner.rndGen()
98         )
99     )
101     // Construct parcel diameters
102     forAll(diameters_, i)
103     {
104         diameters_[i] = parcelPDF_->sample();
105     }
107     // Determine volume of particles to inject
108     this->volumeTotal_ = sum(pow(diameters_, 3))
109         *mathematicalConstant::pi/6.0;
113 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
115 template<class CloudType>
116 Foam::ManualInjection<CloudType>::~ManualInjection()
120 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
122 template<class CloudType>
123 bool Foam::ManualInjection<CloudType>::active() const
125     return true;
129 template<class CloudType>
130 Foam::scalar Foam::ManualInjection<CloudType>::timeEnd() const
132     // Not used
133     return 0.0;
137 template<class CloudType>
138 Foam::vector Foam::ManualInjection<CloudType>::position
140     const label iParcel,
141     const scalar time,
142     const polyMeshInfo& meshInfo
145     vector pos = positions_[iParcel];
146     if (meshInfo.caseIs2d())
147     {
148         if (meshInfo.caseIs2dWedge())
149         {
150             pos.component(meshInfo.emptyComponent()) = 0.0;
151         }
152         else if (meshInfo.caseIs2dSlab())
153         {
154             pos.component(meshInfo.emptyComponent()) =
155                 meshInfo.centrePoint().component(meshInfo.emptyComponent());
156         }
157         else
158         {
159             FatalErrorIn
160             (
161                 "Foam::vector Foam::ManualInjection<CloudType>::position"
162             )   << "Could not determine 2-D case geometry" << nl
163                 << abort(FatalError);
164         }
165     }
167     return pos;
171 template<class CloudType>
172 Foam::vector Foam::ManualInjection<CloudType>::velocity
174     const label,
175     const scalar,
176     const polyMeshInfo& meshInfo
179     vector vel = U0_;
180     if (meshInfo.caseIs2dSlab())
181     {
182         vel.component(meshInfo.emptyComponent()) =
183             meshInfo.centrePoint().component(meshInfo.emptyComponent());
184     }
186     return vel;
190 template<class CloudType>
191 Foam::scalar Foam::ManualInjection<CloudType>::d0
193     const label iParcel,
194     const scalar
195 ) const
197     return diameters_[iParcel];
201 // ************************************************************************* //