IntegralScaleBox.H
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) 2022 OpenCFD Ltd.
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 Class
27  Foam::turbulence::IntegralScaleBox
28 
29 Description
30  Class to describe the integral-scale container being used in the
31  \c turbulentDigitalFilterInletFvPatchField boundary condition.
32 
33 SourceFiles
34  IntegralScaleBox.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef Foam_IntegralScaleBox_H
39 #define Foam_IntegralScaleBox_H
40 
41 #include "Random.H"
42 #include "primitivePatch.H"
43 #include "boundBox.H"
44 #include "coordinateSystem.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 namespace turbulence
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class IntegralScaleBox Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 template<class Type>
58 class IntegralScaleBox
59 {
60  // Private Typedefs
61 
62  //- Type of input set of integral scales
63  typedef typename std::conditional
64  <
65  std::is_same<Type, scalar>::value,
66  vector,
67  tensor
68  >::type TypeL;
69 
70 
71  // Private Enumerations
72 
73  //- Options for the correlation function form
74  enum kernelType : bool
75  {
76  GAUSSIAN = true,
77  EXPONENTIAL = false
78  };
79 
80  //- Names for kernelType
81  static const Enum<kernelType> kernelTypeNames;
82 
83 
84  // Private Data
85 
86  //- Const reference to patch
87  const fvPatch& p_;
88 
89  //- Virtual patch of the integral-scale box
90  autoPtr<primitivePatch> patchPtr_;
91 
92  //- Local coordinate system
94 
95  //- Correlation function form of the kernel
96  const enum kernelType kernelType_;
97 
98  //- Random number generator
99  Random rndGen_;
100 
101  //- Number of faces on integral-scale box inlet plane (<e2> <e3>)
102  const Vector2D<label> n_;
103 
104  //- Uniform cell size on integral-scale box inlet plane [m]
105  Vector2D<scalar> delta_;
106 
107  //- Span of patch bounding box in local coordinate system
108  vector boundingBoxSpan_;
109 
110  //- Min point of patch bounding box in local coordinate system
111  vector boundingBoxMin_;
112 
113  //- Integral scale set in local coordinates
114  const TypeL L_;
115 
116  //- Random-number box ranges to be filtered by corresponding kernels
117  // e.g. for U: (Lxu, Lxv, Lxw, Lyu, Lyv, Lyw, Lzu, Lzv, Lzw)
118  labelList spans_;
119 
120  //- Random-number box
121  // First List: Components of physical variable, e.g. (u,v,w)
122  // Second List: Containing all random numbers in a continuous list
123  // for each component of Cartesian coordinate system, e.g. (x,y,z)
124  // e.g. for U, box_[0] contains sets for u-component in x-y-z order
125  scalarListList box_;
126 
127  //- Kernel coefficients corresponding to L [-]
128  scalarListList kernel_;
129 
130  //- Vertices of virtual patch of the integral-scale box
131  pointField patchPoints_;
132 
133  //- Faces of virtual patch of the integral-scale box
134  faceList patchFaces_;
135 
136  // Forward-stepwise method information
137 
138  //- Flag to enable the forward-stepwise method
139  const bool fsm_;
140 
141  //- One of the two exponential functions in FSM
142  Type C1_;
143 
144  //- One of the two exponential functions in FSM
145  Type C2_;
146 
147  //- Kernel-applied previous-time-step slice field in FSM
148  Field<Type> slice_;
149 
150 
151  // Private Member Functions
152 
153  //- Initialise the local coordinate system with dictionary
154  autoPtr<coordinateSystem> calcCoordinateSystem
155  (
156  const dictionary& dict
157  ) const;
158 
159  //- Initialise the local coordinate system with patch properties
160  void calcCoordinateSystem();
161 
162  //- Initialise the bounding box in local coordinates
163  Vector2D<vector> calcBoundBox() const;
164 
165  //- Initialise the object: delta
166  Vector2D<scalar> calcDelta() const;
167 
168  //- Initialise the object: spans
169  labelList calcSpans() const;
170 
171  //- Return the kernel once per simulation
172  scalarListList calcKernel() const;
173 
174  //- Initialise the object: box
175  scalarListList calcBox();
176 
177  //- Initialise the object: patchPoints
178  pointField calcPatchPoints() const;
179 
180  //- Initialise the object: patchFaces
181  faceList calcPatchFaces() const;
182 
183  //- Initialise the object: virtual patch
184  void calcPatch();
185 
186  //- Convert integral scales in seconds/meters
187  //- to time-step units/integral-scale box inlet cell-size units
188  TypeL convert(const TypeL& L) const;
189 
190  // Forward-stepwise method
191 
192  //- Initialise the object: C1 for scalar fields
193  scalar calcC1(const vector& L) const;
194 
195  //- Initialise the object: C1 for vector fields
196  vector calcC1(const tensor& L) const;
197 
198  //- Initialise the object: C2 for scalar fields
199  scalar calcC2(const vector& L) const;
200 
201  //- Initialise the object: C2 for vector fields
202  vector calcC2(const tensor& L) const;
203 
204  //- Update C1 and C2 coefficients
205  void updateC1C2();
206 
207 
208  //- No copy assignment
209  void operator=(const IntegralScaleBox&) = delete;
210 
211 
212 public:
213 
214  // Constructors
215 
216  //- Construct from patch
217  explicit IntegralScaleBox(const fvPatch& p);
218 
219  //- Copy construct with patch
220  IntegralScaleBox(const fvPatch& p, const IntegralScaleBox& b);
221 
222  //- Construct from patch and dictionary
223  IntegralScaleBox(const fvPatch& p, const dictionary& dict);
224 
225  //- Copy construct
227 
228 
229  // Public Data
230 
231  //- Flag to activate debug statements
232  static int debug;
233 
234 
235  // Member Functions
236 
237  // Access
238 
239  //- Return const reference to integral-scale box inlet patch
240  const primitivePatch& patch()
241  {
242  if (!patchPtr_)
243  {
244  calcPatch();
245  }
246 
247  return *patchPtr_;
248  }
249 
250  //- Return the object: fsm
251  bool fsm() const noexcept
252  {
253  return fsm_;
254  }
255 
256 
257  // Evaluation
258 
259  //- Initialise integral-scale box properties
260  // Avoid constructor level evaluations
261  void initialise();
262 
263  //- Discard current time-step integral-scale box slice (the closest
264  //- to the patch) by shifting from the back to the front
265  void shift();
266 
267  //- Add a new integral-scale box slice to the rear of the box
268  void refill();
269 
270  //- Embed two-point correlations, i.e. L
271  // Apply three-dimensional "valid"-type separable
272  // convolution summation algorithm
273  Field<Type> convolve() const;
274 
275  // Forward-stepwise method
276 
277  //- Apply forward-stepwise correlation for scalar fields
278  void correlate(scalarField& fld);
279 
280  //- Apply forward-stepwise correlation for vector fields
281  void correlate(vectorField& fld);
282 
283 
284  // I-O
285 
286  //- Write integral-scale box settings
287  void write(Ostream&) const;
288 };
289 
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 } // End namespace turbulence
294 } // End namespace Foam
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 #ifdef NoRepository
299  #include "IntegralScaleBox.C"
300 #endif
301 
302 #endif
303 
304 // ************************************************************************* //
dictionary dict
static int debug
Flag to activate debug statements.
Field< Type > convolve() const
Embed two-point correlations, i.e. L.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Class to describe the integral-scale container being used in the turbulentDigitalFilterInletFvPatchFi...
bool fsm() const noexcept
Return the object: fsm.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
A list of faces which address into the list of points.
void refill()
Add a new integral-scale box slice to the rear of the box.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
void shift()
Discard current time-step integral-scale box slice (the closest to the patch) by shifting from the ba...
void initialise()
Initialise integral-scale box properties.
Vector< scalar > vector
Definition: vector.H:57
Random number generator.
Definition: Random.H:55
void correlate(scalarField &fld)
Apply forward-stepwise correlation for scalar fields.
void write(Ostream &) const
Write integral-scale box settings.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
volScalarField & p
IntegralScaleBox(const fvPatch &p)
Construct from patch.
Tensor of scalars, i.e. Tensor<scalar>.
const primitivePatch & patch()
Return const reference to integral-scale box inlet patch.
Namespace for OpenFOAM.