KinematicParcel.C
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-2017 OpenFOAM Foundation
9  Copyright (C) 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 #include "KinematicParcel.H"
30 #include "forceSuSp.H"
31 #include "integrationScheme.H"
32 #include "meshTools.H"
33 #include "cloudSolution.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 template<class ParcelType>
39 
40 
41 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
42 
43 template<class ParcelType>
44 template<class TrackCloudType>
46 (
47  TrackCloudType& cloud,
48  trackingData& td
49 )
50 {
51  tetIndices tetIs = this->currentTetIndices();
52 
53  td.rhoc() = td.rhoInterp().interpolate(this->coordinates(), tetIs);
54 
55  if (td.rhoc() < cloud.constProps().rhoMin())
56  {
57  if (debug)
58  {
60  << "Limiting observed density in cell " << this->cell()
61  << " to " << cloud.constProps().rhoMin() << nl << endl;
62  }
63 
64  td.rhoc() = cloud.constProps().rhoMin();
65  }
66 
67  td.Uc() = td.UInterp().interpolate(this->coordinates(), tetIs);
68 
69  td.muc() = td.muInterp().interpolate(this->coordinates(), tetIs);
70 }
71 
72 
73 template<class ParcelType>
74 template<class TrackCloudType>
76 (
77  TrackCloudType& cloud,
78  trackingData& td,
79  const scalar dt
80 )
81 {
82  td.Uc() = cloud.dispersion().update
83  (
84  dt,
85  this->cell(),
86  U_,
87  td.Uc(),
88  UTurb_,
89  tTurb_
90  );
91 }
92 
93 
94 template<class ParcelType>
95 template<class TrackCloudType>
97 (
98  TrackCloudType& cloud,
99  trackingData& td,
100  const scalar dt
101 )
102 {
103  typename TrackCloudType::parcelType& p =
104  static_cast<typename TrackCloudType::parcelType&>(*this);
105 
106  this->UCorrect_ = Zero;
107 
108  this->UCorrect_ =
109  cloud.dampingModel().velocityCorrection(p, dt);
110 
111  this->UCorrect_ +=
112  cloud.packingModel().velocityCorrection(p, dt);
113 }
114 
115 
116 template<class ParcelType>
117 template<class TrackCloudType>
119 (
120  TrackCloudType& cloud,
121  trackingData& td,
122  const scalar dt
123 )
124 {
125  td.Uc() += cloud.UTrans()[this->cell()]/massCell(td);
126 }
127 
128 
129 template<class ParcelType>
130 template<class TrackCloudType>
132 (
133  TrackCloudType& cloud,
134  trackingData& td,
135  const scalar dt
136 )
137 {
138  // Define local properties at beginning of time step
139  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
140  const scalar np0 = nParticle_;
141  const scalar mass0 = mass();
142 
143  // Reynolds number
144  const scalar Re = this->Re(td);
145 
146 
147  // Sources
148  //~~~~~~~~
149 
150  // Explicit momentum source for particle
151  vector Su = Zero;
152 
153  // Linearised momentum source coefficient
154  scalar Spu = 0.0;
155 
156  // Momentum transfer from the particle to the carrier phase
157  vector dUTrans = Zero;
158 
159 
160  // Motion
161  // ~~~~~~
162 
163  // Calculate new particle velocity
164  this->U_ =
165  calcVelocity(cloud, td, dt, Re, td.muc(), mass0, Su, dUTrans, Spu);
166 
167  this->U_ += this->UCorrect_;
168 
169  // Accumulate carrier phase source terms
170  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
171  if (cloud.solution().coupled())
172  {
173  // Update momentum transfer
174  cloud.UTrans()[this->cell()] += np0*dUTrans;
175 
176  // Update momentum transfer coefficient
177  cloud.UCoeff()[this->cell()] += np0*Spu;
178  }
179 }
180 
181 
182 template<class ParcelType>
183 template<class TrackCloudType>
185 (
186  TrackCloudType& cloud,
187  trackingData& td,
188  const scalar dt,
189  const scalar Re,
190  const scalar mu,
191  const scalar mass,
192  const vector& Su,
193  vector& dUTrans,
194  scalar& Spu
195 ) const
196 {
197  const typename TrackCloudType::parcelType& p =
198  static_cast<const typename TrackCloudType::parcelType&>(*this);
199  typename TrackCloudType::parcelType::trackingData& ttd =
200  static_cast<typename TrackCloudType::parcelType::trackingData&>(td);
201 
202  const typename TrackCloudType::forceType& forces = cloud.forces();
203 
204  // Momentum source due to particle forces
205  const forceSuSp Fcp = forces.calcCoupled(p, ttd, dt, mass, Re, mu);
206  const forceSuSp Fncp = forces.calcNonCoupled(p, ttd, dt, mass, Re, mu);
207  const scalar massEff = forces.massEff(p, ttd, mass);
208 
209  /*
210  // Proper splitting ...
211  // Calculate the integration coefficients
212  const vector acp = (Fcp.Sp()*td.Uc() + Fcp.Su())/massEff;
213  const vector ancp = (Fncp.Sp()*td.Uc() + Fncp.Su() + Su)/massEff;
214  const scalar bcp = Fcp.Sp()/massEff;
215  const scalar bncp = Fncp.Sp()/massEff;
216 
217  // Integrate to find the new parcel velocity
218  const vector deltaUcp =
219  cloud.UIntegrator().partialDelta
220  (
221  U_, dt, acp + ancp, bcp + bncp, acp, bcp
222  );
223  const vector deltaUncp =
224  cloud.UIntegrator().partialDelta
225  (
226  U_, dt, acp + ancp, bcp + bncp, ancp, bncp
227  );
228  const vector deltaT = deltaUcp + deltaUncp;
229  */
230 
231  // Shortcut splitting assuming no implicit non-coupled force ...
232  // Calculate the integration coefficients
233  const vector acp = (Fcp.Sp()*td.Uc() + Fcp.Su())/massEff;
234  const vector ancp = (Fncp.Su() + Su)/massEff;
235  const scalar bcp = Fcp.Sp()/massEff;
236 
237  // Integrate to find the new parcel velocity
238  const vector deltaU = cloud.UIntegrator().delta(U_, dt, acp + ancp, bcp);
239  const vector deltaUncp = ancp*dt;
240  const vector deltaUcp = deltaU - deltaUncp;
241 
242  // Calculate the new velocity and the momentum transfer terms
243  vector Unew = U_ + deltaU;
244 
245  dUTrans -= massEff*deltaUcp;
246 
247  Spu = dt*Fcp.Sp();
248 
249  // Apply correction to velocity and dUTrans for reduced-D cases
250  const polyMesh& mesh = cloud.pMesh();
251  meshTools::constrainDirection(mesh, mesh.solutionD(), Unew);
252  meshTools::constrainDirection(mesh, mesh.solutionD(), dUTrans);
253 
254  return Unew;
255 }
256 
257 
258 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
259 
260 template<class ParcelType>
262 (
263  const KinematicParcel<ParcelType>& p
264 )
265 :
266  ParcelType(p),
267  active_(p.active_),
268  typeId_(p.typeId_),
269  nParticle_(p.nParticle_),
270  d_(p.d_),
271  dTarget_(p.dTarget_),
272  U_(p.U_),
273  rho_(p.rho_),
274  age_(p.age_),
275  tTurb_(p.tTurb_),
276  UTurb_(p.UTurb_),
277  UCorrect_(p.UCorrect_)
278 {}
279 
280 
281 template<class ParcelType>
283 (
284  const KinematicParcel<ParcelType>& p,
285  const polyMesh& mesh
286 )
287 :
288  ParcelType(p, mesh),
289  active_(p.active_),
290  typeId_(p.typeId_),
291  nParticle_(p.nParticle_),
292  d_(p.d_),
293  dTarget_(p.dTarget_),
294  U_(p.U_),
295  rho_(p.rho_),
296  age_(p.age_),
297  tTurb_(p.tTurb_),
298  UTurb_(p.UTurb_),
299  UCorrect_(p.UCorrect_)
300 {}
302 
303 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
304 
305 template<class ParcelType>
306 template<class TrackCloudType>
308 (
309  TrackCloudType& cloud,
310  trackingData& td,
311  const scalar trackTime
312 )
313 {
314  auto& p = static_cast<typename TrackCloudType::parcelType&>(*this);
315  auto& ttd =
316  static_cast<typename TrackCloudType::parcelType::trackingData&>(td);
317 
318  ttd.switchProcessor = false;
319  ttd.keepParticle = true;
320 
321  const cloudSolution& solution = cloud.solution();
322  const scalarField& cellLengthScale = cloud.cellLengthScale();
323  const scalar maxCo = solution.maxCo();
324 
325  while (ttd.keepParticle && !ttd.switchProcessor && p.stepFraction() < 1)
326  {
327  // Cache the current position, cell and step-fraction
328  const point start = p.position();
329  const scalar sfrac = p.stepFraction();
330 
331  // Total displacement over the time-step
332  const vector s = trackTime*U_;
333 
334  // Cell length scale
335  const scalar l = cellLengthScale[p.cell()];
336 
337  // Deviation from the mesh centre for reduced-D cases
338  const vector d = p.deviationFromMeshCentre();
339 
340  // Fraction of the displacement to track in this loop. This is limited
341  // to ensure that the both the time and distance tracked is less than
342  // maxCo times the total value.
343  scalar f = 1 - p.stepFraction();
344  f = min(f, maxCo);
345  f = min(f, maxCo*l/max(SMALL*l, mag(s)));
346  if (p.active())
347  {
348  // Track to the next face
349  p.trackToFace(f*s - d, f);
350  }
351  else
352  {
353  // At present the only thing that sets active_ to false is a stick
354  // wall interaction. We want the position of the particle to remain
355  // the same relative to the face that it is on. The local
356  // coordinates therefore do not change. We still advance in time and
357  // perform the relevant interactions with the fixed particle.
358  p.stepFraction() += f;
359  }
360 
361  const scalar dt = (p.stepFraction() - sfrac)*trackTime;
362 
363  // Avoid problems with extremely small timesteps
364  if (dt > ROOTVSMALL)
365  {
366  // Update cell based properties
367  p.setCellValues(cloud, ttd);
368 
369  p.calcDispersion(cloud, ttd, dt);
370 
371  if (solution.cellValueSourceCorrection())
372  {
373  p.cellValueSourceCorrection(cloud, ttd, dt);
374  }
375 
376  p.calcUCorrection(cloud, ttd, dt);
377 
378  p.calc(cloud, ttd, dt);
379  }
380 
381  p.age() += dt;
382 
383  if (p.active() && p.onFace())
384  {
385  ttd.keepParticle = cloud.functions().postFace(p, ttd);
386  }
387 
388  ttd.keepParticle = cloud.functions().postMove(p, dt, start, ttd);
389 
390  if (p.active() && p.onFace() && ttd.keepParticle)
391  {
392  p.hitFace(s, cloud, ttd);
393  }
394  }
395 
396  return ttd.keepParticle;
397 }
398 
399 
400 template<class ParcelType>
401 template<class TrackCloudType>
403 (
404  TrackCloudType& cloud,
405  trackingData& td
406 )
407 {
408  auto& p = static_cast<typename TrackCloudType::parcelType&>(*this);
409  auto& ttd =
410  static_cast<typename TrackCloudType::parcelType::trackingData&>(td);
411 
412  const polyPatch& pp = p.mesh().boundaryMesh()[p.patch()];
413 
414  // Invoke post-processing model
415  td.keepParticle = cloud.functions().postPatch(p, pp, ttd);
416 
417  if (isA<processorPolyPatch>(pp))
418  {
419  // Skip processor patches
420  return false;
421  }
422  else if (cloud.surfaceFilm().transferParcel(p, pp, td.keepParticle))
423  {
424  // Surface film model consumes the interaction, i.e. all
425  // interactions done
426  return true;
427  }
428  else
429  {
430  // This does not take into account the wall interaction model
431  // Just the polyPatch type. Then, a patch type which has 'rebound'
432  // interaction model will count as escaped parcel while it is not
433  if (!isA<wallPolyPatch>(pp) && !polyPatch::constraintType(pp.type()))
434  {
435  cloud.patchInteraction().addToEscapedParcels(nParticle_*mass());
436  }
437 
438  // Invoke patch interaction model
439  return cloud.patchInteraction().correct(p, pp, td.keepParticle);
440  }
441 }
442 
443 
444 template<class ParcelType>
445 template<class TrackCloudType>
447 (
448  TrackCloudType&,
449  trackingData& td
450 )
451 {
452  td.switchProcessor = true;
453 }
454 
455 
456 template<class ParcelType>
457 template<class TrackCloudType>
459 (
460  TrackCloudType&,
461  trackingData&
462 )
463 {
464  // wall interactions are handled by the generic hitPatch method
465 }
466 
467 
468 template<class ParcelType>
470 {
471  ParcelType::transformProperties(T);
473  U_ = transform(T, U_);
474 }
475 
476 
477 template<class ParcelType>
479 (
480  const vector& separation
481 )
482 {
483  ParcelType::transformProperties(separation);
484 }
485 
486 
487 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
488 
489 #include "KinematicParcelIO.C"
490 
491 // ************************************************************************* //
zeroField Su
Definition: alphaSuSp.H:1
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const interpolation< vector > & UInterp() const
Return const access to the interpolator for continuous phase velocity field.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
KinematicParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
dynamicFvMesh & mesh
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const vector calcVelocity(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar mu, const scalar mass, const vector &Su, vector &dUTrans, scalar &Spu) const
Calculate new particle velocity.
const Mesh & mesh() const noexcept
Return mesh.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
scalar maxCo
int debug
Static debugging option.
scalar rhoc() const
Return the continuous phase density.
scalar muc() const
Return the continuous phase viscosity.
labelList f(nPoints)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionedScalar mu
Atomic mass unit.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:49
PtrList< coordinateSystem > coordinates(solidRegions.size())
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:269
#define WarningInFunction
Report a warning using Foam::Warning.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
const interpolation< scalar > & rhoInterp() const
Return const access to the interpolator for continuous phase density field.
const interpolation< scalar > & muInterp() const
Return const access to the interpolator for continuous phase dynamic viscosity field.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:92
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:680
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Tensor of scalars, i.e. Tensor<scalar>.
void hitProcessorPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
bool hitPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
const vector & Uc() const
Return the continuous phase velocity.
void calcUCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct U following MP-PIC sub-models.
void calcDispersion(TrackCloudType &cloud, trackingData &td, const scalar dt)
Apply dispersion to the carrier phase velocity and update.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127