SprayCloudI.H
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-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class CloudType>
32 inline const Foam::SprayCloud<CloudType>&
34 {
35  return *cloudCopyPtr_;
36 }
37 
38 
39 template<class CloudType>
42 {
43  return atomizationModel_;
44 }
45 
46 
47 template<class CloudType>
50 {
51  return *atomizationModel_;
52 }
53 
54 
55 template<class CloudType>
58 {
59  return breakupModel_;
60 }
61 
62 
63 template<class CloudType>
66 {
67  return *breakupModel_;
68 }
69 
70 
71 template<class CloudType>
73 {
74  return averageParcelMass_;
75 }
76 
77 
78 template<class CloudType>
80 (
81  const scalar fraction
82 ) const
83 {
84  if ((fraction < 0) || (fraction > 1))
85  {
87  << "fraction should be in the range 0 < fraction < 1"
88  << exit(FatalError);
89  }
90 
91  scalar distance = 0.0;
92 
93  const label nParcel = this->size();
94  globalIndex globalParcels(nParcel);
95  const label nParcelSum = globalParcels.size();
96 
97  if (nParcelSum == 0)
98  {
99  return distance;
100  }
101 
102  // lists of parcels mass and distance from initial injection point
103  List<List<scalar>> procMass(Pstream::nProcs());
104  List<List<scalar>> procDist(Pstream::nProcs());
105 
106  List<scalar>& mass = procMass[Pstream::myProcNo()];
107  List<scalar>& dist = procDist[Pstream::myProcNo()];
108 
109  mass.setSize(nParcel);
110  dist.setSize(nParcel);
111 
112  label i = 0;
113  scalar mSum = 0.0;
114  for (const parcelType& p : *this)
115  {
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, Zero);
135  SortableList<scalar> allDist(nParcelSum, Zero);
136  for (const int proci : Pstream::allProcs())
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 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::SprayCloud
Templated base class for spray cloud.
Definition: SprayCloud.H:51
Foam::SortableList::sort
void sort()
Forward (stable) sort the list (if changed after construction).
Definition: SortableList.C:124
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::AtomizationModel
Templated atomization model class.
Definition: SprayCloud.H:49
Foam::globalIndex::localSize
label localSize() const
My local size.
Definition: globalIndexI.H:187
Foam::SprayCloud< Foam::DSMCCloud >::parcelType
Foam::DSMCCloud ::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition: SprayCloud.H:72
Foam::sumOp
Definition: ops.H:213
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::SprayCloud::cloudCopy
const SprayCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: SprayCloudI.H:33
Foam::FatalError
error FatalError
reduce
reduce(hasMovingMesh, orOp< bool >())
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:60
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::BreakupModel
Templated break-up model class.
Definition: SprayCloud.H:50
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::globalIndex::size
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:151
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::SprayCloud::averageParcelMass
scalar averageParcelMass() const
Return const-access to the average parcel mass.
Definition: SprayCloudI.H:72
Foam::SprayCloud::atomization
const AtomizationModel< SprayCloud< CloudType > > & atomization() const
Return const-access to the atomization model.
Definition: SprayCloudI.H:41
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:108
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:309
Foam::globalIndex::offset
label offset(const label proci) const
Start of proci data.
Definition: globalIndexI.H:163
Foam::SprayCloud::breakup
const BreakupModel< SprayCloud< CloudType > > & breakup() const
Return const-access to the breakup model.
Definition: SprayCloudI.H:57
Foam::SprayCloud::penetration
scalar penetration(const scalar fraction) const
Penetration for fraction [0-1] of the current total mass.
Definition: SprayCloudI.H:80