31 template<
class ParcelType>
32 template<
class TrackCloudType>
35 TrackCloudType& cloud,
37 const scalar trackTime
40 typename TrackCloudType::parcelType&
p =
41 static_cast<typename TrackCloudType::parcelType&
>(*this);
43 td.switchProcessor =
false;
44 td.keepParticle =
true;
54 while (td.keepParticle && !td.switchProcessor && p.stepFraction() < 1)
63 const vector d = p.deviationFromMeshCentre();
65 const scalar
f = 1 - p.stepFraction();
66 p.trackToAndHitFace(f*trackTime*Utracking - d, f, cloud, td);
69 return td.keepParticle;
73 template<
class ParcelType>
74 template<
class TrackCloudType>
81 template<
class ParcelType>
82 template<
class TrackCloudType>
89 td.switchProcessor =
true;
93 template<
class ParcelType>
94 template<
class TrackCloudType>
97 TrackCloudType& cloud,
101 const label wppIndex = this->patch();
106 this->
mesh().boundaryMesh()[wppIndex]
113 const scalar deltaT = cloud.pMesh().time().deltaTValue();
117 scalar m = constProps.
mass();
122 scalar U_dot_nw = U_ &
nw;
126 scalar invMagUnfA = 1/
max(
mag(U_dot_nw)*fA, vSmall);
128 cloud.rhoNBF()[wppIndex][wppLocalFace] += invMagUnfA;
130 cloud.rhoMBF()[wppIndex][wppLocalFace] += m*invMagUnfA;
132 cloud.linearKEBF()[wppIndex][wppLocalFace] +=
133 0.5*m*(U_ & U_)*invMagUnfA;
135 cloud.internalEBF()[wppIndex][wppLocalFace] += Ei_*invMagUnfA;
137 cloud.iDofBF()[wppIndex][wppLocalFace] +=
138 constProps.internalDegreesOfFreedom()*invMagUnfA;
140 cloud.momentumBF()[wppIndex][wppLocalFace] += m*Ut*invMagUnfA;
143 scalar preIE = 0.5*m*(U_ & U_) + Ei_;
148 cloud.wallInteraction().correct(*
this);
152 Ut = U_ - U_dot_nw*
nw;
154 invMagUnfA = 1/
max(
mag(U_dot_nw)*fA, vSmall);
156 cloud.rhoNBF()[wppIndex][wppLocalFace] += invMagUnfA;
158 cloud.rhoMBF()[wppIndex][wppLocalFace] += m*invMagUnfA;
160 cloud.linearKEBF()[wppIndex][wppLocalFace] +=
161 0.5*m*(U_ & U_)*invMagUnfA;
163 cloud.internalEBF()[wppIndex][wppLocalFace] += Ei_*invMagUnfA;
165 cloud.iDofBF()[wppIndex][wppLocalFace] +=
166 constProps.internalDegreesOfFreedom()*invMagUnfA;
168 cloud.momentumBF()[wppIndex][wppLocalFace] += m*Ut*invMagUnfA;
171 scalar postIE = 0.5*m*(U_ & U_) + Ei_;
176 scalar deltaQ = cloud.nParticle()*(preIE - postIE)/(deltaT*fA);
178 vector deltaFD = cloud.nParticle()*(preIMom - postIMom)/(deltaT*fA);
180 cloud.qBF()[wppIndex][wppLocalFace] += deltaQ;
182 cloud.fDBF()[wppIndex][wppLocalFace] += deltaFD;
186 template<
class ParcelType>
192 ParcelType::transformProperties(transform);
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A face is a list of labels corresponding to mesh vertices.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
Class to hold DSMC particle constant properties.
scalar mass() const
Return const access to the particle mass [kg].
const vectorField::subField faceAreas() const
Return face areas.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
const scalarField::subField magFaceAreas() const
Return face area magnitudes.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
ParcelType::trackingData trackingData
Use base tracking data.
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
label whichFace(const label l) const
Return label of face in patch from global face label.