Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include <OpenFOAM/polyMesh.H>
00027
00028
00029
00030 template<class ParticleType>
00031 inline Foam::scalar Foam::Particle<ParticleType>::lambda
00032 (
00033 const vector& from,
00034 const vector& to,
00035 const label facei,
00036 const scalar stepFraction
00037 ) const
00038 {
00039 const polyMesh& mesh = cloud_.polyMesh_;
00040 bool movingMesh = mesh.moving();
00041
00042 if (movingMesh)
00043 {
00044 vector Sf = mesh.faceAreas()[facei];
00045 Sf /= mag(Sf);
00046 vector Cf = mesh.faceCentres()[facei];
00047
00048
00049 if (!cloud_.internalFace(facei))
00050 {
00051 label patchi = patch(facei);
00052 const polyPatch& patch = mesh.boundaryMesh()[patchi];
00053
00054
00055 if (isA<wallPolyPatch>(patch))
00056 {
00057 const vector& C = mesh.cellCentres()[celli_];
00058 scalar CCf = mag((C - Cf) & Sf);
00059
00060
00061 const ParticleType& p =
00062 static_cast<const ParticleType&>(*this);
00063 if (CCf > p.wallImpactDistance(Sf))
00064 {
00065 Cf -=p.wallImpactDistance(Sf)*Sf;
00066 }
00067 }
00068 }
00069
00070
00071
00072
00073 const vectorField& oldPoints = mesh.oldPoints();
00074
00075 vector Cf00 = mesh.faces()[facei].centre(oldPoints);
00076 vector Cf0 = Cf00 + stepFraction*(Cf - Cf00);
00077
00078 vector Sf00 = mesh.faces()[facei].normal(oldPoints);
00079
00080
00081 if (mag(Sf00) > SMALL)
00082 {
00083 Sf00 /= mag(Sf00);
00084 }
00085 else
00086 {
00087 Sf00 = Sf;
00088 }
00089
00090 scalar magSfDiff = mag(Sf - Sf00);
00091
00092
00093 if (magSfDiff > SMALL)
00094 {
00095 vector Sf0 = Sf00 + stepFraction*(Sf - Sf00);
00096
00097
00098 vector omega = Sf0 ^ Sf;
00099 scalar magOmega = mag(omega);
00100 omega /= magOmega + SMALL;
00101 vector n0 = omega ^ Sf0;
00102 scalar lam = ((Cf - Cf0) & Sf)/(n0 & Sf);
00103 vector r0 = Cf0 + lam*n0;
00104
00105
00106
00107
00108
00109
00110 vector alpha = from - r0;
00111 vector beta = to - from;
00112 scalar a = beta & (Sf - Sf0);
00113 scalar b = (alpha & (Sf - Sf0)) + (beta & Sf0);
00114 scalar c = alpha & Sf0;
00115
00116 if (mag(a) > SMALL)
00117 {
00118
00119 scalar ap = b/a;
00120 scalar bp = c/a;
00121 scalar cp = ap*ap - 4.0*bp;
00122
00123 if (cp < 0.0)
00124 {
00125
00126 return GREAT;
00127 }
00128 else
00129 {
00130 scalar l1 = -0.5*(ap - ::sqrt(cp));
00131 scalar l2 = -0.5*(ap + ::sqrt(cp));
00132
00133
00134
00135 if (mag(l1) < mag(l2))
00136 {
00137 return l1;
00138 }
00139 else
00140 {
00141 return l2;
00142 }
00143 }
00144 }
00145 else
00146 {
00147
00148 return (-c/b);
00149 }
00150 }
00151 else
00152 {
00153 vector alpha = from - Cf0;
00154 vector beta = to - from - (Cf - Cf0);
00155 scalar lambdaNominator = alpha & Sf;
00156 scalar lambdaDenominator = beta & Sf;
00157
00158
00159 if (mag(lambdaDenominator) < SMALL)
00160 {
00161 if (lambdaDenominator < 0.0)
00162 {
00163 lambdaDenominator = -SMALL;
00164 }
00165 else
00166 {
00167 lambdaDenominator = SMALL;
00168 }
00169 }
00170
00171 return (-lambdaNominator/lambdaDenominator);
00172 }
00173 }
00174 else
00175 {
00176
00177 return lambda(from, to, facei);
00178 }
00179 }
00180
00181
00182 template<class ParticleType>
00183 inline Foam::scalar Foam::Particle<ParticleType>::lambda
00184 (
00185 const vector& from,
00186 const vector& to,
00187 const label facei
00188 ) const
00189 {
00190 const polyMesh& mesh = cloud_.polyMesh_;
00191
00192 vector Sf = mesh.faceAreas()[facei];
00193 Sf /= mag(Sf);
00194 vector Cf = mesh.faceCentres()[facei];
00195
00196
00197 if (!cloud_.internalFace(facei))
00198 {
00199 label patchi = patch(facei);
00200 const polyPatch& patch = mesh.boundaryMesh()[patchi];
00201
00202
00203 if (isA<wallPolyPatch>(patch))
00204 {
00205 const vector& C = mesh.cellCentres()[celli_];
00206 scalar CCf = mag((C - Cf) & Sf);
00207
00208
00209 const ParticleType& p = static_cast<const ParticleType&>(*this);
00210 if (CCf > p.wallImpactDistance(Sf))
00211 {
00212 Cf -=p.wallImpactDistance(Sf)*Sf;
00213 }
00214 }
00215 }
00216
00217 scalar lambdaNominator = (Cf - from) & Sf;
00218 scalar lambdaDenominator = (to - from) & Sf;
00219
00220
00221 if (mag(lambdaDenominator) < SMALL)
00222 {
00223 if (lambdaDenominator < 0.0)
00224 {
00225 lambdaDenominator = -SMALL;
00226 }
00227 else
00228 {
00229 lambdaDenominator = SMALL;
00230 }
00231 }
00232
00233 return lambdaNominator/lambdaDenominator;
00234 }
00235
00236
00237 template<class ParticleType>
00238 inline bool Foam::Particle<ParticleType>::inCell() const
00239 {
00240 DynamicList<label>& faces = cloud_.labels_;
00241 findFaces(position_, faces);
00242
00243 return (!faces.size());
00244 }
00245
00246
00247 template<class ParticleType>
00248 inline bool Foam::Particle<ParticleType>::inCell
00249 (
00250 const vector& position,
00251 const label celli,
00252 const scalar stepFraction
00253 ) const
00254 {
00255 DynamicList<label>& faces = cloud_.labels_;
00256 findFaces(position, celli, stepFraction, faces);
00257
00258 return (!faces.size());
00259 }
00260
00261
00262
00263
00264 template<class ParticleType>
00265 inline Foam::Particle<ParticleType>::trackData::trackData
00266 (
00267 Cloud<ParticleType>& cloud
00268 )
00269 :
00270 cloud_(cloud)
00271 {}
00272
00273
00274 template<class ParticleType>
00275 inline Foam::Cloud<ParticleType>&
00276 Foam::Particle<ParticleType>::trackData::cloud()
00277 {
00278 return cloud_;
00279 }
00280
00281
00282
00283
00284 template<class ParticleType>
00285 inline const Foam::Cloud<ParticleType>&
00286 Foam::Particle<ParticleType>::cloud() const
00287 {
00288 return cloud_;
00289 }
00290
00291
00292 template<class ParticleType>
00293 inline const Foam::vector& Foam::Particle<ParticleType>::position() const
00294 {
00295 return position_;
00296 }
00297
00298
00299 template<class ParticleType>
00300 inline Foam::vector& Foam::Particle<ParticleType>::position()
00301 {
00302 return position_;
00303 }
00304
00305
00306 template<class ParticleType>
00307 inline Foam::label Foam::Particle<ParticleType>::cell() const
00308 {
00309 return celli_;
00310 }
00311
00312
00313 template<class ParticleType>
00314 inline Foam::label& Foam::Particle<ParticleType>::cell()
00315 {
00316 return celli_;
00317 }
00318
00319
00320 template<class ParticleType>
00321 inline Foam::label Foam::Particle<ParticleType>::face() const
00322 {
00323 return facei_;
00324 }
00325
00326
00327 template<class ParticleType>
00328 inline bool Foam::Particle<ParticleType>::onBoundary() const
00329 {
00330 return facei_ != -1 && facei_ >= cloud_.pMesh().nInternalFaces();
00331 }
00332
00333
00334 template<class ParticleType>
00335 inline Foam::scalar& Foam::Particle<ParticleType>::stepFraction()
00336 {
00337 return stepFraction_;
00338 }
00339
00340
00341 template<class ParticleType>
00342 inline Foam::scalar Foam::Particle<ParticleType>::stepFraction() const
00343 {
00344 return stepFraction_;
00345 }
00346
00347
00348 template<class ParticleType>
00349 inline Foam::label Foam::Particle<ParticleType>::origProc() const
00350 {
00351 return origProc_;
00352 }
00353
00354
00355 template<class ParticleType>
00356 inline Foam::label Foam::Particle<ParticleType>::origId() const
00357 {
00358 return origId_;
00359 }
00360
00361
00362 template<class ParticleType>
00363 inline bool Foam::Particle<ParticleType>::softImpact() const
00364 {
00365 return false;
00366 }
00367
00368
00369 template<class ParticleType>
00370 inline Foam::scalar Foam::Particle<ParticleType>::currentTime() const
00371 {
00372 return
00373 cloud_.pMesh().time().value()
00374 + stepFraction_*cloud_.pMesh().time().deltaT().value();
00375 }
00376
00377
00378 template<class ParticleType>
00379 inline Foam::label Foam::Particle<ParticleType>::patch(const label facei) const
00380 {
00381 return cloud_.facePatch(facei);
00382 }
00383
00384
00385 template<class ParticleType>
00386 inline Foam::label Foam::Particle<ParticleType>::patchFace
00387 (
00388 const label patchi,
00389 const label facei
00390 ) const
00391 {
00392 return cloud_.patchFace(patchi, facei);
00393 }
00394
00395
00396 template<class ParticleType>
00397 inline Foam::scalar
00398 Foam::Particle<ParticleType>::wallImpactDistance(const vector&) const
00399 {
00400 return 0.0;
00401 }
00402
00403
00404 template<class ParticleType>
00405 inline Foam::label Foam::Particle<ParticleType>::faceInterpolation() const
00406 {
00407 return facei_;
00408 }
00409
00410
00411