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-2022 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>
72 inline Foam::scalar Foam::SprayCloud<CloudType>::averageParcelMass() const
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  const label nParcel = this->size();
92  const globalIndex globalParcels(nParcel);
93  const label nTotParcel = globalParcels.totalSize();
94 
95  if (nTotParcel == 0)
96  {
97  return 0;
98  }
99 
100  // lists of parcels mass and distance from initial injection point
101  List<scalar> mass(nParcel);
102  List<scalar> dist(nParcel);
103 
104  scalar mTotal = 0;
105  {
106  label i = 0;
107  for (const parcelType& p : *this)
108  {
109  scalar m = p.nParticle()*p.mass();
110  scalar d = mag(p.position() - p.position0());
111  mTotal += m;
112 
113  mass[i] = m;
114  dist[i] = d;
115  ++i;
116  }
117  }
118  // Total mass across all processors
119  reduce(mTotal, sumOp<scalar>());
120 
121  scalar distance = 0;
122  globalParcels.gatherInplace(mass);
123  globalParcels.gatherInplace(dist);
124 
125  if (Pstream::master())
126  {
127  if (nTotParcel == 1)
128  {
129  distance = dist[0];
130  }
131  else
132  {
133  // Distances - sored into ascending order
134  // Masses - leave unsorted
135 
136  const SortList<scalar> sortedDist(dist);
137 
138  const scalar mLimit = fraction*mTotal;
139  const labelList& indices = sortedDist.indices();
140 
141  if (mLimit > (mTotal - mass[indices.last()]))
142  {
143  distance = sortedDist.last();
144  }
145  else
146  {
147  // assuming that 'fraction' is generally closer to 1 than 0,
148  // loop through in reverse distance order
149  const scalar mThreshold = (1.0 - fraction)*mTotal;
150  scalar mCurrent = 0;
151  label i0 = 0;
152 
153  forAllReverse(indices, i)
154  {
155  label indI = indices[i];
156 
157  mCurrent += mass[indI];
158 
159  if (mCurrent > mThreshold)
160  {
161  i0 = i;
162  break;
163  }
164  }
165 
166  if (i0 == indices.size() - 1)
167  {
168  distance = sortedDist.last();
169  }
170  else
171  {
172  // linearly interpolate to determine distance
173  scalar alpha = (mCurrent - mThreshold)/mass[indices[i0]];
174  distance =
175  (
176  sortedDist[i0]
177  + alpha*(sortedDist[i0+1] - sortedDist[i0])
178  );
179  }
180  }
181  }
182  }
183 
184  Pstream::broadcast(distance);
185 
186  return distance;
187 }
188 
189 
190 // ************************************************************************* //
scalar averageParcelMass() const
Return const-access to the average parcel mass.
Definition: SprayCloudI.H:65
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const AtomizationModel< SprayCloud< CloudType > > & atomization() const
Return const-access to the atomization model.
Definition: SprayCloudI.H:34
Templated atomization model class.
Definition: SprayCloud.H:43
scalar penetration(const scalar fraction) const
Penetration for fraction [0-1] of the current total mass.
Definition: SprayCloudI.H:73
const SprayCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: SprayCloudI.H:26
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
const BreakupModel< SprayCloud< CloudType > > & breakup() const
Return const-access to the breakup model.
Definition: SprayCloudI.H:50
Templated base class for spray cloud.
Definition: SprayCloud.H:45
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
Templated break-up model class.
Definition: SprayCloud.H:44
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:437
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].