MaxwellianThermal.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 "MaxwellianThermal.H"
27 #include "constants.H"
28 
29 using namespace Foam::constant;
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  const dictionary& dict,
38 )
39 :
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  const polyMesh& mesh = this->owner().mesh();
60 
61  vector& U = p.U();
62 
63  scalar& Ei = p.Ei();
64 
65  label typeId = p.typeId();
66 
67  const label wppIndex = p.patch(mesh);
68 
69  const polyPatch& wpp = mesh.boundaryMesh()[wppIndex];
70 
71  label wppLocalFace = wpp.whichFace(p.face());
72 
73  const vector nw = p.normal(mesh);
74 
75  // Normal velocity magnitude
76  scalar U_dot_nw = U & nw;
77 
78  // Wall tangential velocity (flow direction)
79  vector Ut = U - U_dot_nw*nw;
80 
81  CloudType& cloud(this->owner());
82 
83  Random& rndGen(cloud.rndGen());
84 
85  while (mag(Ut) < small)
86  {
87  // If the incident velocity is parallel to the face normal, no
88  // tangential direction can be chosen. Add a perturbation to the
89  // incoming velocity and recalculate.
90 
91  U = vector
92  (
93  U.x()*(0.8 + 0.2*rndGen.scalar01()),
94  U.y()*(0.8 + 0.2*rndGen.scalar01()),
95  U.z()*(0.8 + 0.2*rndGen.scalar01())
96  );
97 
98  U_dot_nw = U & nw;
99 
100  Ut = U - U_dot_nw*nw;
101  }
102 
103  // Wall tangential unit vector
104  vector tw1 = Ut/mag(Ut);
105 
106  // Other tangential unit vector
107  vector tw2 = nw^tw1;
108 
109  scalar T = cloud.boundaryT().boundaryField()[wppIndex][wppLocalFace];
110 
111  scalar mass = cloud.constProps(typeId).mass();
112 
113  direction iDof = cloud.constProps(typeId).internalDegreesOfFreedom();
114 
115  U =
116  sqrt(physicoChemical::k.value()*T/mass)
117  *(
118  rndGen.scalarNormal()*tw1
119  + rndGen.scalarNormal()*tw2
120  - sqrt(-2.0*log(max(1 - rndGen.scalar01(), vSmall)))*nw
121  );
122 
123  U += cloud.boundaryU().boundaryField()[wppIndex][wppLocalFace];
124 
125  Ei = cloud.equipartitionInternalEnergy(T, iDof);
126 }
127 
128 
129 // ************************************************************************* //
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
MaxwellianThermal(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
virtual ~MaxwellianThermal()
Destructor.
virtual void correct(typename CloudType::parcelType &p)
Apply wall correction.
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
Collection of constants.
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))