FreeStream.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-2024 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 "FreeStream.H"
27 #include "constants.H"
28 #include "triPointRef.H"
29 #include "tetIndices.H"
30 #include "standardNormal.H"
31 
32 using namespace Foam::constant::mathematical;
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class CloudType>
38 (
39  const dictionary& dict,
41 )
42 :
44  patches_(),
45  moleculeTypeIds_(),
46  numberDensities_(),
47  particleFluxAccumulators_()
48 {
49  // Identify which patches to use
50 
52 
53  forAll(cloud.mesh().boundaryMesh(), p)
54  {
55  const polyPatch& patch = cloud.mesh().boundaryMesh()[p];
56 
57  if (isType<polyPatch>(patch))
58  {
59  patches.append(p);
60  }
61  }
62 
63  patches_.transfer(patches);
64 
65  const dictionary& numberDensitiesDict
66  (
67  this->coeffDict().subDict("numberDensities")
68  );
69 
70  List<word> molecules(numberDensitiesDict.toc());
71 
72  // Initialise the particleFluxAccumulators_
73  particleFluxAccumulators_.setSize(patches_.size());
74 
75  forAll(patches_, p)
76  {
77  const polyPatch& patch = cloud.mesh().boundaryMesh()[patches_[p]];
78 
79  particleFluxAccumulators_[p] = List<Field<scalar>>
80  (
81  molecules.size(),
82  Field<scalar>(patch.size(), 0.0)
83  );
84  }
85 
86  moleculeTypeIds_.setSize(molecules.size());
87 
88  numberDensities_.setSize(molecules.size());
89 
90  forAll(molecules, i)
91  {
92  numberDensities_[i] =
93  numberDensitiesDict.lookup<scalar>(molecules[i]);
94 
95  moleculeTypeIds_[i] = findIndex(cloud.typeIdList(), molecules[i]);
96 
97  if (moleculeTypeIds_[i] == -1)
98  {
100  << "typeId " << molecules[i] << "not defined in cloud." << nl
101  << abort(FatalError);
102  }
103  }
104 
105  numberDensities_ /= cloud.nParticle();
106 }
107 
108 
109 
110 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
111 
112 template<class CloudType>
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
119 template<class CloudType>
121 {
122  CloudType& cloud(this->owner());
123  const polyMesh& mesh(cloud.mesh());
124 
125  forAll(patches_, p)
126  {
127  label patchi = patches_[p];
128 
129  const polyPatch& patch = mesh.boundaryMesh()[patchi];
130  List<Field<scalar>>& pFA = particleFluxAccumulators_[p];
131 
132  forAll(pFA, facei)
133  {
134  pFA[facei].setSize(patch.size(), 0);
135  }
136  }
137 }
138 
139 
140 template<class CloudType>
142 {
143  CloudType& cloud(this->owner());
144 
145  const polyMesh& mesh(cloud.mesh());
146 
147  const scalar deltaT = mesh.time().deltaTValue();
148 
149  randomGenerator& rndGen = cloud.rndGen();
150  distributions::standardNormal& stdNormal = cloud.stdNormal();
151 
152  scalar sqrtPi = sqrt(pi);
153 
154  label particlesInserted = 0;
155 
156  const volScalarField::Boundary& boundaryT
157  (
158  cloud.boundaryT().boundaryField()
159  );
160 
161  const volVectorField::Boundary& boundaryU
162  (
163  cloud.boundaryU().boundaryField()
164  );
165 
166  label nLocateBoundaryHits = 0;
167 
168  forAll(patches_, p)
169  {
170  label patchi = patches_[p];
171 
172  const polyPatch& patch = mesh.boundaryMesh()[patchi];
173 
174  // Add mass to the accumulators. negative face area dotted with the
175  // velocity to point flux into the domain.
176 
177  // Take a reference to the particleFluxAccumulator for this patch
178  List<Field<scalar>>& pFA = particleFluxAccumulators_[p];
179 
180  forAll(pFA, i)
181  {
182  label typeId = moleculeTypeIds_[i];
183 
184  scalar mass = cloud.constProps(typeId).mass();
185 
186  if (min(boundaryT[patchi]) < small)
187  {
189  << "Zero boundary temperature detected, check boundaryT "
190  << "condition." << nl
191  << nl << abort(FatalError);
192  }
193 
194  scalarField mostProbableSpeed
195  (
196  cloud.maxwellianMostProbableSpeed
197  (
198  boundaryT[patchi],
199  mass
200  )
201  );
202 
203  // Dotting boundary velocity with the face unit normal
204  // (which points out of the domain, so it must be
205  // negated), dividing by the most probable speed to form
206  // molecularSpeedRatio * cosTheta
207 
208  scalarField sCosTheta
209  (
210  (boundaryU[patchi] & -patch.faceAreas()/patch.magFaceAreas())
211  / mostProbableSpeed
212  );
213 
214  // From Bird eqn 4.22
215 
216  pFA[i] +=
217  patch.magFaceAreas()*numberDensities_[i]*deltaT
218  *mostProbableSpeed
219  *(
220  exp(-sqr(sCosTheta)) + sqrtPi*sCosTheta*(1 + erf(sCosTheta))
221  )
222  /(2.0*sqrtPi);
223  }
224 
225  forAll(patch, pFI)
226  {
227  // Loop over all faces as the outer loop to avoid calculating
228  // geometrical properties multiple times.
229 
230  const face& f = patch[pFI];
231 
232  label globalFaceIndex = pFI + patch.start();
233 
234  label celli = mesh.faceOwner()[globalFaceIndex];
235 
236  const vector& fC = patch.faceCentres()[pFI];
237 
238  scalar fA = patch.magFaceAreas()[pFI];
239 
241  (
242  mesh,
243  globalFaceIndex,
244  celli
245  );
246 
247  // Cumulative triangle area fractions
248  List<scalar> cTriAFracs(faceTets.size(), 0.0);
249 
250  scalar previousCumulativeSum = 0.0;
251 
252  forAll(faceTets, triI)
253  {
254  const tetIndices& faceTetIs = faceTets[triI];
255 
256  cTriAFracs[triI] =
257  faceTetIs.faceTri(mesh).mag()/fA
258  + previousCumulativeSum;
259 
260  previousCumulativeSum = cTriAFracs[triI];
261  }
262 
263  // Force the last area fraction value to 1.0 to avoid any
264  // rounding/non-flat face errors giving a value < 1.0
265  cTriAFracs.last() = 1.0;
266 
267  // Normal unit vector *negative* so normal is pointing into the
268  // domain
269  vector n = patch.faceAreas()[pFI];
270  n /= -mag(n);
271 
272  // Wall tangential unit vector. Use the direction between the
273  // face centre and the first vertex in the list
274  vector t1 = fC - (mesh.points()[f[0]]);
275  t1 /= mag(t1);
276 
277  // Other tangential unit vector. Rescaling in case face is not
278  // flat and n and t1 aren't perfectly orthogonal
279  vector t2 = n^t1;
280  t2 /= mag(t2);
281 
282  scalar faceTemperature = boundaryT[patchi][pFI];
283 
284  const vector& faceVelocity = boundaryU[patchi][pFI];
285 
286  forAll(pFA, i)
287  {
288  scalar& faceAccumulator = pFA[i][pFI];
289 
290  // Number of whole particles to insert
291  label nI = max(label(faceAccumulator), 0);
292 
293  // Add another particle with a probability proportional to the
294  // remainder of taking the integer part of faceAccumulator
295  if ((faceAccumulator - nI) > rndGen.scalar01())
296  {
297  nI++;
298  }
299 
300  faceAccumulator -= nI;
301 
302  label typeId = moleculeTypeIds_[i];
303 
304  scalar mass = cloud.constProps(typeId).mass();
305 
306  for (label i = 0; i < nI; i++)
307  {
308  // Choose a triangle to insert on, based on their relative
309  // area
310 
311  scalar triSelection = rndGen.scalar01();
312 
313  // Selected triangle
314  label selectedTriI = -1;
315 
316  forAll(cTriAFracs, triI)
317  {
318  selectedTriI = triI;
319 
320  if (cTriAFracs[triI] >= triSelection)
321  {
322  break;
323  }
324  }
325 
326  // Randomly distribute the points on the triangle.
327 
328  const tetIndices& faceTetIs = faceTets[selectedTriI];
329 
330  point p = faceTetIs.faceTri(mesh).randomPoint(rndGen);
331 
332  // Velocity generation
333 
334  scalar mostProbableSpeed
335  (
336  cloud.maxwellianMostProbableSpeed
337  (
338  faceTemperature,
339  mass
340  )
341  );
342 
343  scalar sCosTheta = (faceVelocity & n)/mostProbableSpeed;
344 
345  // Coefficients required for Bird eqn 12.5
346  scalar uNormProbCoeffA =
347  sCosTheta + sqrt(sqr(sCosTheta) + 2.0);
348 
349  scalar uNormProbCoeffB =
350  0.5*
351  (
352  1.0
353  + sCosTheta*(sCosTheta - sqrt(sqr(sCosTheta) + 2.0))
354  );
355 
356  // Equivalent to the QA value in Bird's DSMC3.FOR
357  scalar randomScaling = 3.0;
358 
359  if (sCosTheta < -3)
360  {
361  randomScaling = mag(sCosTheta) + 1;
362  }
363 
364  scalar P = -1;
365 
366  // Normalised candidates for the normal direction velocity
367  // component
368  scalar uNormal;
369  scalar uNormalThermal;
370 
371  // Select a velocity using Bird eqn 12.5
372  do
373  {
374  uNormalThermal =
375  randomScaling*(2.0*rndGen.scalar01() - 1);
376 
377  uNormal = uNormalThermal + sCosTheta;
378 
379  if (uNormal < 0.0)
380  {
381  P = -1;
382  }
383  else
384  {
385  P = 2.0*uNormal/uNormProbCoeffA
386  *exp(uNormProbCoeffB - sqr(uNormalThermal));
387  }
388 
389  } while (P < rndGen.scalar01());
390 
391  vector U =
392  sqrt(physicoChemical::k.value()*faceTemperature/mass)
393  *(stdNormal.sample()*t1 + stdNormal.sample()*t2)
394  + (t1 & faceVelocity)*t1
395  + (t2 & faceVelocity)*t2
396  + mostProbableSpeed*uNormal*n;
397 
398  scalar Ei = cloud.equipartitionInternalEnergy
399  (
400  faceTemperature,
401  cloud.constProps(typeId).internalDegreesOfFreedom()
402  );
403 
404  cloud.addNewParcel
405  (
406  p,
407  celli,
408  nLocateBoundaryHits,
409  U,
410  Ei,
411  typeId
412  );
413 
414  particlesInserted++;
415  }
416  }
417  }
418  }
419 
420  reduce(nLocateBoundaryHits, sumOp<label>());
421  if (nLocateBoundaryHits != 0)
422  {
424  << "Freestream inflow model for cloud " << this->owner().name()
425  << " did not accurately locate " << nLocateBoundaryHits
426  << " particles" << endl;
427  }
428 
429  reduce(particlesInserted, sumOp<label>());
430 
431  Info<< " Particles inserted = "
432  << particlesInserted << endl;
433 }
434 
435 
436 // ************************************************************************* //
label k
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
virtual void topoChange()
Update following mesh change.
Definition: FreeStream.C:120
FreeStream(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Definition: FreeStream.C:38
virtual void inflow()
Introduce particles.
Definition: FreeStream.C:141
virtual ~FreeStream()
Destructor.
Definition: FreeStream.C:113
Generic GeometricBoundaryField class.
Templated inflow boundary model class.
const dictionary & coeffDict() const
Return the coefficients dictionary.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
T & last()
Return the last element of the list.
Definition: UListI.H:128
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
wordList toc() const
Return the table of contents.
Definition: dictionary.C:976
Standard normal distribution. Not selectable.
virtual scalar sample() const
Sample the distribution.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
const Time & time() const
Return time.
static List< tetIndices > faceTetIndices(const polyMesh &mesh, label fI, label cI)
Return the tet decomposition of the given face, with.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const vectorField::subField faceAreas() const
Return face areas.
Definition: polyPatch.C:282
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:276
const scalarField::subField magFaceAreas() const
Return face area magnitudes.
Definition: polyPatch.C:288
Random number generator.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:82
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:123
scalar mag() const
Return scalar magnitude.
Definition: triangleI.H:103
Point randomPoint(randomGenerator &rndGen) const
Return a random point on the triangle from a uniform.
Definition: triangleI.H:217
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const fvPatchList & patches
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
mathematical constants.
dimensionedScalar exp(const dimensionedScalar &ds)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar erf(const dimensionedScalar &ds)
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
static const char nl
Definition: Ostream.H:266
labelList f(nPoints)
dictionary dict
randomGenerator rndGen(653213)
volScalarField & p