turbulentDigitalFilterInletFvPatchVectorField.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) 2019-2020 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::turbulentDigitalFilterFvPatchVectorField
28 
29 Group
30  grpInletBoundaryConditions
31 
32 Description
33  Digital-filter based boundary condition for velocity, i.e. \c U, to generate
34  synthetic turbulence-alike time-series for LES and DES turbulent flow
35  computations from input turbulence statistics.
36 
37  References:
38  \verbatim
39  Digital-filter method-based generator (DFM) (tag:KSJ):
40  Klein, M., Sadiki, A., & Janicka, J. (2003).
41  A digital filter based generation of inflow data for spatially
42  developing direct numerical or large eddy simulations.
43  Journal of computational Physics, 186(2), 652-665.
44  DOI:10.1016/S0021-9991(03)00090-1
45 
46  Forward-stepwise method-based generator (FSM) (tag:XC)
47  Xie, Z. T., & Castro, I. P. (2008).
48  Efficient generation of inflow conditions for
49  large eddy simulation of street-scale flows.
50  Flow, turbulence and combustion, 81(3), 449-470.
51  DOI:10.1007/s10494-008-9151-5
52 
53  Mass-inflow rate correction (tag:KCX):
54  Kim, Y., Castro, I. P., & Xie, Z. T. (2013).
55  Divergence-free turbulence inflow conditions for
56  large-eddy simulations with incompressible flow solvers.
57  Computers & Fluids, 84, 56-68.
58  DOI:10.1016/j.compfluid.2013.06.001
59  \endverbatim
60 
61  In \c DFM or \c FSM, a random number set (mostly white noise), and a group
62  of target statistics (mostly mean flow, Reynolds stress tensor profiles and
63  length-scale sets) are merged into a new number set (stochastic time-series,
64  yet consisting of the statistics) by a chain of mathematical operations
65  whose characteristics are designated by the target statistics, so that the
66  realised statistics of the new sets could match the target.
67 
68  \verbatim
69  Random number sets ---->-|
70  |
71  DFM or FSM ---> New stochastic time-series consisting
72  | turbulence statistics
73  Turbulence statistics ->-|
74  \endverbatim
75 
76  The main difference between \c DFM and \c FSM is that \c FSM replaces
77  the expensive-to-run streamwise convolution summation in \c DFM by a simpler
78  and an almost-equivalent-in-effect numerical procedure in order to reduce
79  computational costs. Accordingly, \c FSM potentially brings computational
80  resource advantages for computations involving relatively large streamwise
81  length-scale sets and small time-steps.
82 
83  Synthetic turbulence is generated on a virtual rectangular structured-mesh
84  plane, which is parallel to the chosen patch, and is mapped onto this patch
85  by the selected mapping method.
86 
87 Usage
88  Example of the boundary condition specification:
89  \verbatim
90  <patchName>
91  {
92  // Mandatory entries (unmodifiable)
93  type turbulentDigitalFilterInlet;
94  n (<nHeight> <nWidth>);
95  L (<L1> <L2> ... <L9>);
96  R uniform (<Rxx> <Rxy> <Rxz> <Ryy> <Ryz> <Rzz>);
97  UMean uniform (1 0 0);
98  Ubulk 10.0;
99 
100  // Optional entries (unmodifiable)
101  fsm false;
102  Gaussian true; // always false for FSM
103  fixSeed true;
104  continuous false;
105  correctFlowRate true;
106  mapMethod nearestCell;
107  perturb 1e-5;
108  C1 -1.5707; //-0.5*PI;
109  C1FSM -0.7854 //-0.25*PI;
110  C2FSM -1.5707; //-0.5*PI;
111 
112  // Optional (inherited) entries
113  ...
114  }
115  \endverbatim
116 
117  where the entries mean:
118  \table
119  Property | Description | Type | Req'd | Dflt
120  type | Type name: turbulentDigitalFilterInlet | word | yes | -
121  n | Number of cells on turbulence generation plane <!--
122  --> | tuple of labels | yes | -
123  L | Integral length-scale set <!--
124  --> (Lxu Lxv Lxw Lyu Lyv Lyw Lzu Lzv Lzw) [m] | tensor | yes | -
125  R | Reynolds stress tensor set <!--
126  --> (xx xy xz yy yz zz) [m2/s2] | symmTensorField | yes | -
127  UMean | Mean velocity profile [m/s] | vectorField | yes | -
128  Ubulk | Characteristic patch-normal bulk flow speed [m/s] <!--
129  --> | scalar | yes | -
130  fsm | Flag to turn on the forward-stepwise method | bool | no | false
131  Gaussian | Autocorrelation function form | bool | no | true
132  fixSeed | Flag to fix random-number generator seed to 1234 <!--
133  --> or generate a new seed based on clock-time per simulation <!--
134  --> | bool | no | true
135  continuous | Flag to write random-number sets at output time, <!--
136  --> and to read them on restart. Otherwise, generate <!--
137  --> new random-number sets of restart | bool | no false
138  correctFlowRate | Flag to correct mass-inflow rate on turbulence <!--
139  --> plane in (only) streamwise direction | bool | no | true
140  mapMethod | Interpolation-to-patch method | word | no | nearestCell
141  perturb | Point perturbation for planarInterpolation mapMethod <!--
142  --> | scalar | no | 1e-5
143  C1 | Model constant shaping autocorrelation function <!--
144  --> (KSJ:Eq. 14) | scalar | no | -0.5*PI
145  C1FSM | Model coefficient in FSM (XC:Eq. 14) | scalar | no | -0.25*PI
146  C2FSM | Model coefficient in FSM (XC:Eq. 14) | scalar | no | -0.5*PI
147  \endtable
148 
149  The inherited entries are elaborated in:
150  - \link fixedValueFvPatchFields.H \endlink
151 
152  Options for the \c fsm entry:
153  \verbatim
154  false | Method due to (KSJ)
155  true | Method due to (XC)
156  \endverbatim
157 
158  Options for the \c Gaussian entry:
159  \verbatim
160  true | Gaussian function
161  false | Exponential function (only option for FSM)
162  \endverbatim
163 
164  Options for the \c mapMethod entry:
165  \verbatim
166  nearestCell | One-to-one direct map, no interpolation
167  planarInterpolation | Bilinear interpolation
168  \endverbatim
169 
170  Patch-profile input is available for two entries:
171  \verbatim
172  R | Reynolds stress tensor
173  UMean | Mean velocity
174  \endverbatim
175  where the input profiles and profile coordinates are located in:
176  \verbatim
177  Coordinates | $FOAM_CASE/constant/boundaryData/<patchName>/points
178  R/UMean | $FOAM_CASE/constant/boundaryData/<patchName>/0/\{R/UMean\}
179  \endverbatim
180 
181  \c points file contains a list of three-dimensional coordinates, and
182  profile data files provide a value corresponding to each coordinate.
183 
184 Note
185  - \c mapMethod=planarInterpolation option needs point coordinates that can
186  form a plane.
187  - \c adjustTimeStep=true option is currently not fully supported.
188  - In order to obtain Reynolds stress tensor information, experiments, RANS
189  simulations or engineering relations can be used.
190  - \c continuous=true means deterministic-statistical consistent restart
191  (relatively more expensive), and \c continuous=false means deterministic
192  discontinuity in synthetic turbulence time-series by keeping statistical
193  consistency (relatively cheaper).
194  - For \c L, the first three entries should always correspond to the
195  length scales in association with the convective (streamwise) mean flow
196  direction.
197  - Streamwise integral length scales are converted to integral time scales
198  by using Taylor's frozen turbulence hypothesis, and \c Ubulk.
199 
200 See also
201  - turbulentDFSEMInletFvPatchVectorField.C
202 
203 SourceFiles
204  turbulentDigitalFilterInletFvPatchVectorField.C
205  turbulentDigitalFilterInletFvPatchVectorFieldTemplates.C
206 
207 \*---------------------------------------------------------------------------*/
208 
209 #ifndef turbulentDigitalFilterInletFvPatchVectorField_H
210 #define turbulentDigitalFilterInletFvPatchVectorField_H
211 
212 #include "fixedValueFvPatchFields.H"
213 #include "Random.H"
214 #include "fieldTypes.H"
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 namespace Foam
219 {
220 
221 class pointToPointPlanarInterpolation;
222 
223 /*---------------------------------------------------------------------------*\
224  Class turbulentDigitalFilterInletFvPatchVectorField Declaration
225 \*---------------------------------------------------------------------------*/
226 
227 class turbulentDigitalFilterInletFvPatchVectorField
228 :
229  public fixedValueFvPatchVectorField
230 {
231  // Private Data
232 
233  //- Bilinear interpolation (for 'mapMethod=planarInterpolation')
234  mutable autoPtr<pointToPointPlanarInterpolation> mapperPtr_;
235 
236  //- Flag to enable the forward-stepwise method
237  const bool fsm_;
238 
239  //- Flag to select correlation function form: Gaussian or exponential
240  const bool Gaussian_;
241 
242  //- Flag to fix the random-number generator seed to 1234 or
243  //- generate a new seed based on clock-time per simulation
244  const bool fixSeed_;
245 
246  //- Flag to write random-number sets at output time, and to read them
247  //- on restart. Otherwise, generate new random-number sets on restart
248  const bool continuous_;
249 
250  //- Flag to correct mass flow rate on turbulence plane
251  const bool correctFlowRate_;
252 
253  //- Internal flag to read R from data files
254  bool interpR_;
255 
256  //- Internal flag to read UMean from data files
257  bool interpUMean_;
258 
259  //- Method for interpolation between a patch and turbulence plane
260  const word mapMethod_;
261 
262  //- Current time index
263  label curTimeIndex_;
264 
265  //- Characteristic patch-normal bulk flow speed [m/s]
266  const scalar Ubulk_;
267 
268  //- Model constant shaping autocorrelation function (KSJ:Eq. 14)
269  const scalar C1_;
270 
271  //- Fraction of perturbation (fraction of bounding box) to add
272  const scalar perturb_;
273 
274  //- First time-step mass/volumetric flow rate
275  scalar flowRate_;
276 
277  //- Random number generator
278  Random rndGen_;
279 
280  //- Number of cells on turbulence plane (<nHeight> <nWidth>) [-]
281  const Tuple2<label, label> n_;
282 
283  //- Uniform mesh size on turbulence plane (reversed) [1/m]
284  Vector2D<scalar> invDelta_;
285 
286  //- Index pairs between patch and turbulence plane
287  //- for the nearest cell mapping
288  const List<Pair<label>> indexPairs_;
289 
290  //- Reynolds stress tensor profile field in global coordinates [m2/s2]
291  const symmTensorField R_;
292 
293  //- Lund-Wu-Squires transformed R field (Cholesky decomp.) [m/s]
294  //- Mapped onto actual mesh patch rather than turbulence plane
295  // (KSJ:Eq. 5)
296  const symmTensorField Lund_;
297 
298  //- Mean inlet velocity profile field in global coordinates [m/s]
299  vectorField UMean_;
300 
301  //- Integral length-scale set per turbulence plane section in local
302  //- coordinates (e1u, e1v, e1w, e2u, e2v, e2w, e3u, e3v, e3w) [m]
303  // Backup of L_ for restart purposes
304  const tensor Lbak_;
305 
306  //- Integral length-scale set in mesh units [cell]
307  const tensor L_;
308 
309  //- One of the two model coefficients in FSM
310  // (XC:The argument of the first exp func in Eq. 14)
311  const scalar C1FSM_;
312 
313  //- One of the two model coefficients in FSM
314  // (XC:The argument of the second exp func in Eq. 14)
315  const scalar C2FSM_;
316 
317  //- One of the two exponential functions in FSM
318  // (XC:The first exponential function in Eq. 14)
319  const List<scalar> coeffs1FSM_;
320 
321  //- One of the two exponential functions in FSM
322  // (XC:The first exponential function in Eq. 14)
323  const List<scalar> coeffs2FSM_;
324 
325  //- Number of cells in random-number box [cell]
326  //- Random-number sets within box are filtered with filterCoeffs_
327  const List<label> szBox_;
328 
329  //- Convenience factors for two-dimensional random-number box [-]
330  const List<label> boxFactors2D_;
331 
332  //- Convenience factors for three-dimensional random-number box [-]
333  const List<label> boxFactors3D_;
334 
335  //- Index to the first elem of last plane of random-number box [-]
336  const List<label> iNextToLastPlane_;
337 
338  //- Random-number sets distributed over three-dimensional box (u, v, w)
339  List<List<scalar>> box_;
340 
341  //- Filter coefficients corresponding to L [-]
342  const List<List<scalar>> filterCoeffs_;
343 
344  //- Filter-applied random-number sets [m/s], i.e. turbulence plane
345  List<List<scalar>> filteredBox_;
346 
347  //- Filter-applied previous-time-step velocity field [m/s] used in FSM
348  vectorField U0_;
349 
350 
351  // Private Member Functions
352 
353  //- Return a reference to the patch mapper object
354  const pointToPointPlanarInterpolation& patchMapper() const;
355 
356  //- Helper function to interpolate values from the boundary data or
357  //- read from dictionary
358  template<class Type>
359  tmp<Field<Type>> interpolateOrRead
360  (
361  const word& fieldName,
362  const dictionary& dict,
363  bool& interpolateField
364  ) const;
365 
366  //- Helper function to interpolate values from the boundary data
367  template<class Type>
368  tmp<Field<Type>> interpolateBoundaryData
369  (
370  const word& fieldName
371  ) const;
372 
373  //- Generate random-number sets obeying the standard normal distribution
374  template<class Form, class Type>
375  Form randomSet(const label len);
376 
377  //- Compute nearest cell index-pairs between turbulence plane and patch
378  List<Pair<label>> indexPairs();
379 
380  //- Check R on patch for mathematical domain errors
381  void checkR() const;
382 
383  //- Compute Lund-Wu-Squires transformation
384  symmTensorField calcLund() const;
385 
386  //- Compute the first time-step mass/
387  //- volumetric flow rate based on UMean
388  scalar calcFlowRate() const;
389 
390  //- Convert integral length scales in meters
391  //- to turbulence plane cell-size units
392  tensor meterToCell(const tensor& L) const;
393 
394  //- Resource allocation functions for the convenience factors
395  List<label> initBox() const;
396  List<label> initFactors2D() const;
397  List<label> initFactors3D() const;
398  List<label> initPlaneFactors() const;
399 
400  //- Compute various convenience factors for random-number box
401  List<List<scalar>> fillBox();
402 
403  //- Compute filter coeffs once per simulation
404  List<List<scalar>> calcFilterCoeffs() const;
405 
406  //- Discard current time-step random-box plane (closest to patch) by
407  //- shifting from the back to the front, and add new plane to the back
408  void shiftRefill();
409 
410  //- Map two-point correlated random-number sets
411  //- on patch based on chosen mapping method
412  void mapFilteredBox(vectorField& U);
413 
414  //- Embed one-point correlations, i.e. R, on patch
415  void onePointCorrs(vectorField& U) const;
416 
417  //- Embed two-point correlations, i.e. L, on box
418  // Three-dimensional "valid"-type separable
419  // convolution summation algorithm
420  // (Based on Song Ho Ahn's two-dimensional "full"-type convolution)
421  void twoPointCorrs();
422 
423  //- Compute coeffs1FSM_ once per simulation
424  List<scalar> calcCoeffs1FSM() const;
425 
426  //- Compute coeffs2FSM_ once per simulation
427  List<scalar> calcCoeffs2FSM() const;
428 
429 
430 public:
431 
432  //- Runtime type information
433  TypeName("turbulentDigitalFilterInlet");
434 
435 
436  // Constructors
437 
438  //- Construct from patch and internal field
440  (
441  const fvPatch&,
443  );
444 
445  //- Construct from patch, internal field and dictionary
447  (
448  const fvPatch&,
450  const dictionary&
451  );
452 
453  //- Construct by mapping given
454  //- turbulentDigitalFilterInletFvPatchVectorField onto a new patch
456  (
458  const fvPatch&,
460  const fvPatchFieldMapper&
461  );
462 
463  //- Construct as copy
465  (
467  );
468 
469  //- Construct and return a clone
470  virtual tmp<fvPatchVectorField> clone() const
471  {
473  (
475  );
476  }
477 
478  //- Construct as copy setting internal field reference
480  (
483  );
484 
485  //- Construct and return a clone setting internal field reference
487  (
489  ) const
490  {
492  (
494  );
495  }
496 
497 
498  //- Destructor
500 
501 
502  // Member Functions
503 
504  // Evaluation functions
505 
506  //- Update the coefficients associated with the patch field
507  virtual void updateCoeffs();
508 
509 
510  // Mapping functions
511 
512  //- Map (and resize as needed) from self given a mapping object
513  virtual void autoMap(const fvPatchFieldMapper& m);
514 
515  //- Reverse map the given fvPatchField onto this fvPatchField
516  virtual void rmap
517  (
518  const fvPatchVectorField& ptf,
519  const labelList& addr
520  );
521 
522 
523  //- Write
524  virtual void write(Ostream&) const;
525 };
526 
527 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
528 
529 } // End namespace Foam
530 
531 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
532 
533 #ifdef NoRepository
535 #endif
536 
537 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
538 
539 #endif
540 
541 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::Tensor< scalar >
L
const vector L(dict.get< vector >("L"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::pointToPointPlanarInterpolation
Interpolates between two sets of unstructured points using 2D Delaunay triangulation....
Definition: pointToPointPlanarInterpolation.H:53
Foam::Vector2D< scalar >
Foam::turbulentDigitalFilterInletFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentDigitalFilterInletFvPatchVectorField.C:1005
Foam::symmTensorField
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
Definition: primitiveFieldsFwd.H:56
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::turbulentDigitalFilterInletFvPatchVectorField::autoMap
virtual void autoMap(const fvPatchFieldMapper &m)
Map (and resize as needed) from self given a mapping object.
Definition: turbulentDigitalFilterInletFvPatchVectorField.C:1057
Foam::turbulentDigitalFilterInletFvPatchVectorField::TypeName
TypeName("turbulentDigitalFilterInlet")
Runtime type information.
Foam::Field< symmTensor >
Foam::turbulentDigitalFilterInletFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentDigitalFilterInletFvPatchVectorField.C:953
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
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
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::turbulentDigitalFilterInletFvPatchVectorField::turbulentDigitalFilterInletFvPatchVectorField
turbulentDigitalFilterInletFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: turbulentDigitalFilterInletFvPatchVectorField.C:646
Foam::turbulentDigitalFilterInletFvPatchVectorField
Definition: turbulentDigitalFilterInletFvPatchVectorField.H:349
Foam::turbulentDigitalFilterInletFvPatchVectorField::rmap
virtual void rmap(const fvPatchVectorField &ptf, const labelList &addr)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: turbulentDigitalFilterInletFvPatchVectorField.C:1066
Random.H
turbulentDigitalFilterInletFvPatchVectorFieldTemplates.C
Foam::fvPatchVectorField
fvPatchField< vector > fvPatchVectorField
Definition: fvPatchFieldsFwd.H:43
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
fixedValueFvPatchFields.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::Tuple2< label, label >
Foam::turbulentDigitalFilterInletFvPatchVectorField::clone
virtual tmp< fvPatchVectorField > clone() const
Construct and return a clone.
Definition: turbulentDigitalFilterInletFvPatchVectorField.H:592
fieldTypes.H
Header files for all the primitive types that Fields are instantiated for.
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::turbulentDigitalFilterInletFvPatchVectorField::~turbulentDigitalFilterInletFvPatchVectorField
virtual ~turbulentDigitalFilterInletFvPatchVectorField()=default
Destructor.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54