30 template<
class CloudType>
44 const typename CloudType::parcelType& p = iter();
50 if (useEquivalentSize_)
52 dEff *=
cbrt(p.nParticle()*volumeFactor_);
55 rMin =
min(dEff, rMin);
61 mag(p.U()) +
mag(p.omega())*dEff/2,
73 template<
class CloudType>
76 typename CloudType::parcelType& p,
78 const WallSiteData<vector>& data,
84 label wPI = patchMap_[data.patchIndex()];
87 scalar Estar = Estar_[wPI];
88 scalar Gstar = Gstar_[wPI];
89 scalar alpha = alpha_[wPI];
92 scalar cohesionEnergyDensity = cohesionEnergyDensity_[wPI];
93 cohesion = cohesion && cohesion_[wPI];
95 vector r_PW = p.position() - site;
97 vector U_PW = p.U() - data.wallData();
99 scalar r_PW_mag =
mag(r_PW);
101 scalar normalOverlapMag =
max(pREff - r_PW_mag, 0.0);
103 vector rHat_PW = r_PW/(r_PW_mag + VSMALL);
105 scalar kN = (4.0/3.0)*
sqrt(pREff)*Estar;
107 scalar etaN = alpha*
sqrt(p.mass()*kN)*
pow025(normalOverlapMag);
111 *(kN*
pow(normalOverlapMag, b) - etaN*(U_PW & rHat_PW));
118 -cohesionEnergyDensity
126 U_PW - (U_PW & rHat_PW)*rHat_PW
127 + (p.omega() ^ (pREff*-rHat_PW));
129 scalar deltaT = this->owner().mesh().time().deltaTValue();
131 vector& tangentialOverlap_PW =
132 p.collisionRecords().matchWallRecord(-r_PW, pREff).collisionData();
134 tangentialOverlap_PW += USlip_PW*deltaT;
136 scalar tangentialOverlapMag =
mag(tangentialOverlap_PW);
138 if (tangentialOverlapMag > VSMALL)
140 scalar kT = 8.0*
sqrt(pREff*normalOverlapMag)*Gstar;
147 if (kT*tangentialOverlapMag > mu*
mag(fN_PW))
152 fT_PW = -mu*
mag(fN_PW)*USlip_PW/
mag(USlip_PW);
154 tangentialOverlap_PW = vector::zero;
159 -kT*tangentialOverlapMag
160 *tangentialOverlap_PW/tangentialOverlapMag
166 p.torque() += (pREff*-rHat_PW) ^ fT_PW;
173 template<
class CloudType>
186 cohesionEnergyDensity_(),
190 collisionResolutionSteps_
194 this->coeffDict().
lookup(
"collisionResolutionSteps")
198 useEquivalentSize_(
Switch(this->coeffDict().
lookup(
"useEquivalentSize")))
200 if (useEquivalentSize_)
205 scalar pNu = this->owner().constProps().poissonsRatio();
207 scalar pE = this->owner().constProps().youngsModulus();
219 if (isA<wallPolyPatch>(bMesh[patchI]))
221 wallPatchIndices.
append(bMesh[patchI].index());
225 label nWallPatches = wallPatchIndices.
size();
227 Estar_.setSize(nWallPatches);
228 Gstar_.setSize(nWallPatches);
229 alpha_.setSize(nWallPatches);
230 b_.setSize(nWallPatches);
231 mu_.setSize(nWallPatches);
232 cohesionEnergyDensity_.setSize(nWallPatches);
233 cohesion_.setSize(nWallPatches);
235 scalar maxEstar = -GREAT;
237 forAll(wallPatchIndices, wPI)
241 this->coeffDict().subDict(bMesh[wallPatchIndices[wPI]].
name())
244 patchMap_[wallPatchIndices[wPI]] = wPI;
250 Estar_[wPI] = 1/((1 -
sqr(pNu))/pE + (1 -
sqr(nu))/E);
252 Gstar_[wPI] = 1/(2*((2 + pNu -
sqr(pNu))/pE + (2 + nu -
sqr(nu))/E));
262 patchCoeffDict.
lookup(
"cohesionEnergyDensity")
265 cohesion_[wPI] = (
mag(cohesionEnergyDensity_[wPI]) > VSMALL);
267 if (Estar_[wPI] > maxEstar)
269 maxEstarIndex_ = wPI;
271 maxEstar = Estar_[wPI];
279 template<
class CloudType>
286 template<
class CloudType>
292 if (useEquivalentSize_)
294 return p.d()/2*
cbrt(p.nParticle()*volumeFactor_);
303 template<
class CloudType>
310 template<
class CloudType>
313 if (!(this->owner().size()))
322 findMinMaxProperties(rMin, rhoMax, UMagMax);
325 scalar minCollisionDeltaT =
328 *
pow(rhoMax/(Estar_[maxEstarIndex_]*
sqrt(UMagMax) + VSMALL), 0.4)
329 /collisionResolutionSteps_;
331 return ceil(this->owner().time().deltaTValue()/minCollisionDeltaT);
335 template<
class CloudType>
345 scalar pREff = this->pREff(p);
347 forAll(flatSitePoints, siteI)
352 flatSitePoints[siteI],
359 forAll(sharpSitePoints, siteI)
366 sharpSitePoints[siteI],
367 sharpSiteData[siteI],
dimensionedScalar sqrt(const dimensionedScalar &ds)
Templated wall interaction class.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
label size() const
Return the number of elements in the PtrList.
vector point
Point is a vector.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Stores the patch ID and templated data to represent a collision with a wall to be passed to the wall ...
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedScalar cbrt(const dimensionedScalar &ds)
void size(const label)
Override size to be inconsistent with allocated storage.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
DSMCCloud< dsmcParcel > CloudType
A list of keyword definitions, which are a keyword followed by any number of values (e...
ParcelType parcelType
Type of parcel the cloud was instantiated for.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual scalar pREff(const typename CloudType::parcelType &p) const
Return the effective radius for a particle for the model.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual label nSubCycles() const
For WallModels that control the timestep, calculate the.
Forces between particles and walls, interacting with a spring, slider, damper model.
virtual bool controlsTimestep() const
Whether the WallModel has a timestep limit that will.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
A list of faces which address into the list of points.
WallLocalSpringSliderDashpot(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Mesh consisting of general polyhedral cells.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Vector< scalar > vector
A scalar version of the templated Vector.
virtual ~WallLocalSpringSliderDashpot()
Destructor.
Templated base class for dsmc cloud.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())