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-------------------------------------------------------------------------------
11License
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
37template<class ParcelType>
39
40
41// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
42
43template<class ParcelType>
44template<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
73template<class ParcelType>
74template<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
94template<class ParcelType>
95template<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
116template<class ParcelType>
117template<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
129template<class ParcelType>
130template<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
182template<class ParcelType>
183template<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
260template<class ParcelType>
262(
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
281template<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{}
301
302
303// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
304
305template<class ParcelType>
306template<class TrackCloudType>
308(
309 TrackCloudType& cloud,
310 trackingData& td,
311 const scalar trackTime
312)
313{
314 typename TrackCloudType::parcelType& p =
315 static_cast<typename TrackCloudType::parcelType&>(*this);
316 typename TrackCloudType::parcelType::trackingData& ttd =
317 static_cast<typename TrackCloudType::parcelType::trackingData&>(td);
318
319 ttd.switchProcessor = false;
320 ttd.keepParticle = true;
321
322 const cloudSolution& solution = cloud.solution();
323 const scalarField& cellLengthScale = cloud.cellLengthScale();
324 const scalar maxCo = solution.maxCo();
325
326 while (ttd.keepParticle && !ttd.switchProcessor && p.stepFraction() < 1)
327 {
328 // Cache the current position, cell and step-fraction
329 const point start = p.position();
330 const scalar sfrac = p.stepFraction();
331
332 // Total displacement over the time-step
333 const vector s = trackTime*U_;
334
335 // Cell length scale
336 const scalar l = cellLengthScale[p.cell()];
337
338 // Deviation from the mesh centre for reduced-D cases
339 const vector d = p.deviationFromMeshCentre();
340
341 // Fraction of the displacement to track in this loop. This is limited
342 // to ensure that the both the time and distance tracked is less than
343 // maxCo times the total value.
344 scalar f = 1 - p.stepFraction();
345 f = min(f, maxCo);
346 f = min(f, maxCo*l/max(SMALL*l, mag(s)));
347 if (p.active())
348 {
349 // Track to the next face
350 p.trackToFace(f*s - d, f);
351 }
352 else
353 {
354 // At present the only thing that sets active_ to false is a stick
355 // wall interaction. We want the position of the particle to remain
356 // the same relative to the face that it is on. The local
357 // coordinates therefore do not change. We still advance in time and
358 // perform the relevant interactions with the fixed particle.
359 p.stepFraction() += f;
360 }
361
362 const scalar dt = (p.stepFraction() - sfrac)*trackTime;
363
364 // Avoid problems with extremely small timesteps
365 if (dt > ROOTVSMALL)
366 {
367 // Update cell based properties
368 p.setCellValues(cloud, ttd);
369
370 p.calcDispersion(cloud, ttd, dt);
371
372 if (solution.cellValueSourceCorrection())
373 {
374 p.cellValueSourceCorrection(cloud, ttd, dt);
375 }
376
377 p.calcUCorrection(cloud, ttd, dt);
378
379 p.calc(cloud, ttd, dt);
380 }
381
382 p.age() += dt;
383
384 if (p.active() && p.onFace())
385 {
386 cloud.functions().postFace(p, ttd.keepParticle);
387 }
388 cloud.functions().postMove(p, dt, start, ttd.keepParticle);
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
400template<class ParcelType>
401template<class TrackCloudType>
403(
404 TrackCloudType& cloud,
405 trackingData& td
406)
407{
408 typename TrackCloudType::parcelType& p =
409 static_cast<typename TrackCloudType::parcelType&>(*this);
410
411 const polyPatch& pp = p.mesh().boundaryMesh()[p.patch()];
412
413 // Invoke post-processing model
414 cloud.functions().postPatch(p, pp, td.keepParticle);
415
416 if (isA<processorPolyPatch>(pp))
417 {
418 // Skip processor patches
419 return false;
420 }
421 else if (cloud.surfaceFilm().transferParcel(p, pp, td.keepParticle))
422 {
423 // Surface film model consumes the interaction, i.e. all
424 // interactions done
425 return true;
426 }
427 else
428 {
429 // This does not take into account the wall interaction model
430 // Just the polyPatch type. Then, a patch type which has 'rebound'
431 // interaction model will count as escaped parcel while it is not
432 if (!isA<wallPolyPatch>(pp) && !polyPatch::constraintType(pp.type()))
433 {
434 cloud.patchInteraction().addToEscapedParcels(nParticle_*mass());
435 }
436
437 // Invoke patch interaction model
438 return cloud.patchInteraction().correct(p, pp, td.keepParticle);
439 }
440}
441
442
443template<class ParcelType>
444template<class TrackCloudType>
446(
447 TrackCloudType&,
448 trackingData& td
449)
450{
451 td.switchProcessor = true;
452}
453
454
455template<class ParcelType>
456template<class TrackCloudType>
458(
459 TrackCloudType&,
461)
462{
463 // wall interactions are handled by the generic hitPatch method
464}
465
466
467template<class ParcelType>
469{
471
472 U_ = transform(T, U_);
473}
474
475
476template<class ParcelType>
478(
479 const vector& separation
480)
481{
483}
484
485
486// * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
487
488#include "KinematicParcelIO.C"
489
490// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
const vector & Uc() const
Return the continuous phase velocity.
scalar muc() const
Return the continuous phase viscosity.
const interpolation< vector > & UInterp() const
scalar rhoc() const
Return the continuous phase density.
const interpolation< scalar > & muInterp() const
const interpolation< scalar > & rhoInterp() const
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
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.
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.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:54
const Switch cellValueSourceCorrection() const
Return const access to the cell value correction flag.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
virtual const word & constraintType() const
Return the constraint type this pointPatch implements.
virtual void move()=0
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:67
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Definition: forceSuSpI.H:61
scalar Sp() const
Return const access to the implicit coefficient [kg/s].
Definition: forceSuSpI.H:67
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.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:1072
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:84
volScalarField & p
const volScalarField & T
const volScalarField & mu
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
scalar maxCo
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))
zeroField Su
Definition: alphaSuSp.H:1
#define WarningInFunction
Report a warning using Foam::Warning.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)