33 template<
class ParcelType>
39 template<
class ParcelType>
40 template<
class TrackCloudType>
43 TrackCloudType&
cloud,
47 tetIndices tetIs = this->currentTetIndices(td.mesh);
49 td.
rhoc() = td.
rhoInterp().interpolate(this->coordinates(), tetIs);
51 if (td.
rhoc() <
cloud.constProps().rhoMin())
56 <<
"Limiting observed density in cell " << this->
cell()
57 <<
" to " <<
cloud.constProps().rhoMin() <<
nl <<
endl;
63 td.
Uc() = td.
UInterp().interpolate(this->coordinates(), tetIs);
65 td.
muc() = td.
muInterp().interpolate(this->coordinates(), tetIs);
69 template<
class ParcelType>
70 template<
class TrackCloudType>
73 TrackCloudType&
cloud,
78 td.
Uc() =
cloud.dispersion().update
90 template<
class ParcelType>
91 template<
class TrackCloudType>
94 TrackCloudType&
cloud,
99 td.
Uc() +=
cloud.UTransRef()[this->
cell()]/massCell(td);
103 template<
class ParcelType>
104 template<
class TrackCloudType>
107 TrackCloudType&
cloud,
114 const scalar np0 = nParticle_;
115 const scalar mass0 = mass();
118 const scalar
Re = this->
Re(td);
139 calcVelocity(
cloud, td, dt,
Re, td.
muc(), mass0,
Su, dUTrans, Spu);
144 if (
cloud.solution().coupled())
147 cloud.UTransRef()[this->
cell()] += np0*dUTrans;
150 cloud.UCoeffRef()[this->
cell()] += np0*Spu;
155 template<
class ParcelType>
156 template<
class TrackCloudType>
159 TrackCloudType&
cloud,
170 const typename TrackCloudType::parcelType&
p =
171 static_cast<const typename TrackCloudType::parcelType&
>(*this);
172 typename TrackCloudType::parcelType::trackingData& ttd =
173 static_cast<typename TrackCloudType::parcelType::trackingData&
>(td);
175 const typename TrackCloudType::forceType& forces =
cloud.forces();
178 const forceSuSp Fcp = forces.calcCoupled(
p, ttd, dt, mass,
Re,
mu);
179 const forceSuSp Fncp = forces.calcNonCoupled(
p, ttd, dt, mass,
Re,
mu);
180 const scalar massEff = forces.massEff(
p, ttd, mass);
206 const vector acp = (Fcp.
Sp()*td.
Uc() + Fcp.
Su())/massEff;
208 const scalar bcp = Fcp.
Sp()/massEff;
211 const vector deltaU =
cloud.UIntegrator().delta(U_, dt, acp + ancp, bcp);
212 const vector deltaUncp = ancp*dt;
213 const vector deltaUcp = deltaU - deltaUncp;
216 vector Unew = U_ + deltaU;
218 dUTrans -= massEff*deltaUcp;
233 template<
class ParcelType>
242 nParticle_(
p.nParticle_),
244 dTarget_(
p.dTarget_),
255 template<
class ParcelType>
256 template<
class TrackCloudType>
259 TrackCloudType&
cloud,
263 typename TrackCloudType::parcelType&
p =
264 static_cast<typename TrackCloudType::parcelType&
>(*this);
265 typename TrackCloudType::parcelType::trackingData& ttd =
266 static_cast<typename TrackCloudType::parcelType::trackingData&
>(td);
268 ttd.keepParticle =
true;
272 const scalar
maxCo =
cloud.solution().maxCo();
277 && ttd.sendToProc == -1
278 &&
p.stepFraction() < ttd.stepFractionRange().second()
281 if (
p.moving() &&
p.onFace())
283 cloud.functions().postFace(
p);
287 const point start =
p.position(td.mesh);
288 const scalar sfrac =
p.stepFraction();
291 const vector s = ttd.trackTime()*U_;
294 const scalar l = cellLengthScale[
p.cell()];
297 const vector d =
p.deviationFromMeshCentre(td.mesh);
302 scalar
f = ttd.stepFractionRange().second() -
p.stepFraction();
308 p.trackToFace(td.mesh,
f*
s - d,
f);
317 p.stepFraction() +=
f;
320 const scalar dt = (
p.stepFraction() - sfrac)*ttd.trackTime();
326 p.setCellValues(
cloud, ttd);
328 p.calcDispersion(
cloud, ttd, dt);
330 if (
cloud.solution().cellValueSourceCorrection())
332 p.cellValueSourceCorrection(
cloud, ttd, dt);
340 cloud.functions().postMove(
p, dt, start, ttd.keepParticle);
342 if (
p.moving() &&
p.onFace() && ttd.keepParticle)
344 cloud.functions().preFace(
p);
350 return ttd.keepParticle;
354 template<
class ParcelType>
360 ParcelType::transformProperties(
transform);
366 template<
class ParcelType>
367 template<
class TrackCloudType>
370 TrackCloudType&
cloud,
374 ParcelType::correctAfterParallelTransfer(
cloud, td);
376 typename TrackCloudType::parcelType&
p =
377 static_cast<typename TrackCloudType::parcelType&
>(*this);
379 const polyPatch& pp = td.mesh.boundaryMesh()[td.sendToPatch];
381 cloud.functions().postPatch(
p, pp);
385 template<
class ParcelType>
386 template<
class TrackCloudType>
389 TrackCloudType&
cloud,
393 typename TrackCloudType::parcelType&
p =
394 static_cast<typename TrackCloudType::parcelType&
>(*this);
396 const polyPatch& pp = td.mesh.boundaryMesh()[
p.patch(td.mesh)];
399 if (
cloud.surfaceFilm().transferParcel(
p, pp, td.keepParticle))
401 cloud.functions().postPatch(
p, pp);
406 if (
cloud.patchInteraction().correct(
p, pp, td.keepParticle))
408 cloud.functions().postPatch(
p, pp);
416 template<
class ParcelType>
417 template<
class TrackCloudType>
420 TrackCloudType&
cloud,
424 ParcelType::hitWedgePatch(
cloud, td);
426 typename TrackCloudType::parcelType&
p =
427 static_cast<typename TrackCloudType::parcelType&
>(*this);
429 const polyPatch& pp = td.mesh.boundaryMesh()[
p.patch(td.mesh)];
431 cloud.functions().postPatch(
p, pp);
435 template<
class ParcelType>
436 template<
class TrackCloudType>
439 TrackCloudType&
cloud,
443 ParcelType::hitSymmetryPlanePatch(
cloud, td);
445 typename TrackCloudType::parcelType&
p =
446 static_cast<typename TrackCloudType::parcelType&
>(*this);
448 const polyPatch& pp = td.mesh.boundaryMesh()[
p.patch(td.mesh)];
450 cloud.functions().postPatch(
p, pp);
454 template<
class ParcelType>
455 template<
class TrackCloudType>
458 TrackCloudType&
cloud,
462 ParcelType::hitSymmetryPatch(
cloud, td);
464 typename TrackCloudType::parcelType&
p =
465 static_cast<typename TrackCloudType::parcelType&
>(*this);
467 const polyPatch& pp = td.mesh.boundaryMesh()[
p.patch(td.mesh)];
469 cloud.functions().postPatch(
p, pp);
473 template<
class ParcelType>
474 template<
class TrackCloudType>
477 TrackCloudType&
cloud,
481 ParcelType::hitCyclicPatch(
cloud, td);
483 typename TrackCloudType::parcelType&
p =
484 static_cast<typename TrackCloudType::parcelType&
>(*this);
486 const polyPatch& pp = td.mesh.boundaryMesh()[
p.patch(td.mesh)];
488 cloud.functions().postPatch(
p, pp);
492 template<
class ParcelType>
493 template<
class TrackCloudType>
496 const vector& displacement,
497 const scalar fraction,
499 TrackCloudType&
cloud,
504 ParcelType::hitNonConformalCyclicPatch
513 if (td.sendToProc == -1 && result)
515 typename TrackCloudType::parcelType&
p =
516 static_cast<typename TrackCloudType::parcelType&
>(*this);
520 cloud.functions().postPatch(
p, pp);
527 template<
class ParcelType>
528 template<
class TrackCloudType>
535 const polyPatch& pp = td.mesh.boundaryMesh()[this->patch(td.mesh)];
538 <<
"Particle " << this->origId() <<
" hit " << pp.type() <<
" patch "
539 << pp.
name() <<
" at " << this->position(td.mesh)
540 <<
" but no interaction model was specified for this patch"
545 template<
class ParcelType>
546 template<
class TrackCloudType>
549 TrackCloudType&
cloud,
553 ParcelType::hitBasicPatch(
cloud, td);
555 typename TrackCloudType::parcelType&
p =
556 static_cast<typename TrackCloudType::parcelType&
>(*this);
558 const polyPatch& pp = td.mesh.boundaryMesh()[
p.patch(td.mesh)];
560 cloud.functions().postPatch(
p, pp);
const vector & Uc() const
Return the continuous phase velocity.
scalar muc() const
Return the continuous phase viscosity.
const interpolation< vector > & UInterp() const
Return const access to the interpolator for continuous.
scalar rhoc() const
Return the continuous phase density.
const interpolation< scalar > & muInterp() const
Return const access to the interpolator for continuous.
const interpolation< scalar > & rhoInterp() const
Return const access to the interpolator for continuous.
Momentum parcel class with rotational motion (as spherical particles only) and one/two-way coupling w...
void calcDispersion(TrackCloudType &cloud, trackingData &td, const scalar dt)
Apply dispersion to the carrier phase velocity and update.
bool move(TrackCloudType &cloud, trackingData &td)
Move the parcel.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
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.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
MomentumParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
bool hitPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
void correctAfterParallelTransfer(TrackCloudType &, trackingData &)
Make changes following a parallel transfer.
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
void hitBasicPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a basic.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
A cell is defined as a list of faces with extra functionality.
A cloud is a collection of lagrangian particles.
friend dimensionSet transform(const dimensionSet &)
Helper container for force Su and Sp terms.
const vector & Su() const
Return const access to the explicit contribution [kg m/s^2].
scalar Sp() const
Return const access to the implicit coefficient [kg/s].
const word & name() const
Return name.
Mesh consisting of general polyhedral cells.
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
A patch is a list of labels that address the faces in the global face list.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar mu
Atomic mass unit.
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionSet transform(const dimensionSet &)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
scalarField Re(const UList< complex > &cf)