solidParticle.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-2020 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 
26 #include "solidParticleCloud.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(solidParticle, 0);
33 }
34 
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
39 (
40  solidParticleCloud& cloud,
41  trackingData& td,
42  const scalar trackTime
43 )
44 {
45  td.switchProcessor = false;
46  td.keepParticle = true;
47 
48  while (td.keepParticle && !td.switchProcessor && stepFraction() < 1)
49  {
50  if (debug)
51  {
52  Info<< "Time = " << mesh().time().timeName()
53  << " trackTime = " << trackTime
54  << " stepFraction() = " << stepFraction() << endl;
55  }
56 
57 
58  const scalar sfrac = stepFraction();
59 
60  const scalar f = 1 - stepFraction();
61  trackToAndHitFace(f*trackTime*U_, f, cloud, td);
62 
63  const scalar dt = (stepFraction() - sfrac)*trackTime;
64 
65  const tetIndices tetIs = this->currentTetIndices();
66  scalar rhoc = td.rhoInterp().interpolate(this->coordinates(), tetIs);
67  vector Uc = td.UInterp().interpolate(this->coordinates(), tetIs);
68  scalar nuc = td.nuInterp().interpolate(this->coordinates(), tetIs);
69 
70  scalar rhop = cloud.rhop();
71  scalar magUr = mag(Uc - U_);
72 
73  scalar ReFunc = 1.0;
74  scalar Re = magUr*d_/nuc;
75 
76  if (Re > 0.01)
77  {
78  ReFunc += 0.15*pow(Re, 0.687);
79  }
80 
81  scalar Dc = (24.0*nuc/d_)*ReFunc*(3.0/4.0)*(rhoc/(d_*rhop));
82 
83  U_ = (U_ + dt*(Dc*Uc + (1.0 - rhoc/rhop)*td.g()))/(1.0 + dt*Dc);
84  }
85 
86  return td.keepParticle;
87 }
88 
89 
91 {
92  return false;
93 }
94 
95 
97 (
99  trackingData& td
100 )
101 {
102  td.switchProcessor = true;
103 }
104 
105 
107 {
108  const vector nw = normal();
109 
110  const scalar Un = U_ & nw;
111  const vector Ut = U_ - Un*nw;
112 
113  if (Un > 0)
114  {
115  U_ -= (1.0 + cloud.e())*Un*nw;
116  }
117 
118  U_ -= cloud.mu()*Ut;
119 }
120 
121 
123 {
125  U_ = transform.transform(U_);
126 }
127 
128 
129 // ************************************************************************* //
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool hitPatch(solidParticleCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
Definition: solidParticle.C:90
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
Definition: particle.C:1048
dynamicFvMesh & mesh
bool switchProcessor
Flag to switch processor.
Definition: particle.H:109
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
A Cloud of solid particles.
const interpolationCellPoint< scalar > & rhoInterp() const
const interpolationCellPoint< scalar > & nuInterp() const
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
const interpolationCellPoint< vector > & UInterp() const
Class used to pass tracking data to the trackToFace function.
Definition: solidParticle.H:87
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void hitProcessorPatch(solidParticleCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
Definition: solidParticle.C:97
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:112
bool move(solidParticleCloud &, trackingData &, const scalar)
Move.
Definition: solidParticle.C:39
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Type transform(const Type &) const
Transform the given type.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
void hitWallPatch(solidParticleCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
label nw
Definition: createFields.H:12
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477