parallel running fixes
[OpenFOAM-1.6.x.git] / src / lagrangian / intermediate / submodels / Kinematic / InjectionModel / PatchInjection / PatchInjection.C
blob4e88b90bfebf896aab6bdad356bda3f6dfd301e3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2009-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 "PatchInjection.H"
28 #include "DataEntry.H"
29 #include "pdf.H"
31 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
33 template<class CloudType>
34 Foam::label Foam::PatchInjection<CloudType>::parcelsToInject
36     const scalar time0,
37     const scalar time1
38 ) const
40     if ((time0 >= 0.0) && (time0 < duration_))
41     {
42         return round(fraction_*(time1 - time0)*parcelsPerSecond_);
43     }
44     else
45     {
46         return 0;
47     }
51 template<class CloudType>
52 Foam::scalar Foam::PatchInjection<CloudType>::volumeToInject
54     const scalar time0,
55     const scalar time1
56 ) const
58     if ((time0 >= 0.0) && (time0 < duration_))
59     {
60         return fraction_*volumeFlowRate_().integrate(time0, time1);
61     }
62     else
63     {
64         return 0.0;
65     }
69 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
71 template<class CloudType>
72 Foam::PatchInjection<CloudType>::PatchInjection
74     const dictionary& dict,
75     CloudType& owner
78     InjectionModel<CloudType>(dict, owner, typeName),
79     patchName_(this->coeffDict().lookup("patchName")),
80     duration_(readScalar(this->coeffDict().lookup("duration"))),
81     parcelsPerSecond_
82     (
83         readScalar(this->coeffDict().lookup("parcelsPerSecond"))
84     ),
85     U0_(this->coeffDict().lookup("U0")),
86     volumeFlowRate_
87     (
88         DataEntry<scalar>::New
89         (
90             "volumeFlowRate",
91             this->coeffDict()
92         )
93     ),
94     parcelPDF_
95     (
96         pdf::New
97         (
98             this->coeffDict().subDict("parcelPDF"),
99             owner.rndGen()
100         )
101     ),
102     cellOwners_(),
103     fraction_(1.0)
105     label patchId = owner.mesh().boundaryMesh().findPatchID(patchName_);
107     if (patchId < 0)
108     {
109         FatalErrorIn
110         (
111             "PatchInjection<CloudType>::PatchInjection"
112             "("
113                 "const dictionary&, "
114                 "CloudType&"
115             ")"
116         )   << "Requested patch " << patchName_ << " not found" << nl
117             << "Available patches are: " << owner.mesh().boundaryMesh().names()
118             << nl << exit(FatalError);
119     }
121     const polyPatch& patch = owner.mesh().boundaryMesh()[patchId];
123     cellOwners_ = patch.faceCells();
125     label patchSize = cellOwners_.size();
126     label totalPatchSize = patchSize;
127     reduce(totalPatchSize, sumOp<label>());
128     fraction_ = scalar(patchSize)/totalPatchSize;
130     // Set total volume/mass to inject
131     this->volumeTotal_ = fraction_*volumeFlowRate_().integrate(0.0, duration_);
132     this->massTotal_ *= fraction_;
136 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
138 template<class CloudType>
139 Foam::PatchInjection<CloudType>::~PatchInjection()
143 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
145 template<class CloudType>
146 bool Foam::PatchInjection<CloudType>::active() const
148     return true;
152 template<class CloudType>
153 Foam::scalar Foam::PatchInjection<CloudType>::timeEnd() const
155     return this->SOI_ + duration_;
159 template<class CloudType>
160 void Foam::PatchInjection<CloudType>::setPositionAndCell
162     const label,
163     const label,
164     const scalar,
165     vector& position,
166     label& cellOwner
169     if (cellOwners_.size() > 0)
170     {
171         label cellI = this->owner().rndGen().integer(0, cellOwners_.size() - 1);
173         cellOwner = cellOwners_[cellI];
174         position = this->owner().mesh().C()[cellOwner];
175     }
176     else
177     {
178         cellOwner = -1;
179         // dummy position
180         position = pTraits<vector>::max;
181     }
185 template<class CloudType>
186 void Foam::PatchInjection<CloudType>::setProperties
188     const label,
189     const label,
190     const scalar,
191     typename CloudType::parcelType& parcel
194     // set particle velocity
195     parcel.U() = U0_;
197     // set particle diameter
198     parcel.d() = parcelPDF_->sample();
202 template<class CloudType>
203 bool Foam::PatchInjection<CloudType>::fullyDescribed() const
205     return false;
209 template<class CloudType>
210 bool Foam::PatchInjection<CloudType>::validInjection(const label)
212     return true;
216 // ************************************************************************* //