solidParticle.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "solidParticleCloud.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTemplateTypeNameAndDebug(Cloud<solidParticle>, 0);
35 }
36 
37 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38 
40 (
42  trackingData& td,
43  const scalar trackTime
44 )
45 {
46  td.switchProcessor = false;
47  td.keepParticle = true;
48 
49  while (td.keepParticle && !td.switchProcessor && stepFraction() < 1)
50  {
51  if (debug)
52  {
53  Pout<< "Time = " << mesh().time().timeName()
54  << " trackTime = " << trackTime
55  << " steptFraction() = " << stepFraction() << endl;
56  }
57 
58 
59  const scalar sfrac = stepFraction();
60 
61  const scalar f = 1 - stepFraction();
62  trackToAndHitFace(f*trackTime*U_, f, cloud, td);
63 
64  const scalar dt = (stepFraction() - sfrac)*trackTime;
65 
66  const tetIndices tetIs = this->currentTetIndices();
67  scalar rhoc = td.rhoInterp().interpolate(this->coordinates(), tetIs);
68  vector Uc = td.UInterp().interpolate(this->coordinates(), tetIs);
69  scalar nuc = td.nuInterp().interpolate(this->coordinates(), tetIs);
70 
71  scalar rhop = cloud.rhop();
72  scalar magUr = mag(Uc - U_);
73 
74  scalar ReFunc = 1.0;
75  scalar Re = magUr*d_/nuc;
76 
77  if (Re > 0.01)
78  {
79  ReFunc += 0.15*pow(Re, 0.687);
80  }
81 
82  scalar Dc = (24.0*nuc/d_)*ReFunc*(3.0/4.0)*(rhoc/(d_*rhop));
83 
84  U_ = (U_ + dt*(Dc*Uc + (1.0 - rhoc/rhop)*td.g()))/(1.0 + dt*Dc);
85  }
86 
87  return td.keepParticle;
88 }
89 
90 
92 {
93  return false;
94 }
95 
96 
98 (
100  trackingData& td
101 )
102 {
103  td.switchProcessor = true;
104 }
105 
106 
108 {
109  const vector nw = normal();
110 
111  scalar Un = U_ & nw;
112  vector Ut = U_ - Un*nw;
113 
114  if (Un > 0)
115  {
116  U_ -= (1.0 + cloud.e())*Un*nw;
117  }
118 
119  U_ -= cloud.mu()*Ut;
120 }
121 
122 
124 {
126  U_ = transform(T, U_);
127 }
128 
129 
131 {
132  particle::transformProperties(separation);
133 }
134 
135 
136 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::Tensor< scalar >
Foam::solidParticleCloud
A Cloud of solid particles.
Definition: solidParticleCloud.H:58
Foam::solidParticle::trackingData::g
const vector & g() const
Definition: solidParticleI.H:100
Foam::solidParticle::hitWallPatch
void hitWallPatch(solidParticleCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
Definition: solidParticle.C:107
Foam::defineTemplateTypeNameAndDebug
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::solidParticle::transformProperties
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: solidParticle.C:123
Foam::solidParticle::trackingData::UInterp
const interpolationCellPoint< vector > & UInterp() const
Definition: solidParticleI.H:88
Foam::solidParticle::trackingData::nuInterp
const interpolationCellPoint< scalar > & nuInterp() const
Definition: solidParticleI.H:95
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::particle::trackingData::keepParticle
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:105
Foam::particle::trackingData::switchProcessor
bool switchProcessor
Flag to switch processor.
Definition: particle.H:102
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::interpolationCellPoint::interpolate
Type interpolate(const cellPointWeight &cpw) const
Interpolate field for the given cellPointWeight.
Definition: interpolationCellPointI.H:32
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::solidParticle::trackingData
Class used to pass tracking data to the trackToFace function.
Definition: solidParticle.H:83
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
solidParticleCloud.H
f
labelList f(nPoints)
Foam::particle::transformProperties
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:1072
Foam::Vector< scalar >
Foam::solidParticle::move
bool move(solidParticleCloud &, trackingData &, const scalar)
Move.
Definition: solidParticle.C:40
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::solidParticle::trackingData::rhoInterp
const interpolationCellPoint< scalar > & rhoInterp() const
Definition: solidParticleI.H:81
Foam::solidParticle::hitPatch
bool hitPatch(solidParticleCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting a patch.
Definition: solidParticle.C:91
Foam::solidParticle::hitProcessorPatch
void hitProcessorPatch(solidParticleCloud &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
Definition: solidParticle.C:98