WallLocalSpringSliderDashpot.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class CloudType>
35(
36 scalar& rMin,
37 scalar& rhoMax,
38 scalar& UMagMax
39) const
40{
41 rMin = VGREAT;
42 rhoMax = -VGREAT;
43 UMagMax = -VGREAT;
44
45 for (const typename CloudType::parcelType& p : this->owner())
46 {
47 // Finding minimum diameter to avoid excessive arithmetic
48
49 scalar dEff = p.d();
50
51 if (useEquivalentSize_)
52 {
53 dEff *= cbrt(p.nParticle()*volumeFactor_);
54 }
55
56 rMin = min(dEff, rMin);
57
58 rhoMax = max(p.rho(), rhoMax);
59
60 UMagMax = max
61 (
62 mag(p.U()) + mag(p.omega())*dEff/2,
63 UMagMax
64 );
65 }
66
67 // Transform the minimum diameter into minimum radius
68 // rMin = dMin/2
69
70 rMin /= 2.0;
71}
72
73
74template<class CloudType>
76(
77 typename CloudType::parcelType& p,
78 const point& site,
79 const WallSiteData<vector>& data,
80 scalar pREff,
81 bool cohesion
82) const
83{
84 // wall patch index
85 label wPI = patchMap_[data.patchIndex()];
86
87 // data for this patch
88 scalar Estar = Estar_[wPI];
89 scalar Gstar = Gstar_[wPI];
90 scalar alpha = alpha_[wPI];
91 scalar b = b_[wPI];
92 scalar mu = mu_[wPI];
93 scalar cohesionEnergyDensity = cohesionEnergyDensity_[wPI];
94 cohesion = cohesion && cohesion_[wPI];
95
96 vector r_PW = p.position() - site;
97
98 vector U_PW = p.U() - data.wallData();
99
100 scalar r_PW_mag = mag(r_PW);
101
102 scalar normalOverlapMag = max(pREff - r_PW_mag, 0.0);
103
104 vector rHat_PW = r_PW/(r_PW_mag + VSMALL);
105
106 scalar kN = (4.0/3.0)*sqrt(pREff)*Estar;
107
108 scalar etaN = alpha*sqrt(p.mass()*kN)*pow025(normalOverlapMag);
109
110 vector fN_PW =
111 rHat_PW
112 *(kN*pow(normalOverlapMag, b) - etaN*(U_PW & rHat_PW));
113
114 // Cohesion force, energy density multiplied by the area of wall/particle
115 // overlap
116 if (cohesion)
117 {
118 fN_PW +=
119 -cohesionEnergyDensity
120 *mathematical::pi*(sqr(pREff) - sqr(r_PW_mag))
121 *rHat_PW;
122 }
123
124 p.f() += fN_PW;
125
126 vector USlip_PW =
127 U_PW - (U_PW & rHat_PW)*rHat_PW
128 + (p.omega() ^ (pREff*-rHat_PW));
129
130 scalar deltaT = this->owner().mesh().time().deltaTValue();
131
132 vector& tangentialOverlap_PW =
133 p.collisionRecords().matchWallRecord(-r_PW, pREff).collisionData();
134
135 tangentialOverlap_PW += USlip_PW*deltaT;
136
137 scalar tangentialOverlapMag = mag(tangentialOverlap_PW);
138
139 if (tangentialOverlapMag > VSMALL)
140 {
141 scalar kT = 8.0*sqrt(pREff*normalOverlapMag)*Gstar;
142
143 scalar etaT = etaN;
144
145 // Tangential force
146 vector fT_PW;
147
148 if (kT*tangentialOverlapMag > mu*mag(fN_PW))
149 {
150 // Tangential force greater than sliding friction,
151 // particle slips
152
153 fT_PW = -mu*mag(fN_PW)*USlip_PW/mag(USlip_PW);
154
155 tangentialOverlap_PW = Zero;
156 }
157 else
158 {
159 fT_PW =
160 -kT*tangentialOverlapMag
161 *tangentialOverlap_PW/tangentialOverlapMag
162 - etaT*USlip_PW;
163 }
164
165 p.f() += fT_PW;
166
167 p.torque() += (pREff*-rHat_PW) ^ fT_PW;
168 }
169}
170
171
172// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173
174template<class CloudType>
176(
177 const dictionary& dict,
179)
180:
181 WallModel<CloudType>(dict, cloud, typeName),
182 Estar_(),
183 Gstar_(),
184 alpha_(),
185 b_(),
186 mu_(),
187 cohesionEnergyDensity_(),
188 cohesion_(),
189 patchMap_(),
190 maxEstarIndex_(-1),
191 collisionResolutionSteps_
192 (
193 this->coeffDict().getScalar("collisionResolutionSteps")
194 ),
195 volumeFactor_(1.0),
196 useEquivalentSize_(Switch(this->coeffDict().lookup("useEquivalentSize")))
197{
198 if (useEquivalentSize_)
199 {
200 this->coeffDict().readEntry("volumeFactor", volumeFactor_);
201 }
202
203 scalar pNu = this->owner().constProps().poissonsRatio();
204
205 scalar pE = this->owner().constProps().youngsModulus();
206
207 const polyMesh& mesh = cloud.mesh();
208
210
211 patchMap_.setSize(bMesh.size(), -1);
212
213 DynamicList<label> wallPatchIndices;
214
215 for (const polyPatch& pp : bMesh)
216 {
217 if (isA<wallPolyPatch>(pp))
218 {
219 wallPatchIndices.append(pp.index());
220 }
221 }
222
223 label nWallPatches = wallPatchIndices.size();
224
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);
232
233 scalar maxEstar = -GREAT;
234
235 forAll(wallPatchIndices, wPI)
236 {
237 const dictionary& patchCoeffDict
238 (
239 this->coeffDict().subDict(bMesh[wallPatchIndices[wPI]].name())
240 );
241
242 patchMap_[wallPatchIndices[wPI]] = wPI;
243
244 scalar nu = patchCoeffDict.get<scalar>("poissonsRatio");
245
246 scalar E = patchCoeffDict.get<scalar>("youngsModulus");
247
248 Estar_[wPI] = 1/((1 - sqr(pNu))/pE + (1 - sqr(nu))/E);
249
250 Gstar_[wPI] = 1/(2*((2 + pNu - sqr(pNu))/pE + (2 + nu - sqr(nu))/E));
251
252 alpha_[wPI] = patchCoeffDict.get<scalar>("alpha");
253
254 b_[wPI] = patchCoeffDict.get<scalar>("b");
255
256 mu_[wPI] = patchCoeffDict.get<scalar>("mu");
257
258 cohesionEnergyDensity_[wPI] =
259 patchCoeffDict.get<scalar>("cohesionEnergyDensity");
260
261 cohesion_[wPI] = (mag(cohesionEnergyDensity_[wPI]) > VSMALL);
262
263 if (Estar_[wPI] > maxEstar)
264 {
265 maxEstarIndex_ = wPI;
266
267 maxEstar = Estar_[wPI];
268 }
269 }
270}
271
272
273// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
274
275template<class CloudType>
277(
278 const typename CloudType::parcelType& p
279) const
280{
281 if (useEquivalentSize_)
282 {
283 return p.d()/2*cbrt(p.nParticle()*volumeFactor_);
284 }
285
286 return p.d()/2;
287}
288
289
290template<class CloudType>
292{
293 return true;
294}
295
296
297template<class CloudType>
299{
300 if (!(this->owner().size()))
301 {
302 return 1;
303 }
304
305 scalar rMin;
306 scalar rhoMax;
307 scalar UMagMax;
308
309 findMinMaxProperties(rMin, rhoMax, UMagMax);
310
311 // Note: pi^(7/5)*(5/4)^(2/5) = 5.429675
312 scalar minCollisionDeltaT =
313 5.429675
314 *rMin
315 *pow(rhoMax/(Estar_[maxEstarIndex_]*sqrt(UMagMax) + VSMALL), 0.4)
316 /collisionResolutionSteps_;
317
318 return ceil(this->owner().time().deltaTValue()/minCollisionDeltaT);
319}
320
321
322template<class CloudType>
324(
325 typename CloudType::parcelType& p,
326 const List<point>& flatSitePoints,
327 const List<WallSiteData<vector>>& flatSiteData,
328 const List<point>& sharpSitePoints,
329 const List<WallSiteData<vector>>& sharpSiteData
330) const
331{
332 scalar pREff = this->pREff(p);
333
334 forAll(flatSitePoints, siteI)
335 {
336 evaluateWall
337 (
338 p,
339 flatSitePoints[siteI],
340 flatSiteData[siteI],
341 pREff,
342 true
343 );
344 }
345
346 forAll(sharpSitePoints, siteI)
347 {
348 // Treating sharp sites like flat sites, except suppress cohesion
349
350 evaluateWall
351 (
352 p,
353 sharpSitePoints[siteI],
354 sharpSiteData[siteI],
355 pREff,
356 false
357 );
358 }
359}
360
361
362// ************************************************************************* //
Y[inertIndex] max(0.0)
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:98
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A list of faces which address into the list of points.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Forces between particles and walls, interacting with a spring, slider, damper model.
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.
virtual bool controlsTimestep() const
Whether the WallModel has a timestep limit that will.
Templated wall interaction class.
Definition: WallModel.H:56
const dictionary & coeffDict() const
Return the coefficients dictionary.
Definition: WallModel.C:80
const CloudType & owner() const
Return the owner cloud object.
Definition: WallModel.C:57
Stores the patch ID and templated data to represent a collision with a wall to be passed to the wall ...
Definition: WallSiteData.H:67
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
const volScalarField & mu
dynamicFvMesh & mesh
const dimensionedScalar rhoMax
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & nu
volScalarField & alpha
dictionary dict
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333