1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2009-2011 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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "PatchInjection.H"
27 #include "DataEntry.H"
28 #include "distributionModel.H"
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 template<class CloudType>
33 Foam::PatchInjection<CloudType>::PatchInjection
35 const dictionary& dict,
39 InjectionModel<CloudType>(dict, owner, typeName),
40 patchName_(this->coeffDict().lookup("patchName")),
41 patchId_(owner.mesh().boundaryMesh().findPatchID(patchName_)),
42 duration_(readScalar(this->coeffDict().lookup("duration"))),
45 readScalar(this->coeffDict().lookup("parcelsPerSecond"))
47 U0_(this->coeffDict().lookup("U0")),
50 DataEntry<scalar>::New("flowRateProfile", this->coeffDict())
54 distributionModels::distributionModel::New
56 this->coeffDict().subDict("sizeDistribution"),
67 "PatchInjection<CloudType>::PatchInjection"
72 ) << "Requested patch " << patchName_ << " not found" << nl
73 << "Available patches are: " << owner.mesh().boundaryMesh().names()
74 << nl << exit(FatalError);
77 const polyPatch& patch = owner.mesh().boundaryMesh()[patchId_];
79 cellOwners_ = patch.faceCells();
81 label patchSize = cellOwners_.size();
82 label totalPatchSize = patchSize;
83 reduce(totalPatchSize, sumOp<label>());
84 fraction_ = scalar(patchSize)/totalPatchSize;
86 // Set total volume/mass to inject
87 this->volumeTotal_ = fraction_*flowRateProfile_().integrate(0.0, duration_);
88 this->massTotal_ *= fraction_;
92 template<class CloudType>
93 Foam::PatchInjection<CloudType>::PatchInjection
95 const PatchInjection<CloudType>& im
98 InjectionModel<CloudType>(im),
99 patchName_(im.patchName_),
100 patchId_(im.patchId_),
101 duration_(im.duration_),
102 parcelsPerSecond_(im.parcelsPerSecond_),
104 flowRateProfile_(im.flowRateProfile_().clone().ptr()),
105 sizeDistribution_(im.sizeDistribution_().clone().ptr()),
106 cellOwners_(im.cellOwners_),
107 fraction_(im.fraction_)
111 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
113 template<class CloudType>
114 Foam::PatchInjection<CloudType>::~PatchInjection()
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 template<class CloudType>
121 Foam::scalar Foam::PatchInjection<CloudType>::timeEnd() const
123 return this->SOI_ + duration_;
127 template<class CloudType>
128 Foam::label Foam::PatchInjection<CloudType>::parcelsToInject
134 if ((time0 >= 0.0) && (time0 < duration_))
136 scalar nParcels =fraction_*(time1 - time0)*parcelsPerSecond_;
138 cachedRandom& rnd = this->owner().rndGen();
140 label nParcelsToInject = floor(nParcels);
142 // Inject an additional parcel with a probability based on the
143 // remainder after the floor function
148 nParcels - scalar(nParcelsToInject)
149 > rnd.position(scalar(0), scalar(1))
156 return nParcelsToInject;
165 template<class CloudType>
166 Foam::scalar Foam::PatchInjection<CloudType>::volumeToInject
172 if ((time0 >= 0.0) && (time0 < duration_))
174 return fraction_*flowRateProfile_().integrate(time0, time1);
183 template<class CloudType>
184 void Foam::PatchInjection<CloudType>::setPositionAndCell
195 if (cellOwners_.size() > 0)
197 cachedRandom& rnd = this->owner().rndGen();
199 label cellI = rnd.position<label>(0, cellOwners_.size() - 1);
201 cellOwner = cellOwners_[cellI];
203 // The position is between the face and cell centre, which could be
204 // in any tet of the decomposed cell, so arbitrarily choose the first
205 // face of the cell as the tetFace and the first point after the base
206 // point on the face as the tetPt. The tracking will pick the cell
207 // consistent with the motion in the firsttracking step.
208 tetFaceI = this->owner().mesh().cells()[cellOwner][0];
211 // position perturbed between cell and patch face centres
212 const vector& pc = this->owner().mesh().C()[cellOwner];
214 this->owner().mesh().Cf().boundaryField()[patchId_][cellI];
216 const scalar a = rnd.sample01<scalar>();
217 const vector d = pf - pc;
218 position = pc + 0.5*a*d;
226 position = pTraits<vector>::max;
231 template<class CloudType>
232 void Foam::PatchInjection<CloudType>::setProperties
237 typename CloudType::parcelType& parcel
240 // set particle velocity
243 // set particle diameter
244 parcel.d() = sizeDistribution_->sample();
248 template<class CloudType>
249 bool Foam::PatchInjection<CloudType>::fullyDescribed() const
255 template<class CloudType>
256 bool Foam::PatchInjection<CloudType>::validInjection(const label)
262 // ************************************************************************* //