MaxwellianThermal.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 "MaxwellianThermal.H"
29 #include "constants.H"
30 
31 using namespace Foam::constant;
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class CloudType>
37 (
38  const dictionary& dict,
40 )
41 :
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
47 
48 template<class CloudType>
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 template<class CloudType>
57 (
58  typename CloudType::parcelType& p
59 )
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();
68 
69  const polyPatch& wpp = p.mesh().boundaryMesh()[wppIndex];
70 
71  label wppLocalFace = wpp.whichFace(p.face());
72 
73  const vector nw = p.normal();
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  U = vector
91  (
92  U.x()*(0.8 + 0.2*rndGen.sample01<scalar>()),
93  U.y()*(0.8 + 0.2*rndGen.sample01<scalar>()),
94  U.z()*(0.8 + 0.2*rndGen.sample01<scalar>())
95  );
96 
97  U_dot_nw = U & nw;
98 
99  Ut = U - U_dot_nw*nw;
100  }
101 
102  // Wall tangential unit vector
103  vector tw1 = Ut/mag(Ut);
104 
105  // Other tangential unit vector
106  vector tw2 = nw^tw1;
107 
108  scalar T = cloud.boundaryT().boundaryField()[wppIndex][wppLocalFace];
109 
110  scalar mass = cloud.constProps(typeId).mass();
111 
112  direction iDof = cloud.constProps(typeId).internalDegreesOfFreedom();
113 
114  U =
115  sqrt(physicoChemical::k.value()*T/mass)
116  *(
117  rndGen.GaussNormal<scalar>()*tw1
118  + rndGen.GaussNormal<scalar>()*tw2
119  - sqrt(-2.0*log(max(1 - rndGen.sample01<scalar>(), VSMALL)))*nw
120  );
121 
122  U += cloud.boundaryU().boundaryField()[wppIndex][wppLocalFace];
123 
124  Ei = cloud.equipartitionInternalEnergy(T, iDof);
125 }
126 
127 
128 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::constant::physicoChemical::k
const dimensionedScalar k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::MaxwellianThermal::~MaxwellianThermal
virtual ~MaxwellianThermal()
Destructor.
Definition: MaxwellianThermal.C:49
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
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::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
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
constants.H
Foam::MaxwellianThermal::correct
virtual void correct(typename CloudType::parcelType &p)
Apply wall correction.
Definition: MaxwellianThermal.C:57
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
U
U
Definition: pEqn.H:72
MaxwellianThermal.H
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)
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::MaxwellianThermal::MaxwellianThermal
MaxwellianThermal(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Definition: MaxwellianThermal.C:37