GradientDispersionRAS.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-2025 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 #include "standardNormal.H"
30 #include "meshTools.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  const dictionary& dict,
38  CloudType& owner
39 )
40 :
42  gradkPtr_(nullptr),
43  ownGradK_(false)
44 {}
45 
46 
47 template<class CloudType>
49 (
51 )
52 :
54  gradkPtr_(dm.gradkPtr_),
55  ownGradK_(dm.ownGradK_)
56 {
57  dm.ownGradK_ = false;
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
63 template<class CloudType>
65 {
66  cacheFields(false);
67 }
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
72 template<class CloudType>
74 {
76 
77  if (store)
78  {
79  gradkPtr_ = fvc::grad(*this->kPtr_).ptr();
80  ownGradK_ = true;
81  }
82  else
83  {
84  if (ownGradK_)
85  {
86  deleteDemandDrivenData(gradkPtr_);
87  gradkPtr_ = nullptr;
88  ownGradK_ = false;
89  }
90  }
91 }
92 
93 
94 template<class CloudType>
96 (
97  const scalar dt,
98  const label celli,
99  const vector& U,
100  const vector& Uc,
101  vector& UTurb,
102  scalar& tTurb
103 )
104 {
105  const polyMesh& mesh = this->owner().mesh();
106 
107  distributions::standardNormal& stdNormal = this->owner().stdNormal();
108 
109  const scalar cps = 0.16432;
110 
111  const scalar k = this->kPtr_->primitiveField()[celli];
112  const scalar epsilon =
113  this->epsilonPtr_->primitiveField()[celli] + rootVSmall;
114  const vector& gradk = this->gradkPtr_->primitiveField()[celli];
115 
116  const scalar UrelMag = mag(U - Uc - UTurb);
117 
118  const scalar tTurbLoc =
119  min(k/epsilon, cps*pow(k, 1.5)/epsilon/(UrelMag + small));
120 
121 
122  // Parcel is perturbed by the turbulence
123  if (dt < tTurbLoc)
124  {
125  tTurb += dt;
126 
127  if (tTurb > tTurbLoc)
128  {
129  tTurb = 0.0;
130 
131  const scalar sigma = sqrt(2.0*k/3.0);
132  const vector dir = -gradk/(mag(gradk) + small);
133 
134  scalar fac = 0.0;
135 
136  // In 2D calculations the -grad(k) is always
137  // away from the axis of symmetry
138  // This creates a 'hole' in the spray and to
139  // prevent this we let fac be both negative/positive
140  if (this->owner().mesh().nSolutionD() == 2)
141  {
142  fac = stdNormal.sample();
143  }
144  else
145  {
146  fac = mag(stdNormal.sample());
147  }
148 
149  UTurb = sigma*fac*dir;
150 
152  }
153  }
154  else
155  {
156  tTurb = great;
157  UTurb = Zero;
158  }
159 
160  return Uc + UTurb;
161 }
162 
163 
164 // ************************************************************************* //
label k
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
Base class for particle dispersion models based on RAS turbulence.
virtual void cacheFields(const bool store)
Cache carrier fields.
The velocity is perturbed in the direction of -grad(k), with a Gaussian random number distribution wi...
virtual ~GradientDispersionRAS()
Destructor.
virtual void cacheFields(const bool store)
Cache carrier fields.
GradientDispersionRAS(const dictionary &dict, CloudType &owner)
Construct from components.
virtual vector update(const scalar dt, const label celli, const vector &U, const vector &Uc, vector &UTurb, scalar &tTurb)
Update (disperse particles)
bool ownGradK_
Take ownership of the grad(k)
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Standard normal distribution. Not selectable.
virtual scalar sample() const
Sample the distribution.
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1029
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Template functions to aid in the implementation of demand driven data.
const scalar epsilon
Calculate the gradient of the given field.
U
Definition: pEqn.H:72
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:671
static const zero Zero
Definition: zero.H:97
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
void deleteDemandDrivenData(DataType *&dataPtr)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict