WallSpringSliderDashpot.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class CloudType>
32 (
33  scalar& rMin,
34  scalar& rhoMax,
35  scalar& UMagMax
36 ) const
37 {
38  rMin = vGreat;
39  rhoMax = -vGreat;
40  UMagMax = -vGreat;
41 
42  forAllConstIter(typename CloudType, this->owner(), iter)
43  {
44  const typename CloudType::parcelType& p = iter();
45 
46  // Finding minimum diameter to avoid excessive arithmetic
47 
48  scalar dEff = p.d();
49 
50  if (useEquivalentSize_)
51  {
52  dEff *= cbrt(p.nParticle()*volumeFactor_);
53  }
54 
55  rMin = min(dEff, rMin);
56 
57  rhoMax = max(p.rho(), rhoMax);
58 
59  UMagMax = max
60  (
61  mag(p.U()) + mag(p.omega())*dEff/2,
62  UMagMax
63  );
64  }
65 
66  // Transform the minimum diameter into minimum radius
67  // rMin = dMin/2
68 
69  rMin /= 2.0;
70 }
71 
72 
73 template<class CloudType>
75 (
76  typename CloudType::parcelType& p,
77  const point& site,
78  const WallSiteData<vector>& data,
79  scalar pREff,
80  scalar kN,
81  bool cohesion
82 ) const
83 {
84  vector r_PW = p.position() - site;
85 
86  vector U_PW = p.U() - data.wallData();
87 
88  scalar r_PW_mag = mag(r_PW);
89 
90  scalar normalOverlapMag = max(pREff - r_PW_mag, 0.0);
91 
92  vector rHat_PW = r_PW/(r_PW_mag + vSmall);
93 
94  scalar etaN = alpha_*sqrt(p.mass()*kN)*pow025(normalOverlapMag);
95 
96  vector fN_PW =
97  rHat_PW
98  *(kN*pow(normalOverlapMag, b_) - etaN*(U_PW & rHat_PW));
99 
100  // Cohesion force, energy density multiplied by the area of wall/particle
101  // overlap
102  if (cohesion)
103  {
104  fN_PW +=
105  -cohesionEnergyDensity_
106  *mathematical::pi*(sqr(pREff) - sqr(r_PW_mag))
107  *rHat_PW;
108  }
109 
110  p.f() += fN_PW;
111 
112  vector USlip_PW =
113  U_PW - (U_PW & rHat_PW)*rHat_PW
114  + (p.omega() ^ (pREff*-rHat_PW));
115 
116  scalar deltaT = this->owner().mesh().time().deltaTValue();
117 
118  vector& tangentialOverlap_PW =
119  p.collisionRecords().matchWallRecord(-r_PW, pREff).collisionData();
120 
121  tangentialOverlap_PW += USlip_PW*deltaT;
122 
123  scalar tangentialOverlapMag = mag(tangentialOverlap_PW);
124 
125  if (tangentialOverlapMag > vSmall)
126  {
127  scalar kT = 8.0*sqrt(pREff*normalOverlapMag)*Gstar_;
128 
129  scalar etaT = etaN;
130 
131  // Tangential force
132  vector fT_PW;
133 
134  if (kT*tangentialOverlapMag > mu_*mag(fN_PW))
135  {
136  // Tangential force greater than sliding friction,
137  // particle slips
138 
139  fT_PW = -mu_*mag(fN_PW)*USlip_PW/mag(USlip_PW);
140 
141  tangentialOverlap_PW = Zero;
142  }
143  else
144  {
145  fT_PW = - kT*tangentialOverlap_PW - etaT*USlip_PW;
146  }
147 
148  p.f() += fT_PW;
149 
150  p.torque() += (pREff*-rHat_PW) ^ fT_PW;
151  }
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
156 
157 template<class CloudType>
159 (
160  const dictionary& dict,
161  CloudType& cloud
162 )
163 :
164  WallModel<CloudType>(dict, cloud, typeName),
165  Estar_(),
166  Gstar_(),
167  alpha_(this->coeffDict().template lookup<scalar>("alpha")),
168  b_(this->coeffDict().template lookup<scalar>("b")),
169  mu_(this->coeffDict().template lookup<scalar>("mu")),
170  cohesionEnergyDensity_
171  (
172  this->coeffDict().template lookup<scalar>("cohesionEnergyDensity")
173  ),
174  cohesion_(false),
175  collisionResolutionSteps_
176  (
177  this->coeffDict().template lookup<scalar>("collisionResolutionSteps")
178  ),
179  volumeFactor_(1.0),
180  useEquivalentSize_(Switch(this->coeffDict().lookup("useEquivalentSize")))
181 {
182  if (useEquivalentSize_)
183  {
184  volumeFactor_ =
185  this->coeffDict().template lookup<scalar>("volumeFactor");
186  }
187 
188  scalar nu = this->coeffDict().template lookup<scalar>("poissonsRatio");
189 
190  scalar E = this->coeffDict().template lookup<scalar>("youngsModulus");
191 
192  scalar pNu = this->owner().constProps().poissonsRatio();
193 
194  scalar pE = this->owner().constProps().youngsModulus();
195 
196  Estar_ = 1/((1 - sqr(pNu))/pE + (1 - sqr(nu))/E);
197 
198  Gstar_ = 1/(2*((2 + pNu - sqr(pNu))/pE + (2 + nu - sqr(nu))/E));
199 
200  cohesion_ = (mag(cohesionEnergyDensity_) > vSmall);
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
205 
206 template<class CloudType>
208 {}
209 
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
213 template<class CloudType>
215 (
216  const typename CloudType::parcelType& p
217 ) const
218 {
219  if (useEquivalentSize_)
220  {
221  return p.d()/2*cbrt(p.nParticle()*volumeFactor_);
222  }
223  else
224  {
225  return p.d()/2;
226  }
227 }
228 
229 
230 template<class CloudType>
232 {
233  return true;
234 }
235 
236 
237 template<class CloudType>
239 {
240  if (!(this->owner().size()))
241  {
242  return 1;
243  }
244 
245  scalar rMin;
246  scalar rhoMax;
247  scalar UMagMax;
248 
249  findMinMaxProperties(rMin, rhoMax, UMagMax);
250 
251  // Note: pi^(7/5)*(5/4)^(2/5) = 5.429675
252  scalar minCollisionDeltaT =
253  5.429675
254  *rMin
255  *pow(rhoMax/(Estar_*sqrt(UMagMax) + vSmall), 0.4)
256  /collisionResolutionSteps_;
257 
258  return ceil(this->owner().time().deltaTValue()/minCollisionDeltaT);
259 }
260 
261 
262 template<class CloudType>
264 (
265  typename CloudType::parcelType& p,
266  const List<point>& flatSitePoints,
267  const List<WallSiteData<vector>>& flatSiteData,
268  const List<point>& sharpSitePoints,
269  const List<WallSiteData<vector>>& sharpSiteData
270 ) const
271 {
272  scalar pREff = this->pREff(p);
273 
274  scalar kN = (4.0/3.0)*sqrt(pREff)*Estar_;
275 
276  forAll(flatSitePoints, siteI)
277  {
278  evaluateWall
279  (
280  p,
281  flatSitePoints[siteI],
282  flatSiteData[siteI],
283  pREff,
284  kN,
285  cohesion_
286  );
287  }
288 
289  forAll(sharpSitePoints, siteI)
290  {
291  // Treating sharp sites like flat sites, except suppress cohesion
292 
293  evaluateWall
294  (
295  p,
296  sharpSitePoints[siteI],
297  sharpSiteData[siteI],
298  pREff,
299  kN,
300  false
301  );
302  }
303 }
304 
305 
306 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
DSMCCloud< dsmcParcel > CloudType
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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...
Definition: HashTable.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
Definition: Switch.H:60
Forces between particles and walls, interacting with a spring, slider, damper model.
virtual label nSubCycles() const
For WallModels that control the timestep, calculate the.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual ~WallSpringSliderDashpot()
Destructor.
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar cbrt(const dimensionedScalar &ds)
static const zero Zero
Definition: zero.H:97
Templated wall interaction class.
Definition: PairCollision.H:51
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Stores the patch ID and templated data to represent a collision with a wall to be passed to the wall ...
Definition: WallSiteData.H:50
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:215
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WallSpringSliderDashpot(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
vector point
Point is a vector.
Definition: point.H:41
dimensioned< scalar > mag(const dimensioned< Type > &)
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
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.