PatchFlowRateInjection.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "distributionModel.H"
32#include "surfaceFields.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class CloudType>
38(
39 const dictionary& dict,
40 CloudType& owner,
41 const word& modelName
42)
43:
44 InjectionModel<CloudType>(dict, owner, modelName,typeName),
45 patchInjectionBase(owner.mesh(), this->coeffDict().getWord("patch")),
46 phiName_(this->coeffDict().template getOrDefault<word>("phi", "phi")),
47 rhoName_(this->coeffDict().template getOrDefault<word>("rho", "rho")),
48 duration_(this->coeffDict().getScalar("duration")),
49 concentration_
50 (
51 Function1<scalar>::New
52 (
53 "concentration",
54 this->coeffDict(),
55 &owner.mesh()
56 )
57 ),
58 parcelConcentration_
59 (
60 this->coeffDict().getScalar("parcelConcentration")
61 ),
62 sizeDistribution_
63 (
65 (
66 this->coeffDict().subDict("sizeDistribution"),
67 owner.rndGen()
68 )
69 )
70{
71 // Convert from user time to reduce the number of time conversion calls
72 const Time& time = owner.db().time();
73 duration_ = time.userTimeToTime(duration_);
74 concentration_->userTimeToTime(time);
75
77
78 // Re-initialise total mass/volume to inject to zero
79 // - will be reset during each injection
80 this->volumeTotal_ = 0.0;
81 this->massTotal_ = 0.0;
82}
83
84
85template<class CloudType>
87(
89)
90:
93 phiName_(im.phiName_),
94 rhoName_(im.rhoName_),
95 duration_(im.duration_),
96 concentration_(im.concentration_.clone()),
97 parcelConcentration_(im.parcelConcentration_),
98 sizeDistribution_(im.sizeDistribution_.clone())
99{}
100
101
102// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
103
104template<class CloudType>
106{}
107
108
109// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110
111template<class CloudType>
113{
114 patchInjectionBase::updateMesh(this->owner().mesh());
115}
116
117
118template<class CloudType>
120{
121 return this->SOI_ + duration_;
122}
123
124
125template<class CloudType>
127{
128 const polyMesh& mesh = this->owner().mesh();
129
130 const surfaceScalarField& phi =
132
133 const scalarField& phip = phi.boundaryField()[patchId_];
134
135 scalar flowRateIn = 0.0;
137 {
138 flowRateIn = max(0.0, -sum(phip));
139 }
140 else
141 {
142 const volScalarField& rho =
144 const scalarField& rhop = rho.boundaryField()[patchId_];
145
146 flowRateIn = max(0.0, -sum(phip/rhop));
147 }
148
149 reduce(flowRateIn, sumOp<scalar>());
150
151 return flowRateIn;
152}
153
154
155template<class CloudType>
157(
158 const scalar time0,
159 const scalar time1
160)
161{
162 if ((time0 >= 0.0) && (time0 < duration_))
163 {
164 scalar dt = time1 - time0;
165
166 scalar c = concentration_->value(0.5*(time0 + time1));
167
168 scalar nParcels = parcelConcentration_*c*flowRate()*dt;
169
170 Random& rnd = this->owner().rndGen();
171
172 label nParcelsToInject = floor(nParcels);
173
174 // Inject an additional parcel with a probability based on the
175 // remainder after the floor function
176 if
177 (
178 nParcelsToInject > 0
179 && (
180 nParcels - scalar(nParcelsToInject)
181 > rnd.globalPosition(scalar(0), scalar(1))
182 )
183 )
184 {
185 ++nParcelsToInject;
186 }
187
188 return nParcelsToInject;
189 }
190
191 return 0;
192}
193
194
195template<class CloudType>
197(
198 const scalar time0,
199 const scalar time1
200)
201{
202 scalar volume = 0.0;
203
204 if ((time0 >= 0.0) && (time0 < duration_))
205 {
206 scalar c = concentration_->value(0.5*(time0 + time1));
207
208 volume = c*(time1 - time0)*flowRate();
209 }
210
211 this->volumeTotal_ = volume;
212 this->massTotal_ = volume*this->owner().constProps().rho0();
213
214 return volume;
215}
216
217
218template<class CloudType>
220(
221 const label,
222 const label,
223 const scalar,
224 vector& position,
225 label& cellOwner,
226 label& tetFacei,
227 label& tetPti
228)
229{
231 (
232 this->owner().mesh(),
233 this->owner().rndGen(),
234 position,
235 cellOwner,
236 tetFacei,
237 tetPti
238 );
239}
240
241
242template<class CloudType>
244(
245 const label,
246 const label,
247 const scalar,
248 typename CloudType::parcelType& parcel
249)
250{
251 // Set particle velocity to carrier velocity
252 parcel.U() = this->owner().U()[parcel.cell()];
253
254 // Set particle diameter
255 parcel.d() = sizeDistribution_->sample();
256}
257
258
259template<class CloudType>
261{
262 return false;
263}
264
265
266template<class CloudType>
268{
269 return true;
270}
271
272
273// ************************************************************************* //
surfaceScalarField & phi
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
const dimensionSet & dimensions() const
Return dimensions.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
Templated injection model class.
scalar massTotal_
Total mass to inject [kg].
Patch injection, by using patch flow rate to determine concentration and velocity.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual ~PatchFlowRateInjection()
Destructor.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual scalar flowRate() const
Return the total volumetric flow rate across the patch [m3/s].
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Random number generator.
Definition: Random.H:60
Type globalPosition(const Type &start, const Type &end)
Return a sample on the interval [start,end].
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Type & value() const
Return const reference to value.
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
void updateMesh()
Update for new mesh topology.
const Time & time() const noexcept
Return time registry.
const Type & lookupObject(const word &name, const bool recursive=false) const
label setPositionAndCell(const fvMesh &mesh, const scalar fraction01, Random &rnd, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
Foam::surfaceFields.
Random rndGen
Definition: createFields.H:23