SprayCloudI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "SortableList.H"
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class CloudType>
31 inline const Foam::SprayCloud<CloudType>&
33 {
34  return cloudCopyPtr_();
35 }
36 
37 
38 template<class CloudType>
41 {
42  return atomisationModel_;
43 }
44 
45 
46 template<class CloudType>
49 {
50  return atomisationModel_();
51 }
52 
53 
54 template<class CloudType>
57 {
58  return breakupModel_;
59 }
60 
61 
62 template<class CloudType>
65 {
66  return breakupModel_();
67 }
68 
69 
70 template<class CloudType>
72 {
73  return averageParcelMass_;
74 }
75 
76 
77 template<class CloudType>
79 (
80  const scalar fraction
81 ) const
82 {
83  if ((fraction < 0) || (fraction > 1))
84  {
86  << "fraction should be in the range 0 < fraction < 1"
87  << exit(FatalError);
88  }
89 
90  scalar distance = 0.0;
91 
92  const label nParcel = this->size();
93  globalIndex globalParcels(nParcel);
94  const label nParcelSum = globalParcels.size();
95 
96  if (nParcelSum == 0)
97  {
98  return distance;
99  }
100 
101  // lists of parcels mass and distance from initial injection point
102  List<List<scalar>> procMass(Pstream::nProcs());
103  List<List<scalar>> procDist(Pstream::nProcs());
104 
105  List<scalar>& mass = procMass[Pstream::myProcNo()];
106  List<scalar>& dist = procDist[Pstream::myProcNo()];
107 
108  mass.setSize(nParcel);
109  dist.setSize(nParcel);
110 
111  label i = 0;
112  scalar mSum = 0.0;
113  forAllConstIter(typename SprayCloud<CloudType>, *this, iter)
114  {
115  const parcelType& p = iter();
116  scalar m = p.nParticle()*p.mass();
117  scalar d = mag(p.position() - p.position0());
118  mSum += m;
119 
120  mass[i] = m;
121  dist[i] = d;
122 
123  i++;
124  }
125 
126  // calculate total mass across all processors
127  reduce(mSum, sumOp<scalar>());
128  Pstream::gatherList(procMass);
129  Pstream::gatherList(procDist);
130 
131  if (Pstream::master())
132  {
133  // flatten the mass lists
134  List<scalar> allMass(nParcelSum, 0.0);
135  SortableList<scalar> allDist(nParcelSum, 0.0);
136  for (label proci = 0; proci < Pstream::nProcs(); proci++)
137  {
139  (
140  allMass,
141  globalParcels.localSize(proci),
142  globalParcels.offset(proci)
143  ) = procMass[proci];
144 
145  // flatten the distance list
147  (
148  allDist,
149  globalParcels.localSize(proci),
150  globalParcels.offset(proci)
151  ) = procDist[proci];
152  }
153 
154  // sort allDist distances into ascending order
155  // note: allMass masses are left unsorted
156  allDist.sort();
157 
158  if (nParcelSum > 1)
159  {
160  const scalar mLimit = fraction*mSum;
161  const labelList& indices = allDist.indices();
162 
163  if (mLimit > (mSum - allMass[indices.last()]))
164  {
165  distance = allDist.last();
166  }
167  else
168  {
169  // assuming that 'fraction' is generally closer to 1 than 0,
170  // loop through in reverse distance order
171  const scalar mThreshold = (1.0 - fraction)*mSum;
172  scalar mCurrent = 0.0;
173  label i0 = 0;
174 
175  forAllReverse(indices, i)
176  {
177  label indI = indices[i];
178 
179  mCurrent += allMass[indI];
180 
181  if (mCurrent > mThreshold)
182  {
183  i0 = i;
184  break;
185  }
186  }
187 
188  if (i0 == indices.size() - 1)
189  {
190  distance = allDist.last();
191  }
192  else
193  {
194  // linearly interpolate to determine distance
195  scalar alpha = (mCurrent - mThreshold)/allMass[indices[i0]];
196  distance =
197  allDist[i0] + alpha*(allDist[i0+1] - allDist[i0]);
198  }
199  }
200  }
201  else
202  {
203  distance = allDist.first();
204  }
205  }
206 
207  Pstream::scatter(distance);
208 
209  return distance;
210 }
211 
212 
213 // ************************************************************************* //
scalar averageParcelMass() const
Return const-access to the average parcel mass.
Definition: SprayCloudI.H:71
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:112
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:80
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const AtomisationModel< SprayCloud< CloudType > > & atomisation() const
Return const-access to the atomisation model.
Definition: SprayCloudI.H:40
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
scalar penetration(const scalar fraction) const
Penetration for fraction [0-1] of the current total mass.
Definition: SprayCloudI.H:79
const SprayCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: SprayCloudI.H:32
Templated atomisation model class.
Definition: SprayCloud.H:52
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
T & first()
Return the first element of the list.
Definition: UListI.H:114
label offset(const label proci) const
Start of proci data.
Definition: globalIndexI.H:48
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A List obtained as a section of another List.
Definition: SubList.H:53
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const BreakupModel< SprayCloud< CloudType > > & breakup() const
Return const-access to the breakup model.
Definition: SprayCloudI.H:56
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:96
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
Templated base class for spray cloud.
Definition: SprayCloud.H:69
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void setSize(const label)
Reset size of List.
Definition: List.C:281
dimensioned< scalar > mag(const dimensioned< Type > &)
Templated break-up model class.
Definition: SprayCloud.H:55
volScalarField & p
T & last()
Return the last element of the list.
Definition: UListI.H:128
label localSize() const
My local size.
Definition: globalIndexI.H:60