ORourkeCollision.C
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-2018 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 "ORourkeCollision.H"
27 #include "SLGThermo.H"
28 #include "CompactListList.H"
29 #include "mathematicalConstants.H"
30 
31 using namespace Foam::constant::mathematical;
32 
33 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
34 
35 template<class CloudType>
37 (
38  typename CloudType::parcelType::trackingData& td,
39  const scalar dt
40 )
41 {
42  // Create the occupancy list for the cells
43  labelList occupancy(this->owner().mesh().nCells(), 0);
44  forAllIter(typename CloudType, this->owner(), iter)
45  {
46  occupancy[iter().cell()]++;
47  }
48 
49  // Initialize the sizes of the lists of parcels in each cell
50  CompactListList<parcelType*> pInCell(occupancy);
51 
52  // Reset the occupancy to use as a counter
53  occupancy = 0;
54 
55  // Set the parcel pointer lists for each cell
56  forAllIter(typename CloudType, this->owner(), iter)
57  {
58  pInCell(iter().cell(), occupancy[iter().cell()]++) = &iter();
59  }
60 
61  for (label celli=0; celli<this->owner().mesh().nCells(); celli++)
62  {
63  UList<parcelType*> pInCelli(pInCell[celli]);
64 
65  if (pInCelli.size() >= 2)
66  {
67  forAll(pInCelli, i)
68  {
69  for (label j=i+1; j<pInCelli.size(); j++)
70  {
71  parcelType& p1 = *pInCelli[i];
72  parcelType& p2 = *pInCelli[j];
73 
74  scalar m1 = p1.nParticle()*p1.mass();
75  scalar m2 = p2.nParticle()*p2.mass();
76 
77  bool massChanged = collideParcels(dt, p1, p2, m1, m2);
78 
79  if (massChanged)
80  {
81  if (m1 > rootVSmall)
82  {
83  const scalarField X(liquids_.X(p1.Y()));
84  p1.setCellValues(this->owner(), td);
85  p1.rho() = liquids_.rho(td.pc(), p1.T(), X);
86  p1.Cp() = liquids_.Cp(td.pc(), p1.T(), X);
87  p1.sigma() = liquids_.sigma(td.pc(), p1.T(), X);
88  p1.mu() = liquids_.mu(td.pc(), p1.T(), X);
89  p1.d() = cbrt(6.0*m1/(p1.nParticle()*p1.rho()*pi));
90  }
91 
92  if (m2 > rootVSmall)
93  {
94  const scalarField X(liquids_.X(p2.Y()));
95  p2.setCellValues(this->owner(), td);
96  p2.rho() = liquids_.rho(td.pc(), p2.T(), X);
97  p2.Cp() = liquids_.Cp(td.pc(), p2.T(), X);
98  p2.sigma() = liquids_.sigma(td.pc(), p2.T(), X);
99  p2.mu() = liquids_.mu(td.pc(), p2.T(), X);
100  p2.d() = cbrt(6.0*m2/(p2.nParticle()*p2.rho()*pi));
101  }
102  }
103  }
104  }
105  }
106  }
107 
108  // Remove coalesced parcels that fall below minimum mass threshold
109  forAllIter(typename CloudType, this->owner(), iter)
110  {
111  parcelType& p = iter();
112  scalar mass = p.nParticle()*p.mass();
113 
114  if (mass < this->owner().constProps().minParcelMass())
115  {
116  this->owner().deleteParticle(p);
117  }
118  }
119 }
120 
121 
122 template<class CloudType>
124 (
125  const scalar dt,
126  parcelType& p1,
127  parcelType& p2,
128  scalar& m1,
129  scalar& m2
130 )
131 {
132  // Return if parcel masses are ~0
133  if ((m1 < rootVSmall) || (m2 < rootVSmall))
134  {
135  return false;
136  }
137 
138  const scalar Vc = this->owner().mesh().V()[p1.cell()];
139  const scalar d1 = p1.d();
140  const scalar d2 = p2.d();
141 
142  scalar magUrel = mag(p1.U() - p2.U());
143  scalar sumD = d1 + d2;
144  scalar nu0 = 0.25*constant::mathematical::pi*sqr(sumD)*magUrel*dt/Vc;
145  scalar nMin = min(p1.nParticle(), p2.nParticle());
146  scalar nu = nMin*nu0;
147  scalar collProb = exp(-nu);
148  scalar xx = this->owner().rndGen().template sample01<scalar>();
149 
150  // Collision occurs
151  if (xx > collProb)
152  {
153  if (d1 > d2)
154  {
155  return collideSorted(dt, p1, p2, m1, m2);
156  }
157  else
158  {
159  return collideSorted(dt, p2, p1, m2, m1);
160  }
161  }
162  else
163  {
164  return false;
165  }
166 }
167 
168 
169 template<class CloudType>
171 (
172  const scalar dt,
173  parcelType& p1,
174  parcelType& p2,
175  scalar& m1,
176  scalar& m2
177 )
178 {
179  const scalar nP1 = p1.nParticle();
180  const scalar nP2 = p2.nParticle();
181 
182  const scalar sigma1 = p1.sigma();
183  const scalar sigma2 = p2.sigma();
184 
185  const scalar d1 = p1.d();
186  const scalar d2 = p2.d();
187 
188  const scalar T1 = p1.T();
189  const scalar T2 = p2.T();
190 
191  const scalar rho1 = p1.rho();
192  const scalar rho2 = p2.rho();
193 
194  const vector& U1 = p1.U();
195  const vector& U2 = p2.U();
196 
197 
198  vector URel = U1 - U2;
199  scalar magURel = mag(URel);
200 
201  scalar mTot = m1 + m2;
202 
203  scalar gamma = d1/max(rootVSmall, d2);
204  scalar f = pow3(gamma) + 2.7*gamma - 2.4*sqr(gamma);
205 
206  // Mass-averaged temperature
207  scalar Tave = (T1*m1 + T2*m2)/mTot;
208 
209  // Interpolate to find average surface tension
210  scalar sigmaAve = sigma1;
211  if (mag(T2 - T1) > small)
212  {
213  sigmaAve += (sigma2 - sigma1)*(Tave - T1)/(T2 - T1);
214  }
215 
216  scalar Vtot = m1/rho1 + m2/rho2;
217  scalar rhoAve = mTot/Vtot;
218 
219  scalar dAve = sqrt(d1*d2);
220  scalar WeColl = 0.5*rhoAve*sqr(magURel)*dAve/max(rootVSmall, sigmaAve);
221 
222  scalar coalesceProb = min(1.0, 2.4*f/max(rootVSmall, WeColl));
223 
224  scalar prob = this->owner().rndGen().template sample01<scalar>();
225 
226  // Coalescence
227  if (coalescence_ && prob < coalesceProb)
228  {
229  // Number of the droplets that coalesce
230  scalar nProb = prob*nP2/nP1;
231 
232  // Conservation of mass, momentum and energy
233  scalar m1Org = m1;
234  scalar m2Org = m2;
235 
236  scalar dm = nP1*nProb*m2/scalar(nP2);
237 
238  m1 += dm;
239  m2 -= dm;
240 
241  p1.T() = (Tave*mTot - m2*T2)/m1;
242 
243  p1.U() = (m1*U1 + (1.0 - m2/m2Org)*m2*U2)/m1;
244 
245  p1.Y() = (m1Org*p1.Y() + dm*p2.Y())/m1;
246 
247  p2.nParticle() = m2/(rho2*p2.volume());
248 
249  return true;
250  }
251  // Grazing collision (no coalescence)
252  else
253  {
254  scalar gf = sqrt(prob) - sqrt(coalesceProb);
255  scalar denom = 1.0 - sqrt(coalesceProb);
256  if (denom < 1.0e-5)
257  {
258  denom = 1.0;
259  }
260  gf /= denom;
261 
262  // If gf negative, this means that coalescence is turned off
263  // and these parcels should have coalesced
264  gf = max(0.0, gf);
265 
266  // gf -> 1 => v1p -> U1 ...
267  // gf -> 0 => v1p -> momentum/mTot
268 
269  vector mr = m1*U1 + m2*U2;
270  vector v1p = (mr + m2*gf*URel)/mTot;
271  vector v2p = (mr - m1*gf*URel)/mTot;
272 
273  if (nP1 < nP2)
274  {
275  p1.U() = v1p;
276  p2.U() = (nP1*v2p + (nP2 - nP1)*U2)/nP2;
277  }
278  else
279  {
280  p1.U() = (nP2*v1p + (nP1 - nP2)*U1)/nP1;
281  p2.U() = v2p;
282  }
283 
284  return false;
285  }
286 }
287 
288 
289 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
290 
291 template<class CloudType>
293 (
294  const dictionary& dict,
295  CloudType& owner,
296  const word& modelName
297 )
298 :
299  StochasticCollisionModel<CloudType>(dict, owner, modelName),
300  liquids_
301  (
302  owner.db().template lookupObject<SLGThermo>("SLGThermo").liquids()
303  ),
304  coalescence_(this->coeffDict().lookup("coalescence"))
305 {}
306 
307 
308 template<class CloudType>
310 (
312 )
313 :
315  liquids_(cm.liquids_),
316  coalescence_(cm.coalescence_)
317 {}
318 
319 
320 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
321 
322 template<class CloudType>
324 {}
325 
326 
327 // ************************************************************************* //
ORourkeCollision(const dictionary &dict, CloudType &cloud, const word &modelName=typeName)
Construct from dictionary.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const scalar mTot
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void collide(typename CloudType::parcelType::trackingData &td, const scalar dt)
Main collision routine.
virtual ~ORourkeCollision()
Destructor.
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
mathematical constants.
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar cbrt(const dimensionedScalar &ds)
virtual bool collideParcels(const scalar dt, parcelType &p1, parcelType &p2, scalar &m1, scalar &m2)
Collide parcels and return true if mass has changed.
volVectorField & U1
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
A packed storage unstructured matrix of objects of type <T> using an offset table for access...
CloudType::parcelType parcelType
Convenience typedef to the cloud&#39;s parcel type.
virtual bool collideSorted(const scalar dt, parcelType &p1, parcelType &p2, scalar &m1, scalar &m2)
labelList f(nPoints)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const liquidMixtureProperties & liquids_
Collision model by P.J. O&#39;Rourke.
dimensionedScalar pow3(const dimensionedScalar &ds)
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Templated stochastic collision model class.
volVectorField & U2
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar & rho2
Definition: createFields.H:40
volScalarField & p
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
const dimensionedScalar & rho1
Definition: createFields.H:39
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
volScalarField & nu
Switch coalescence_
Coalescence activation switch.