Loading...
Searching...
No Matches
WallSpringSliderDashpot.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-2017 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>
34void Foam::WallSpringSliderDashpot<CloudType>::findMinMaxProperties
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>
75void Foam::WallSpringSliderDashpot<CloudType>::evaluateWall
76(
77 typename CloudType::parcelType& p,
78 const point& site,
79 const WallSiteData<vector>& data,
80 scalar pREff,
81 scalar kN,
82 bool cohesion
83) const
84{
85 vector r_PW = p.position() - site;
86
87 vector U_PW = p.U() - data.wallData();
88
89 scalar r_PW_mag = mag(r_PW);
90
91 scalar normalOverlapMag = max(pREff - r_PW_mag, 0.0);
92
93 vector rHat_PW = r_PW/(r_PW_mag + VSMALL);
94
95 scalar etaN = alpha_*sqrt(p.mass()*kN)*pow025(normalOverlapMag);
96
97 vector fN_PW =
98 rHat_PW
99 *(kN*pow(normalOverlapMag, b_) - etaN*(U_PW & rHat_PW));
100
101 // Cohesion force, energy density multiplied by the area of wall/particle
102 // overlap
103 if (cohesion)
104 {
105 fN_PW +=
106 -cohesionEnergyDensity_
107 *mathematical::pi*(sqr(pREff) - sqr(r_PW_mag))
108 *rHat_PW;
109 }
110
111 p.f() += fN_PW;
112
113 vector USlip_PW =
114 U_PW - (U_PW & rHat_PW)*rHat_PW
115 + (p.omega() ^ (pREff*-rHat_PW));
116
117 scalar deltaT = this->owner().mesh().time().deltaTValue();
118
119 vector& tangentialOverlap_PW =
120 p.collisionRecords().matchWallRecord(-r_PW, pREff).collisionData();
121
122 tangentialOverlap_PW += USlip_PW*deltaT;
123
124 scalar tangentialOverlapMag = mag(tangentialOverlap_PW);
125
126 if (tangentialOverlapMag > VSMALL)
127 {
128 scalar kT = 8.0*sqrt(pREff*normalOverlapMag)*Gstar_;
129
130 scalar etaT = etaN;
131
132 // Tangential force
133 vector fT_PW;
134
135 if (kT*tangentialOverlapMag > mu_*mag(fN_PW))
136 {
137 // Tangential force greater than sliding friction,
138 // particle slips
139
140 fT_PW = -mu_*mag(fN_PW)*USlip_PW/mag(USlip_PW);
141
142 tangentialOverlap_PW = Zero;
143 }
144 else
145 {
146 fT_PW = - kT*tangentialOverlap_PW - etaT*USlip_PW;
147 }
148
149 p.f() += fT_PW;
150
151 p.torque() += (pREff*-rHat_PW) ^ fT_PW;
153}
154
155
156// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
157
158template<class CloudType>
160(
161 const dictionary& dict,
162 CloudType& cloud
163)
164:
165 WallModel<CloudType>(dict, cloud, typeName),
166 Estar_(),
167 Gstar_(),
168 alpha_(this->coeffDict().getScalar("alpha")),
169 b_(this->coeffDict().getScalar("b")),
170 mu_(this->coeffDict().getScalar("mu")),
171 cohesionEnergyDensity_
172 (
173 this->coeffDict().getScalar("cohesionEnergyDensity")
174 ),
175 cohesion_(false),
176 collisionResolutionSteps_
177 (
178 this->coeffDict().getScalar("collisionResolutionSteps")
179 ),
180 volumeFactor_(1.0),
181 useEquivalentSize_(Switch(this->coeffDict().lookup("useEquivalentSize")))
182{
183 if (useEquivalentSize_)
184 {
185 this->coeffDict().readEntry("volumeFactor", volumeFactor_);
186 }
187
188 scalar nu = this->coeffDict().getScalar("poissonsRatio");
189
190 scalar E = this->coeffDict().getScalar("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
206template<class CloudType>
208{}
209
210
211// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212
213template<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_);
223
224 return p.d()/2;
225}
226
227
228template<class CloudType>
230{
231 return true;
232}
233
234
235template<class CloudType>
237{
238 if (!(this->owner().size()))
239 {
240 return 1;
241 }
242
243 scalar rMin;
244 scalar rhoMax;
245 scalar UMagMax;
246
247 findMinMaxProperties(rMin, rhoMax, UMagMax);
248
249 // Note: pi^(7/5)*(5/4)^(2/5) = 5.429675
250 scalar minCollisionDeltaT =
251 5.429675
252 *rMin
253 *pow(rhoMax/(Estar_*sqrt(UMagMax) + VSMALL), 0.4)
254 /collisionResolutionSteps_;
255
256 return ceil(this->owner().time().deltaTValue()/minCollisionDeltaT);
257}
258
259
260template<class CloudType>
261void Foam::WallSpringSliderDashpot<CloudType>::evaluateWall
262(
263 typename CloudType::parcelType& p,
264 const List<point>& flatSitePoints,
265 const List<WallSiteData<vector>>& flatSiteData,
266 const List<point>& sharpSitePoints,
267 const List<WallSiteData<vector>>& sharpSiteData
268) const
269{
270 scalar pREff = this->pREff(p);
271
272 scalar kN = (4.0/3.0)*sqrt(pREff)*Estar_;
273
274 forAll(flatSitePoints, siteI)
275 {
276 evaluateWall
277 (
278 p,
279 flatSitePoints[siteI],
280 flatSiteData[siteI],
281 pREff,
282 kN,
283 cohesion_
284 );
285 }
286
287 forAll(sharpSitePoints, siteI)
288 {
289 // Treating sharp sites like flat sites, except suppress cohesion
290
291 evaluateWall
292 (
293 p,
294 sharpSitePoints[siteI],
295 sharpSiteData[siteI],
296 pREff,
297 kN,
298 false
299 );
300 }
301}
302
303
304// ************************************************************************* //
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition DSMCCloudI.H:91
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
WallModel(const dictionary &dict, CloudType &owner, const word &type)
Construct from components.
Definition WallModel.C:27
const dictionary & coeffDict() const
Return the coefficients dictionary.
Definition WallModel.C:73
const dictionary & dict() const
Return the dictionary.
Definition WallModel.C:65
const CloudType & owner() const
Return the owner cloud object.
Definition WallModel.C:50
Stores the patch ID and templated data to represent a collision with a wall to be passed to the wall ...
virtual ~WallSpringSliderDashpot()
Destructor.
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.
WallSpringSliderDashpot(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
virtual bool controlsTimestep() const
Whether the WallModel has a timestep limit that will.
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
scalar getScalar(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<scalar>(const word&, keyType::option).
Lookup type of boundary radiation properties.
Definition lookup.H:60
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
const dimensionedScalar rhoMax
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
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:26
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & nu
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299