1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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
38 if ((0.0 >= time0) && (0.0 < time1))
40 return positions_.size();
49 template<class CloudType>
50 Foam::scalar Foam::ManualInjection<CloudType>::volumeToInject
56 // All parcels introduced at SOI
57 if ((0.0 >= time0) && (0.0 < time1))
59 return this->volumeTotal_;
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 template<class CloudType>
71 Foam::ManualInjection<CloudType>::ManualInjection
73 const dictionary& dict,
77 InjectionModel<CloudType>(dict, owner, typeName),
78 positionsFile_(this->coeffDict().lookup("positionsFile")),
84 owner.db().time().constant(),
90 diameters_(positions_.size()),
91 U0_(this->coeffDict().lookup("U0")),
96 this->coeffDict().subDict("parcelPDF"),
101 // Construct parcel diameters
102 forAll(diameters_, i)
104 diameters_[i] = parcelPDF_->sample();
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
129 template<class CloudType>
130 Foam::scalar Foam::ManualInjection<CloudType>::timeEnd() const
137 template<class CloudType>
138 Foam::vector Foam::ManualInjection<CloudType>::position
142 const polyMeshInfo& meshInfo
145 vector pos = positions_[iParcel];
146 if (meshInfo.caseIs2d())
148 if (meshInfo.caseIs2dWedge())
150 pos.component(meshInfo.emptyComponent()) = 0.0;
152 else if (meshInfo.caseIs2dSlab())
154 pos.component(meshInfo.emptyComponent()) =
155 meshInfo.centrePoint().component(meshInfo.emptyComponent());
161 "Foam::vector Foam::ManualInjection<CloudType>::position"
162 ) << "Could not determine 2-D case geometry" << nl
163 << abort(FatalError);
171 template<class CloudType>
172 Foam::vector Foam::ManualInjection<CloudType>::velocity
176 const polyMeshInfo& meshInfo
180 if (meshInfo.caseIs2dSlab())
182 vel.component(meshInfo.emptyComponent()) =
183 meshInfo.centrePoint().component(meshInfo.emptyComponent());
190 template<class CloudType>
191 Foam::scalar Foam::ManualInjection<CloudType>::d0
197 return diameters_[iParcel];
201 // ************************************************************************* //