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-2022 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,
35 )
36 :
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  const polyMesh& mesh = this->owner().mesh();
61 
62  vector& U = p.U();
63 
64  scalar& Ei = p.Ei();
65 
66  label typeId = p.typeId();
67 
68  const label wppIndex = p.patch(mesh);
69 
70  const polyPatch& wpp = mesh.boundaryMesh()[wppIndex];
71 
72  label wppLocalFace = wpp.whichFace(p.face());
73 
74  const vector nw = p.normal(mesh);
75 
76  // Normal velocity magnitude
77  scalar U_dot_nw = U & nw;
78 
79  CloudType& cloud(this->owner());
80 
81  Random& rndGen(cloud.rndGen());
82 
83  if (diffuseFraction_ > rndGen.scalar01())
84  {
85  // Diffuse reflection
86 
87  // Wall tangential velocity (flow direction)
88  vector Ut = U - U_dot_nw*nw;
89 
90  while (mag(Ut) < small)
91  {
92  // If the incident velocity is parallel to the face normal, no
93  // tangential direction can be chosen. Add a perturbation to the
94  // incoming velocity and recalculate.
95 
96  U = vector
97  (
98  U.x()*(0.8 + 0.2*rndGen.scalar01()),
99  U.y()*(0.8 + 0.2*rndGen.scalar01()),
100  U.z()*(0.8 + 0.2*rndGen.scalar01())
101  );
102 
103  U_dot_nw = U & nw;
104 
105  Ut = U - U_dot_nw*nw;
106  }
107 
108  // Wall tangential unit vector
109  vector tw1 = Ut/mag(Ut);
110 
111  // Other tangential unit vector
112  vector tw2 = nw^tw1;
113 
114  scalar T = cloud.boundaryT().boundaryField()[wppIndex][wppLocalFace];
115 
116  scalar mass = cloud.constProps(typeId).mass();
117 
118  direction iDof = cloud.constProps(typeId).internalDegreesOfFreedom();
119 
120  U =
121  sqrt(physicoChemical::k.value()*T/mass)
122  *(
123  rndGen.scalarNormal()*tw1
124  + rndGen.scalarNormal()*tw2
125  - sqrt(-2.0*log(max(1 - rndGen.scalar01(), vSmall)))*nw
126  );
127 
128  U += cloud.boundaryU().boundaryField()[wppIndex][wppLocalFace];
129 
130  Ei = cloud.equipartitionInternalEnergy(T, iDof);
131  }
132  else
133  {
134  // Specular reflection
135 
136  if (U_dot_nw > 0.0)
137  {
138  U -= 2.0*U_dot_nw*nw;
139  }
140  }
141 
142 }
143 
144 
145 // ************************************************************************* //
label k
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
virtual ~MixedDiffuseSpecular()
Destructor.
virtual void correct(typename CloudType::parcelType &p)
Apply wall correction.
MixedDiffuseSpecular(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Random number generator.
Definition: Random.H:58
scalar scalar01()
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:56
scalar scalarNormal()
Advance the state and return a scalar sample from a normal.
Definition: Random.C:31
Templated wall interaction model class.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:360
U
Definition: pEqn.H:72
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)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
uint8_t direction
Definition: direction.H:45
dictionary dict
volScalarField & p
Random rndGen(label(0))