MixedDiffuseSpecular.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 
26 #include "MixedDiffuseSpecular.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class CloudType>
32 (
33  const dictionary& dict,
34  CloudType& cloud
35 )
36 :
37  WallInteractionModel<CloudType>(dict, cloud, typeName),
38  diffuseFraction_
39  (
40  this->coeffDict().template lookup<scalar>("diffuseFraction")
41  )
42 {}
43 
44 
45 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
46 
47 template<class CloudType>
49 {}
50 
51 
52 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
53 
54 template<class CloudType>
56 (
57  typename CloudType::parcelType& p
58 )
59 {
60  vector& U = p.U();
61 
62  scalar& Ei = p.Ei();
63 
64  label typeId = p.typeId();
65 
66  const label wppIndex = p.patch();
67 
68  const polyPatch& wpp = p.mesh().boundaryMesh()[wppIndex];
69 
70  label wppLocalFace = wpp.whichFace(p.face());
71 
72  const vector nw = p.normal();
73 
74  // Normal velocity magnitude
75  scalar U_dot_nw = U & nw;
76 
77  CloudType& cloud(this->owner());
78 
79  Random& rndGen(cloud.rndGen());
80 
81  if (diffuseFraction_ > rndGen.scalar01())
82  {
83  // Diffuse reflection
84 
85  // Wall tangential velocity (flow direction)
86  vector Ut = U - U_dot_nw*nw;
87 
88  while (mag(Ut) < small)
89  {
90  // If the incident velocity is parallel to the face normal, no
91  // tangential direction can be chosen. Add a perturbation to the
92  // incoming velocity and recalculate.
93 
94  U = vector
95  (
96  U.x()*(0.8 + 0.2*rndGen.scalar01()),
97  U.y()*(0.8 + 0.2*rndGen.scalar01()),
98  U.z()*(0.8 + 0.2*rndGen.scalar01())
99  );
100 
101  U_dot_nw = U & nw;
102 
103  Ut = U - U_dot_nw*nw;
104  }
105 
106  // Wall tangential unit vector
107  vector tw1 = Ut/mag(Ut);
108 
109  // Other tangential unit vector
110  vector tw2 = nw^tw1;
111 
112  scalar T = cloud.boundaryT().boundaryField()[wppIndex][wppLocalFace];
113 
114  scalar mass = cloud.constProps(typeId).mass();
115 
116  direction iDof = cloud.constProps(typeId).internalDegreesOfFreedom();
117 
118  U =
119  sqrt(physicoChemical::k.value()*T/mass)
120  *(
121  rndGen.scalarNormal()*tw1
122  + rndGen.scalarNormal()*tw2
123  - sqrt(-2.0*log(max(1 - rndGen.scalar01(), vSmall)))*nw
124  );
125 
126  U += cloud.boundaryU().boundaryField()[wppIndex][wppLocalFace];
127 
128  Ei = cloud.equipartitionInternalEnergy(T, iDof);
129  }
130  else
131  {
132  // Specular reflection
133 
134  if (U_dot_nw > 0.0)
135  {
136  U -= 2.0*U_dot_nw*nw;
137  }
138  }
139 
140 }
141 
142 
143 // ************************************************************************* //
const volScalarField & boundaryT() const
Return macroscopic temperature.
Definition: DSMCCloudI.H:193
dictionary dict
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
dimensionedScalar log(const dimensionedScalar &ds)
virtual void correct(typename CloudType::parcelType &p)
Apply wall correction.
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 > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
uint8_t direction
Definition: direction.H:45
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Cmpt & z() const
Definition: VectorI.H:87
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label k
Boltzmann constant.
const Cmpt & y() const
Definition: VectorI.H:81
Random rndGen(label(0))
virtual ~MixedDiffuseSpecular()
Destructor.
const volVectorField & boundaryU() const
Return macroscopic velocity.
Definition: DSMCCloudI.H:201
Random & rndGen()
Return references to the random object.
Definition: DSMCCloudI.H:121
const Cmpt & x() const
Definition: VectorI.H:75
Random number generator.
Definition: Random.H:57
scalar equipartitionInternalEnergy(scalar temperature, direction internalDegreesOfFreedom)
Generate a random internal energy, sampled from the.
Definition: DSMCCloud.C:1029
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:215
const volScalarField & T
U
Definition: pEqn.H:72
dimensioned< scalar > mag(const dimensioned< Type > &)
Templated wall interaction model class.
Definition: DSMCCloud.H:59
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:95
label nw
Definition: createFields.H:25
MixedDiffuseSpecular(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:383