GradientDispersionRAS.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 | Copyright (C) 2011-2016 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 "GradientDispersionRAS.H"
27 #include "demandDrivenData.H"
28 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CloudType>
34 (
35  const dictionary& dict,
36  CloudType& owner
37 )
38 :
40  gradkPtr_(NULL),
41  ownGradK_(false)
42 {}
43 
44 
45 template<class CloudType>
47 (
49 )
50 :
52  gradkPtr_(dm.gradkPtr_),
53  ownGradK_(dm.ownGradK_)
54 {
55  dm.ownGradK_ = false;
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
61 template<class CloudType>
63 {
64  cacheFields(false);
65 }
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class CloudType>
72 {
74 
75  if (store)
76  {
77  gradkPtr_ = fvc::grad(*this->kPtr_).ptr();
78  ownGradK_ = true;
79  }
80  else
81  {
82  if (ownGradK_)
83  {
84  deleteDemandDrivenData(gradkPtr_);
85  gradkPtr_ = NULL;
86  ownGradK_ = false;
87  }
88  }
89 }
90 
91 
92 template<class CloudType>
94 (
95  const scalar dt,
96  const label celli,
97  const vector& U,
98  const vector& Uc,
99  vector& UTurb,
100  scalar& tTurb
101 )
102 {
103  cachedRandom& rnd = this->owner().rndGen();
104 
105  const scalar cps = 0.16432;
106 
107  const scalar k = this->kPtr_->primitiveField()[celli];
108  const scalar epsilon =
109  this->epsilonPtr_->primitiveField()[celli] + ROOTVSMALL;
110  const vector& gradk = this->gradkPtr_->primitiveField()[celli];
111 
112  const scalar UrelMag = mag(U - Uc - UTurb);
113 
114  const scalar tTurbLoc =
115  min(k/epsilon, cps*pow(k, 1.5)/epsilon/(UrelMag + SMALL));
116 
117 
118  // Parcel is perturbed by the turbulence
119  if (dt < tTurbLoc)
120  {
121  tTurb += dt;
122 
123  if (tTurb > tTurbLoc)
124  {
125  tTurb = 0.0;
126 
127  const scalar sigma = sqrt(2.0*k/3.0);
128  const vector dir = -gradk/(mag(gradk) + SMALL);
129 
130  scalar fac = 0.0;
131 
132  // In 2D calculations the -grad(k) is always
133  // away from the axis of symmetry
134  // This creates a 'hole' in the spray and to
135  // prevent this we let fac be both negative/positive
136  if (this->owner().mesh().nSolutionD() == 2)
137  {
138  fac = rnd.GaussNormal<scalar>();
139  }
140  else
141  {
142  fac = mag(rnd.GaussNormal<scalar>());
143  }
144 
145  UTurb = sigma*fac*dir;
146  }
147  }
148  else
149  {
150  tTurb = GREAT;
151  UTurb = Zero;
152  }
153 
154  return Uc + UTurb;
155 }
156 
157 
158 // ************************************************************************* //
The velocity is perturbed in the direction of -grad(k), with a Gaussian random number distribution wi...
virtual ~GradientDispersionRAS()
Destructor.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
bool ownGradK_
Take ownership of the grad(k)
Random number generator.
Definition: cachedRandom.H:63
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:812
label k
Boltzmann constant.
Type GaussNormal()
Return a sample whose components are normally distributed.
Calculate the gradient of the given field.
Random & rndGen()
Return refernce to the random object.
Definition: DSMCCloudI.H:121
static const zero Zero
Definition: zero.H:91
const volVectorField * gradkPtr_
Gradient of k.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Template functions to aid in the implementation of demand driven data.
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
scalar epsilon
GradientDispersionRAS(const dictionary &dict, CloudType &owner)
Construct from components.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual vector update(const scalar dt, const label celli, const vector &U, const vector &Uc, vector &UTurb, scalar &tTurb)
Update (disperse particles)
Base class for particle dispersion models based on RAS turbulence.
void deleteDemandDrivenData(DataPtr &dataPtr)
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
virtual void cacheFields(const bool store)
Cache carrier fields.