MixedDiffuseSpecular.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 "MixedDiffuseSpecular.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 (
35  const dictionary& dict,
37 )
38 :
40  diffuseFraction_(this->coeffDict().getScalar("diffuseFraction"))
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45 
46 template<class CloudType>
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
53 template<class CloudType>
55 (
56  typename CloudType::parcelType& p
57 )
58 {
59  vector& U = p.U();
60 
61  scalar& Ei = p.Ei();
62 
63  label typeId = p.typeId();
64 
65  const label wppIndex = p.patch();
66 
67  const polyPatch& wpp = p.mesh().boundaryMesh()[wppIndex];
68 
69  label wppLocalFace = wpp.whichFace(p.face());
70 
71  const vector nw = p.normal();
72 
73  // Normal velocity magnitude
74  scalar U_dot_nw = U & nw;
75 
76  CloudType& cloud(this->owner());
77 
78  Random& rndGen = cloud.rndGen();
79 
80  if (diffuseFraction_ > rndGen.sample01<scalar>())
81  {
82  // Diffuse reflection
83 
84  // Wall tangential velocity (flow direction)
85  vector Ut = U - U_dot_nw*nw;
86 
87  while (mag(Ut) < SMALL)
88  {
89  // If the incident velocity is parallel to the face normal, no
90  // tangential direction can be chosen. Add a perturbation to the
91  // incoming velocity and recalculate.
92  U = vector
93  (
94  U.x()*(0.8 + 0.2*rndGen.sample01<scalar>()),
95  U.y()*(0.8 + 0.2*rndGen.sample01<scalar>()),
96  U.z()*(0.8 + 0.2*rndGen.sample01<scalar>())
97  );
98 
99  U_dot_nw = U & nw;
100 
101  Ut = U - U_dot_nw*nw;
102  }
103 
104  // Wall tangential unit vector
105  vector tw1 = Ut/mag(Ut);
106 
107  // Other tangential unit vector
108  vector tw2 = nw^tw1;
109 
110  scalar T = cloud.boundaryT().boundaryField()[wppIndex][wppLocalFace];
111 
112  scalar mass = cloud.constProps(typeId).mass();
113 
114  direction iDof = cloud.constProps(typeId).internalDegreesOfFreedom();
115 
116  U =
117  sqrt(physicoChemical::k.value()*T/mass)
118  *(
119  rndGen.GaussNormal<scalar>()*tw1
120  + rndGen.GaussNormal<scalar>()*tw2
121  - sqrt(-2.0*log(max(1 - rndGen.sample01<scalar>(), VSMALL)))*nw
122  );
123 
124  U += cloud.boundaryU().boundaryField()[wppIndex][wppLocalFace];
125 
126  Ei = cloud.equipartitionInternalEnergy(T, iDof);
127  }
128  else
129  {
130  // Specular reflection
131 
132  if (U_dot_nw > 0.0)
133  {
134  U -= 2.0*U_dot_nw*nw;
135  }
136  }
137 
138 }
139 
140 
141 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
MixedDiffuseSpecular.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::WallInteractionModel
Templated wall interaction model class.
Definition: DSMCCloud.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::MixedDiffuseSpecular::correct
virtual void correct(typename CloudType::parcelType &p)
Apply wall correction.
Definition: MixedDiffuseSpecular.C:55
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::MixedDiffuseSpecular::~MixedDiffuseSpecular
virtual ~MixedDiffuseSpecular()
Destructor.
Definition: MixedDiffuseSpecular.C:47
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
U
U
Definition: pEqn.H:72
Foam::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:448
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::direction
uint8_t direction
Definition: direction.H:52
rndGen
Random rndGen
Definition: createFields.H:23
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::MixedDiffuseSpecular::MixedDiffuseSpecular
MixedDiffuseSpecular(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Definition: MixedDiffuseSpecular.C:34