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-2024 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 #include "standardNormal.H"
29 
30 using namespace Foam::constant;
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  const dictionary& dict,
39 )
40 :
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  // Wall tangential velocity (flow direction)
80  vector Ut = U - U_dot_nw*nw;
81 
82  CloudType& cloud(this->owner());
83 
84  randomGenerator& rndGen = cloud.rndGen();
85  distributions::standardNormal& stdNormal = cloud.stdNormal();
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 
93  U = vector
94  (
95  U.x()*(0.8 + 0.2*rndGen.scalar01()),
96  U.y()*(0.8 + 0.2*rndGen.scalar01()),
97  U.z()*(0.8 + 0.2*rndGen.scalar01())
98  );
99 
100  U_dot_nw = U & nw;
101 
102  Ut = U - U_dot_nw*nw;
103  }
104 
105  // Wall tangential unit vector
106  vector tw1 = Ut/mag(Ut);
107 
108  // Other tangential unit vector
109  vector tw2 = nw^tw1;
110 
111  scalar T = cloud.boundaryT().boundaryField()[wppIndex][wppLocalFace];
112 
113  scalar mass = cloud.constProps(typeId).mass();
114 
115  direction iDof = cloud.constProps(typeId).internalDegreesOfFreedom();
116 
117  U =
118  sqrt(physicoChemical::k.value()*T/mass)
119  *(
120  stdNormal.sample()*tw1
121  + stdNormal.sample()*tw2
122  - sqrt(-2.0*log(max(1 - rndGen.scalar01(), vSmall)))*nw
123  );
124 
125  U += cloud.boundaryU().boundaryField()[wppIndex][wppLocalFace];
126 
127  Ei = cloud.equipartitionInternalEnergy(T, iDof);
128 }
129 
130 
131 // ************************************************************************* //
label k
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
MaxwellianThermal(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
virtual ~MaxwellianThermal()
Destructor.
virtual void correct(typename CloudType::parcelType &p)
Apply wall correction.
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:162
Standard normal distribution. Not selectable.
virtual scalar sample() const
Sample the distribution.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
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
Random number generator.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
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
randomGenerator rndGen(653213)
volScalarField & p