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