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-2022 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  const scalar fraction
74 ) const
75 {
76  if ((fraction < 0) || (fraction > 1))
77  {
79  << "fraction should be in the range 0 < fraction < 1"
80  << exit(FatalError);
81  }
82 
83  scalar distance = 0.0;
84 
85  const label nParcel = this->size();
86  globalIndex globalParcels(nParcel);
87  const label nParcelSum = globalParcels.size();
88 
89  if (nParcelSum == 0)
90  {
91  return distance;
92  }
93 
94  // lists of parcels mass and distance from initial injection point
95  List<List<scalar>> procMass(Pstream::nProcs());
96  List<List<scalar>> procDist(Pstream::nProcs());
97 
98  List<scalar>& mass = procMass[Pstream::myProcNo()];
99  List<scalar>& dist = procDist[Pstream::myProcNo()];
100 
101  mass.setSize(nParcel);
102  dist.setSize(nParcel);
103 
104  label i = 0;
105  scalar mSum = 0.0;
106  forAllConstIter(typename SprayCloud<CloudType>, *this, iter)
107  {
108  const parcelType& p = iter();
109  scalar m = p.nParticle()*p.mass();
110  scalar d = mag(p.position() - p.position0());
111  mSum += m;
112 
113  mass[i] = m;
114  dist[i] = d;
115 
116  i++;
117  }
118 
119  // calculate total mass across all processors
120  reduce(mSum, sumOp<scalar>());
121  Pstream::gatherList(procMass);
122  Pstream::gatherList(procDist);
123 
124  if (Pstream::master())
125  {
126  // flatten the mass lists
127  List<scalar> allMass(nParcelSum, 0.0);
128  SortableList<scalar> allDist(nParcelSum, 0.0);
129  for (label proci = 0; proci < Pstream::nProcs(); proci++)
130  {
132  (
133  allMass,
134  globalParcels.localSize(proci),
135  globalParcels.offset(proci)
136  ) = procMass[proci];
137 
138  // flatten the distance list
140  (
141  allDist,
142  globalParcels.localSize(proci),
143  globalParcels.offset(proci)
144  ) = procDist[proci];
145  }
146 
147  // sort allDist distances into ascending order
148  // note: allMass masses are left unsorted
149  allDist.sort();
150 
151  if (nParcelSum > 1)
152  {
153  const scalar mLimit = fraction*mSum;
154  const labelList& indices = allDist.indices();
155 
156  if (mLimit > (mSum - allMass[indices.last()]))
157  {
158  distance = allDist.last();
159  }
160  else
161  {
162  // assuming that 'fraction' is generally closer to 1 than 0,
163  // loop through in reverse distance order
164  const scalar mThreshold = (1.0 - fraction)*mSum;
165  scalar mCurrent = 0.0;
166  label i0 = 0;
167 
168  forAllReverse(indices, i)
169  {
170  label indI = indices[i];
171 
172  mCurrent += allMass[indI];
173 
174  if (mCurrent > mThreshold)
175  {
176  i0 = i;
177  break;
178  }
179  }
180 
181  if (i0 == indices.size() - 1)
182  {
183  distance = allDist.last();
184  }
185  else
186  {
187  // linearly interpolate to determine distance
188  scalar alpha = (mCurrent - mThreshold)/allMass[indices[i0]];
189  distance =
190  allDist[i0] + alpha*(allDist[i0+1] - allDist[i0]);
191  }
192  }
193  }
194  else
195  {
196  distance = allDist.first();
197  }
198  }
199 
200  Pstream::scatter(distance);
201 
202  return distance;
203 }
204 
205 
206 // ************************************************************************* //
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:112
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:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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:72
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
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