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);
57 rhoMax =
max(p.rho(), rhoMax);
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 =
Zero;
159 -kT*tangentialOverlapMag
160 *tangentialOverlap_PW/tangentialOverlapMag
166 p.torque() += (pREff*-rHat_PW) ^ fT_PW;
173 template<
class CloudType>
186 cohesionEnergyDensity_(),
190 collisionResolutionSteps_
192 this->coeffDict().
template lookup<scalar>(
"collisionResolutionSteps")
195 useEquivalentSize_(
Switch(this->coeffDict().
lookup(
"useEquivalentSize")))
197 if (useEquivalentSize_)
200 this->coeffDict().template lookup<scalar>(
"volumeFactor");
203 scalar pNu = this->owner().constProps().poissonsRatio();
205 scalar pE = this->owner().constProps().youngsModulus();
217 if (isA<wallPolyPatch>(bMesh[
patchi]))
219 wallPatchIndices.
append(bMesh[patchi].index());
223 label nWallPatches = wallPatchIndices.
size();
225 Estar_.setSize(nWallPatches);
226 Gstar_.setSize(nWallPatches);
227 alpha_.setSize(nWallPatches);
228 b_.setSize(nWallPatches);
229 mu_.setSize(nWallPatches);
230 cohesionEnergyDensity_.setSize(nWallPatches);
231 cohesion_.setSize(nWallPatches);
233 scalar maxEstar = -great;
235 forAll(wallPatchIndices, wPI)
239 this->coeffDict().subDict(bMesh[wallPatchIndices[wPI]].
name())
242 patchMap_[wallPatchIndices[wPI]] = wPI;
244 scalar nu = patchCoeffDict.template lookup<scalar>(
"poissonsRatio");
246 scalar E = patchCoeffDict.template lookup<scalar>(
"youngsModulus");
248 Estar_[wPI] = 1/((1 -
sqr(pNu))/pE + (1 -
sqr(nu))/E);
250 Gstar_[wPI] = 1/(2*((2 + pNu -
sqr(pNu))/pE + (2 + nu -
sqr(nu))/E));
252 alpha_[wPI] = patchCoeffDict.template lookup<scalar>(
"alpha");
254 b_[wPI] = patchCoeffDict.template lookup<scalar>(
"b");
256 mu_[wPI] = patchCoeffDict.template lookup<scalar>(
"mu");
258 cohesionEnergyDensity_[wPI] =
259 patchCoeffDict.
lookup<scalar>(
"cohesionEnergyDensity");
261 cohesion_[wPI] = (
mag(cohesionEnergyDensity_[wPI]) > vSmall);
263 if (Estar_[wPI] > maxEstar)
265 maxEstarIndex_ = wPI;
267 maxEstar = Estar_[wPI];
275 template<
class CloudType>
282 template<
class CloudType>
288 if (useEquivalentSize_)
290 return p.d()/2*
cbrt(p.nParticle()*volumeFactor_);
299 template<
class CloudType>
306 template<
class CloudType>
309 if (!(this->owner().size()))
318 findMinMaxProperties(rMin, rhoMax, UMagMax);
321 scalar minCollisionDeltaT =
324 *
pow(rhoMax/(Estar_[maxEstarIndex_]*
sqrt(UMagMax) + vSmall), 0.4)
325 /collisionResolutionSteps_;
327 return ceil(this->owner().time().deltaTValue()/minCollisionDeltaT);
331 template<
class CloudType>
341 scalar pREff = this->pREff(p);
343 forAll(flatSitePoints, siteI)
348 flatSitePoints[siteI],
355 forAll(sharpSitePoints, siteI)
362 sharpSitePoints[siteI],
363 sharpSiteData[siteI],
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
WallLocalSpringSliderDashpot(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
DSMCCloud< dsmcParcel > CloudType
virtual label nSubCycles() const
For WallModels that control the timestep, calculate the.
virtual scalar pREff(const typename CloudType::parcelType &p) const
Return the effective radius for a particle for the model.
A list of keyword definitions, which are a keyword followed by any number of values (e...
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Forces between particles and walls, interacting with a spring, slider, damper model.
dimensionedScalar pow025(const dimensionedScalar &ds)
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/any.
Vector< scalar > vector
A scalar version of the templated Vector.
virtual bool controlsTimestep() const
Whether the WallModel has a timestep limit that will.
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar cbrt(const dimensionedScalar &ds)
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Templated wall interaction class.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Stores the patch ID and templated data to represent a collision with a wall to be passed to the wall ...
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
label size() const
Return the number of elements in the UPtrList.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
vector point
Point is a vector.
const polyBoundaryMesh & bMesh
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Templated base class for dsmc cloud.
virtual ~WallLocalSpringSliderDashpot()
Destructor.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.