source: palm/trunk/SOURCE/synthetic_turbulence_generator_mod.f90 @ 4843

Last change on this file since 4843 was 4843, checked in by raasch, 3 years ago

local namelist parameter added to switch off the module although the respective module namelist appears in the namelist file, further copyright updates

  • Property svn:keywords set to Id
File size: 113.6 KB
Line 
1!> @synthetic_turbulence_generator_mod.f90
2!--------------------------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms of the GNU General
6! Public License as published by the Free Software Foundation, either version 3 of the License, or
7! (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the
10! implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General
11! Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with PALM. If not, see
14! <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2021 Leibniz Universitaet Hannover
17!--------------------------------------------------------------------------------------------------!
18!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: synthetic_turbulence_generator_mod.f90 4843 2021-01-15 15:22:11Z raasch $
27! local namelist parameter added to switch off the module although the respective module namelist
28! appears in the namelist file
29!
30! 4842 2021-01-14 10:42:28Z raasch
31! reading of namelist file and actions in case of namelist errors revised so that statement labels
32! and goto statements are not required any more
33!
34! 4828 2021-01-05 11:21:41Z Giersch
35! Implementation of downward facing USM and LSM surfaces
36!
37! 4647 2020-08-24 16:36:18Z suehring
38! Change default value of synthetic turbulence adjustment as well as compute_velocity_seeds_local
39! By default, the random-seed computation is now distributed among several cores. Especially for
40! large length scales this is significantly faster.
41!
42! 4640 2020-08-11 16:28:32Z suehring
43! - to avoid that the correction term in r11/r22 computation becomes unrealistically high, limit
44!   Obukhov length (term is not valid for near neutral conditions)
45! - to avoid unrealistically large perturbations, change computation of r21 so that this resembles
46!   the vertical transport of horizontal momentum
47!
48! 4629 2020-07-29 09:37:56Z raasch
49! support for MPI Fortran77 interface (mpif.h) removed
50!
51! 4603 2020-07-14 16:08:30Z suehring
52! Bugfix in initialization from ASCII file - x-length scales at the bottom boundary were not
53! initialized properly
54!
55! 4566 2020-06-16 10:11:51Z suehring
56! - revise parametrization for reynolds-stress components, turbulent length- and time scales
57! - revise computation of velocity disturbances to be consistent to Xie and Castro (2008)
58! - change default value of time interval to adjust turbulence parametrization
59! - bugfix in computation of amplitude-tensor (vertical flux of horizontal momentum)
60!
61! 4562 2020-06-12 08:38:47Z raasch
62! Parts of r4559 re-formatted
63!
64! 4559 2020-06-11 08:51:48Z raasch
65! File re-formatted to follow the PALM coding standard
66!
67! 4535 2020-05-15 12:07:23Z raasch
68! Bugfix for restart data format query
69!
70! 4495 2020-04-13 20:11:20Z raasch
71! Restart data handling with MPI-IO added
72!
73! 4481 2020-03-31 18:55:54Z maronga
74! Bugfix: cpp-directives for serial mode added, dummy statements to prevent compile errors added
75!
76! 4442 2020-03-04 19:21:13Z suehring
77! Set back turbulent length scale to 8 x grid spacing in the parametrized mode
78! (was accidantly changed).
79!
80! 4441 2020-03-04 19:20:35Z suehring
81! Correct misplaced preprocessor directive
82!
83! 4438 2020-03-03 20:49:28Z suehring
84! Performance optimizations in velocity-seed calculation:
85!  - Random number array is only defined and computed locally (except for normalization to zero mean
86!    and unit variance)
87!  - Parallel random number generator is applied independent on the 2D random numbers in other
88!    routines
89!  - Option to decide wheter velocity seeds are computed locally without any further communication
90!    or are computed by all processes along the communicator
91!
92! 4346 2019-12-18 11:55:56Z motisi
93! Introduction of wall_flags_total_0, which currently sets bits based on static topography
94! information used in wall_flags_static_0
95!
96! 4335 2019-12-12 16:39:05Z suehring
97! Commentation of last commit
98!
99! 4332 2019-12-10 19:44:12Z suehring
100! Limit initial velocity seeds in restart runs, if not the seed calculation may become unstable.
101! Further, minor bugfix in initial velocity seed calculation.
102!
103! 4329 2019-12-10 15:46:36Z motisi
104! Renamed wall_flags_0 to wall_flags_static_0
105!
106! 4309 2019-11-26 18:49:59Z suehring
107! Computation of velocity seeds optimized. This implies that random numbers are computed now using
108! the parallel random number generator. Random numbers are now only computed and normalized locally,
109! while distributed over all mpi ranks afterwards, instead of computing random numbers on a global
110! array.
111! Further, the number of calls for the time-consuming velocity-seed generation is reduced - now the
112! left and right, as well as the north and south boundary share the same velocity-seed matrices.
113!
114! 4182 2019-08-22 15:20:23Z scharf
115! Corrected "Former revisions" section
116!
117! 4148 2019-08-08 11:26:00Z suehring
118! Remove unused variable
119!
120! 4144 2019-08-06 09:11:47Z raasch
121! Relational operators .EQ., .NE., etc. replaced by ==, /=, etc.
122!
123! 4071 2019-07-03 20:02:00Z suehring
124! Bugfix, initialize mean_inflow_profiles in case turbulence and inflow information is not read from
125! file.
126!
127! 4022 2019-06-12 11:52:39Z suehring
128! Several bugfixes and improvements
129! - Revise bias correction of the imposed perturbations (correction via volume flow can create
130!   instabilities in case the mean volume flow is close to zero)
131! - Introduce lower limits in calculation of coefficient matrix, else the calculation may become
132!   numerically unstable
133! - Impose perturbations every timestep, even though no new set of perturbations is generated in
134!   case dt_stg_call /= dt_3d
135! - Implement a gradual decrease of Reynolds stress and length scales above ABL height (within 1
136!   length scale above ABL depth to 1/10) rather than a discontinuous decrease
137! - Bugfix in non-nested case: use ABL height for parametrized turbulence
138!
139! 3987 2019-05-22 09:52:13Z kanani
140! Introduce alternative switch for debug output during timestepping
141!
142! 3938 2019-04-29 16:06:25Z suehring
143! Remove unused variables
144!
145! 3937 2019-04-29 15:09:07Z suehring
146! Minor bugfix in case of a very early restart where mc_factor is sill not present.
147! Some modification and fixing of potential bugs in the calculation of scaling parameters used for
148! synthetic turbulence parametrization.
149!
150! 3909 2019-04-17 09:13:25Z suehring
151! Minor bugfix for last commit
152!
153! 3900 2019-04-16 15:17:43Z suehring
154! Missing re-calculation of perturbation seeds in case of restarts
155!
156! 3891 2019-04-12 17:52:01Z suehring
157! Bugfix in initialization in case of restart runs.
158!
159! 3885 2019-04-11 11:29:34Z kanani
160! Changes related to global restructuring of location messages and introduction of additional debug
161! messages
162!
163!
164! Removed unused variables
165!
166! 3719 2019-02-06 13:10:18Z kanani
167! Removed log_point measurement from stg_init, since this part is counted to log_point(2)
168! 'initialization' already. Moved other log_points to calls of the subroutines in time_integration
169! for better overview.
170!
171! 2259 2017-06-08 09:09:11Z gronemeier
172! Initial revision
173!
174! Authors:
175! --------
176! @author Tobias Gronemeier, Matthias Suehring, Atsushi Inagaki, Micha Gryschka, Christoph Knigge
177!
178!
179! Description:
180! ------------
181!> The module generates turbulence at the inflow boundary based on a method by Xie and Castro (2008)
182!> utilizing a Lund rotation (Lund, 1998) and a mass-flux correction by Kim et al. (2013).
183!> The turbulence is correlated based on length scales in y- and z-direction and a time scale for
184!> each velocity component. The profiles of length and time scales, mean u, v, w, e and pt, and all
185!> components of the Reynolds stress tensor can be either read from file STG_PROFILES, or will be
186!> parametrized within the boundary layer.
187!>
188!> @todo Enable cyclic_fill
189!>       Implement turbulence generation for e and pt
190!> @note <Enter notes on the module>
191!> @bug  Height information from input file is not used. Profiles from input must match with current
192!>       PALM grid.
193!>       In case of restart, velocity seeds differ from precursor run if a11, a22, or a33 are zero.
194!--------------------------------------------------------------------------------------------------!
195 MODULE synthetic_turbulence_generator_mod
196
197
198    USE arrays_3d,                                                                                 &
199        ONLY:  dzw,                                                                                &
200               ddzw,                                                                               &
201               drho_air,                                                                           &
202               mean_inflow_profiles,                                                               &
203               pt,                                                                                 &
204               pt_init,                                                                            &
205               q,                                                                                  &
206               q_init,                                                                             &
207               u,                                                                                  &
208               u_init,                                                                             &
209               v,                                                                                  &
210               v_init,                                                                             &
211               w,                                                                                  &
212               zu,                                                                                 &
213               zw
214
215    USE basic_constants_and_equations_mod,                                                         &
216        ONLY:  g,                                                                                  &
217               kappa,                                                                              &
218               pi
219
220    USE control_parameters,                                                                        &
221        ONLY:  bc_lr,                                                                              &
222               bc_ns,                                                                              &
223               child_domain,                                                                       &
224               coupling_char,                                                                      &
225               debug_output_timestep,                                                              &
226               dt_3d,                                                                              &
227               e_init,                                                                             &
228               humidity,                                                                           &
229               initializing_actions,                                                               &
230               intermediate_timestep_count,                                                        &
231               intermediate_timestep_count_max,                                                    &
232               length,                                                                             &
233               message_string,                                                                     &
234               nesting_offline,                                                                    &
235               neutral,                                                                            &
236               num_mean_inflow_profiles,                                                           &
237               random_generator,                                                                   &
238               rans_mode,                                                                          &
239               restart_data_format_output,                                                         &
240               restart_string,                                                                     &
241               syn_turb_gen,                                                                       &
242               time_since_reference_point,                                                         &
243               turbulent_inflow
244
245    USE cpulog,                                                                                    &
246        ONLY:  cpu_log,                                                                            &
247               log_point_s
248
249    USE grid_variables,                                                                            &
250        ONLY:  ddx,                                                                                &
251               ddy,                                                                                &
252               dx,                                                                                 &
253               dy
254
255    USE indices,                                                                                   &
256        ONLY:  nbgp,                                                                               &
257               nz,                                                                                 &
258               nzb,                                                                                &
259               nzt,                                                                                &
260               nx,                                                                                 &
261               nxl,                                                                                &
262               nxlu,                                                                               &
263               nxr,                                                                                &
264               ny,                                                                                 &
265               nys,                                                                                &
266               nysv,                                                                               &
267               nyn,                                                                                &
268               wall_flags_total_0
269
270    USE kinds
271
272#if defined( __parallel )
273    USE MPI
274#endif
275
276    USE nesting_offl_mod,                                                                          &
277        ONLY:  nesting_offl_calc_zi,                                                               &
278               zi_ribulk
279
280    USE pegrid,                                                                                    &
281        ONLY:  comm1dx,                                                                            &
282               comm1dy,                                                                            &
283               comm2d,                                                                             &
284               ierr,                                                                               &
285               myidx,                                                                              &
286               myidy,                                                                              &
287               pdims
288
289    USE pmc_interface,                                                                             &
290        ONLY : rans_mode_parent
291
292    USE random_generator_parallel,                                                                 &
293        ONLY:  init_parallel_random_generator,                                                     &
294               random_dummy,                                                                       &
295               random_number_parallel,                                                             &
296               random_seed_parallel
297
298    USE restart_data_mpi_io_mod,                                                                   &
299        ONLY:  rrd_mpi_io,                                                                         &
300               wrd_mpi_io
301
302    USE transpose_indices,                                                                         &
303        ONLY:  nzb_x,                                                                              &
304               nzt_x
305
306    USE surface_mod,                                                                               &
307        ONLY:  surf_def_h,                                                                         &
308               surf_lsm_h,                                                                         &
309               surf_usm_h
310
311    IMPLICIT NONE
312
313    INTEGER(iwp) ::  id_stg_left        !< left lateral boundary core id in case of turbulence generator
314    INTEGER(iwp) ::  id_stg_north       !< north lateral boundary core id in case of turbulence generator
315    INTEGER(iwp) ::  id_stg_right       !< right lateral boundary core id in case of turbulence generator
316    INTEGER(iwp) ::  id_stg_south       !< south lateral boundary core id in case of turbulence generator
317    INTEGER(iwp) ::  k_zi               !< vertical index of boundary-layer top
318    INTEGER(iwp) ::  mergp_limit = 1000 !< maximum length scale (in gp)
319    INTEGER(iwp) ::  mergp_x            !< maximum length scale in x (in gp)
320    INTEGER(iwp) ::  mergp_xy           !< maximum horizontal length scale (in gp)
321    INTEGER(iwp) ::  mergp_y            !< maximum length scale in y (in gp)
322    INTEGER(iwp) ::  mergp_z            !< maximum length scale in z (in gp)
323    INTEGER(iwp) ::  nzb_x_stg          !< lower bound of z coordinate (required for transposing z on PEs along x)
324    INTEGER(iwp) ::  nzt_x_stg          !< upper bound of z coordinate (required for transposing z on PEs along x)
325    INTEGER(iwp) ::  nzb_y_stg          !< lower bound of z coordinate (required for transposing z on PEs along y)
326    INTEGER(iwp) ::  nzt_y_stg          !< upper bound of z coordinate (required for transposing z on PEs along y)
327#if defined( __parallel )
328    INTEGER(iwp) ::  stg_type_xz        !< MPI type for full z range
329    INTEGER(iwp) ::  stg_type_xz_small  !< MPI type for small z range
330    INTEGER(iwp) ::  stg_type_yz        !< MPI type for full z range
331    INTEGER(iwp) ::  stg_type_yz_small  !< MPI type for small z range
332#endif
333
334    INTEGER(iwp), DIMENSION(3) ::  nr_non_topo_xz = 0  !< number of non-topography grid points at xz cross-sections,
335                                                       !< required for bias correction of imposed perturbations
336    INTEGER(iwp), DIMENSION(3) ::  nr_non_topo_yz = 0  !< number of non-topography grid points at yz cross-sections,
337                                                       !< required for bias correction of imposed perturbations
338
339#if defined( __parallel )
340    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  displs_xz      !< displacement for MPI_GATHERV
341    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  recv_count_xz  !< receive count for MPI_GATHERV
342    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  displs_yz      !< displacement for MPI_GATHERV
343    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  recv_count_yz  !< receive count for MPI_GATHERV
344#endif
345    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nux            !< length scale of u in x direction (in gp)
346    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nuy            !< length scale of u in y direction (in gp)
347    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nuz            !< length scale of u in z direction (in gp)
348    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nvx            !< length scale of v in x direction (in gp)
349    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nvy            !< length scale of v in y direction (in gp)
350    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nvz            !< length scale of v in z direction (in gp)
351    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nwx            !< length scale of w in x direction (in gp)
352    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nwy            !< length scale of w in y direction (in gp)
353    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nwz            !< length scale of w in z direction (in gp)
354
355    INTEGER(isp), DIMENSION(:), ALLOCATABLE   ::  id_rand_xz   !< initial random IDs at xz inflow boundary
356    INTEGER(isp), DIMENSION(:), ALLOCATABLE   ::  id_rand_yz   !< initial random IDs at yz inflow boundary
357    INTEGER(isp), DIMENSION(:,:), ALLOCATABLE ::  seq_rand_xz  !< initial random seeds at xz inflow boundary
358    INTEGER(isp), DIMENSION(:,:), ALLOCATABLE ::  seq_rand_yz  !< initial random seeds at yz inflow boundary
359
360
361    LOGICAL ::  adjustment_step               = .FALSE.  !< control flag indicating that time and lenght scales have been updated and
362                                                         !< no time correlation to the timestep before should be considered
363    LOGICAL ::  compute_velocity_seeds_local  = .FALSE.   !< switch to decide whether velocity seeds are computed locally
364                                                         !< or if computation is distributed over several processes
365    LOGICAL ::  parametrize_inflow_turbulence = .FALSE.  !< flag indicating that inflow turbulence is either read from file
366                                                         !< (.FALSE.) or if it parametrized
367    LOGICAL ::  use_syn_turb_gen              = .FALSE.  !< switch to use synthetic turbulence generator
368    LOGICAL ::  velocity_seed_initialized     = .FALSE.  !< true after first call of stg_main
369
370
371    REAL(wp) ::  blend                     !< value to create gradually and smooth blending of Reynolds stress and length
372                                           !< scales above the boundary layer
373    REAL(wp) ::  blend_coeff = -9.3_wp     !< coefficient used to ensure that blending functions decreases to 1/10 after
374                                           !< one length scale above ABL top
375    REAL(wp) ::  d_l                       !< blend_coeff/length_scale
376    REAL(wp) ::  d_nxy                     !< inverse of the total number of xy-grid points
377    REAL(wp) ::  dt_stg_adjust = 1800.0_wp !< time interval for adjusting turbulence statistics
378    REAL(wp) ::  dt_stg_call = 0.0_wp      !< time interval for calling synthetic turbulence generator
379    REAL(wp) ::  scale_l                   !< scaling parameter used for turbulence parametrization - Obukhov length
380    REAL(wp) ::  scale_us                  !< scaling parameter used for turbulence parametrization - friction velocity
381    REAL(wp) ::  scale_wm                  !< scaling parameter used for turbulence parametrization - momentum scale
382    REAL(wp) ::  time_stg_adjust = 0.0_wp  !< time counter for adjusting turbulence information
383    REAL(wp) ::  time_stg_call = 0.0_wp    !< time counter for calling generator
384
385    REAL(wp), DIMENSION(3) ::  mc_factor = 1.0_wp  !< correction factor for the u,v,w-components to maintain original mass flux
386
387
388    REAL(wp),DIMENSION(:), ALLOCATABLE ::  r11  !< Reynolds parameter
389    REAL(wp),DIMENSION(:), ALLOCATABLE ::  r21  !< Reynolds parameter
390    REAL(wp),DIMENSION(:), ALLOCATABLE ::  r22  !< Reynolds parameter
391    REAL(wp),DIMENSION(:), ALLOCATABLE ::  r31  !< Reynolds parameter
392    REAL(wp),DIMENSION(:), ALLOCATABLE ::  r32  !< Reynolds parameter
393    REAL(wp),DIMENSION(:), ALLOCATABLE ::  r33  !< Reynolds parameter
394
395    REAL(wp), DIMENSION(:), ALLOCATABLE ::  a11  !< coefficient for Lund rotation
396    REAL(wp), DIMENSION(:), ALLOCATABLE ::  a21  !< coefficient for Lund rotation
397    REAL(wp), DIMENSION(:), ALLOCATABLE ::  a22  !< coefficient for Lund rotation
398    REAL(wp), DIMENSION(:), ALLOCATABLE ::  a31  !< coefficient for Lund rotation
399    REAL(wp), DIMENSION(:), ALLOCATABLE ::  a32  !< coefficient for Lund rotation
400    REAL(wp), DIMENSION(:), ALLOCATABLE ::  a33  !< coefficient for Lund rotation
401    REAL(wp), DIMENSION(:), ALLOCATABLE ::  tu   !< Lagrangian time scale of u
402    REAL(wp), DIMENSION(:), ALLOCATABLE ::  tv   !< Lagrangian time scale of v
403    REAL(wp), DIMENSION(:), ALLOCATABLE ::  tw   !< Lagrangian time scale of w
404
405    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  bux     !< filter function for u in x direction
406    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  buy     !< filter function for u in y direction
407    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  buz     !< filter function for u in z direction
408    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  bvx     !< filter function for v in x direction
409    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  bvy     !< filter function for v in y direction
410    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  bvz     !< filter function for v in z direction
411    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  bwx     !< filter function for w in y direction
412    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  bwy     !< filter function for w in y direction
413    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  bwz     !< filter function for w in z direction
414    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fu_xz   !< velocity seed for u at xz plane
415    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fuo_xz  !< velocity seed for u at xz plane with new random number
416    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fu_yz   !< velocity seed for u at yz plane
417    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fuo_yz  !< velocity seed for u at yz plane with new random number
418    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fv_xz   !< velocity seed for v at xz plane
419    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fvo_xz  !< velocity seed for v at xz plane with new random number
420    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fv_yz   !< velocity seed for v at yz plane
421    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fvo_yz  !< velocity seed for v at yz plane with new random number
422    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fw_xz   !< velocity seed for w at xz plane
423    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fwo_xz  !< velocity seed for w at xz plane with new random number
424    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fw_yz   !< velocity seed for w at yz plane
425    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fwo_yz  !< velocity seed for w at yz plane with new random number
426
427    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  dist_xz !< disturbances for parallel/crosswind/vertical component at north/south boundary
428    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  dist_yz !< disturbances for parallel/crosswind/vertical component  at north/south boundary
429
430!
431!-- PALM interfaces:
432!-- Adjust time and lenght scales, Reynolds stress, and filter functions
433    INTERFACE stg_adjust
434       MODULE PROCEDURE stg_adjust
435    END INTERFACE stg_adjust
436!
437!-- Input parameter checks to be done in check_parameters
438    INTERFACE stg_check_parameters
439       MODULE PROCEDURE stg_check_parameters
440    END INTERFACE stg_check_parameters
441
442!
443!-- Calculate filter functions
444    INTERFACE stg_filter_func
445       MODULE PROCEDURE stg_filter_func
446    END INTERFACE stg_filter_func
447
448!
449!-- Generate velocity seeds at south and north domain boundary
450    INTERFACE stg_generate_seed_xz
451       MODULE PROCEDURE stg_generate_seed_xz
452    END INTERFACE stg_generate_seed_xz
453!
454!-- Generate velocity seeds at left and/or right domain boundary
455    INTERFACE stg_generate_seed_yz
456       MODULE PROCEDURE stg_generate_seed_yz
457    END INTERFACE stg_generate_seed_yz
458
459!
460!-- Output of information to the header file
461    INTERFACE stg_header
462       MODULE PROCEDURE stg_header
463    END INTERFACE stg_header
464
465!
466!-- Initialization actions
467    INTERFACE stg_init
468       MODULE PROCEDURE stg_init
469    END INTERFACE stg_init
470
471!
472!-- Main procedure of synth. turb. gen.
473    INTERFACE stg_main
474       MODULE PROCEDURE stg_main
475    END INTERFACE stg_main
476
477!
478!-- Reading of NAMELIST parameters
479    INTERFACE stg_parin
480       MODULE PROCEDURE stg_parin
481    END INTERFACE stg_parin
482
483!
484!-- Reading of parameters for restart runs
485    INTERFACE stg_rrd_global
486       MODULE PROCEDURE stg_rrd_global_ftn
487       MODULE PROCEDURE stg_rrd_global_mpi
488    END INTERFACE stg_rrd_global
489
490!
491!-- Writing of binary output for restart runs
492    INTERFACE stg_wrd_global
493       MODULE PROCEDURE stg_wrd_global
494    END INTERFACE stg_wrd_global
495
496    SAVE
497
498    PRIVATE
499
500!
501!-- Public interfaces
502    PUBLIC  stg_adjust,                                                                            &
503            stg_check_parameters,                                                                  &
504            stg_header,                                                                            &
505            stg_init,                                                                              &
506            stg_main,                                                                              &
507            stg_parin,                                                                             &
508            stg_rrd_global,                                                                        &
509            stg_wrd_global
510
511!
512!-- Public variables
513    PUBLIC  dt_stg_call,                                                                           &
514            dt_stg_adjust,                                                                         &
515            id_stg_left,                                                                           &
516            id_stg_north,                                                                          &
517            id_stg_right,                                                                          &
518            id_stg_south,                                                                          &
519            parametrize_inflow_turbulence,                                                         &
520            time_stg_adjust,                                                                       &
521            time_stg_call,                                                                         &
522            use_syn_turb_gen
523
524
525 CONTAINS
526
527
528!--------------------------------------------------------------------------------------------------!
529! Description:
530! ------------
531!> Check parameters routine for synthetic turbulence generator
532!--------------------------------------------------------------------------------------------------!
533 SUBROUTINE stg_check_parameters
534
535    IF ( .NOT. use_syn_turb_gen  .AND.  .NOT. rans_mode  .AND.                                     &
536          nesting_offline )  THEN
537       message_string = 'Synthetic turbulence generator is required ' //                           &
538                        'if offline nesting is applied and PALM operates in LES mode.'
539       CALL message( 'stg_check_parameters', 'PA0520', 0, 0, 0, 6, 0 )
540    ENDIF
541
542    IF ( .NOT. use_syn_turb_gen  .AND.  child_domain                                               &
543         .AND. rans_mode_parent  .AND.  .NOT. rans_mode )  THEN
544       message_string = 'Synthetic turbulence generator is required when nesting is applied ' //   &
545                        'and parent operates in RANS-mode but current child in LES mode.'
546       CALL message( 'stg_check_parameters', 'PA0524', 1, 2, 0, 6, 0 )
547    ENDIF
548
549    IF ( use_syn_turb_gen )  THEN
550
551       IF ( child_domain  .AND.  .NOT. rans_mode  .AND.  .NOT. rans_mode_parent )  THEN
552          message_string = 'Using synthetic turbulence generator is not allowed in LES-LES nesting.'
553          CALL message( 'stg_check_parameters', 'PA0620', 1, 2, 0, 6, 0 )
554
555       ENDIF
556
557       IF ( child_domain  .AND.  rans_mode  .AND.  rans_mode_parent )  THEN
558          message_string = 'Using synthetic turbulence generator is not allowed in RANS-RANS nesting.'
559          CALL message( 'stg_check_parameters', 'PA0621', 1, 2, 0, 6, 0 )
560
561       ENDIF
562
563       IF ( .NOT. nesting_offline  .AND.  .NOT. child_domain )  THEN
564
565          IF ( INDEX( initializing_actions, 'set_constant_profiles' ) == 0                         &
566          .AND.  INDEX( initializing_actions, 'read_restart_data' ) == 0 )  THEN
567             message_string = 'Using synthetic turbulence generator requires ' //                  &
568                              '%initializing_actions = "set_constant_profiles" or ' //             &
569                              ' "read_restart_data", if not offline nesting is applied.'
570             CALL message( 'stg_check_parameters', 'PA0015', 1, 2, 0, 6, 0 )
571          ENDIF
572
573          IF ( bc_lr /= 'dirichlet/radiation' )  THEN
574             message_string = 'Using synthetic turbulence generator requires &bc_lr = ' //         &
575                              ' "dirichlet/radiation", if not offline nesting is applied.'
576             CALL message( 'stg_check_parameters', 'PA0035', 1, 2, 0, 6, 0 )
577          ENDIF
578          IF ( bc_ns /= 'cyclic' )  THEN
579             message_string = 'Using synthetic turbulence generator requires &bc_ns = ' //         &
580                              ' "cyclic", if not offline nesting is applied.'
581             CALL message( 'stg_check_parameters', 'PA0037', 1, 2, 0, 6, 0 )
582          ENDIF
583
584       ENDIF
585
586       IF ( turbulent_inflow )  THEN
587          message_string = 'Using synthetic turbulence generator in combination ' //               &
588                           '&with turbulent_inflow = .T. is not allowed'
589          CALL message( 'stg_check_parameters', 'PA0039', 1, 2, 0, 6, 0 )
590       ENDIF
591!
592!--    Synthetic turbulence generator requires the parallel random generator
593       IF ( random_generator /= 'random-parallel' )  THEN
594          message_string = 'Using synthetic turbulence generator requires random_generator = ' //  &
595                           'random-parallel.'
596          CALL message( 'stg_check_parameters', 'PA0421', 1, 2, 0, 6, 0 )
597       ENDIF
598
599    ENDIF
600
601 END SUBROUTINE stg_check_parameters
602
603
604!--------------------------------------------------------------------------------------------------!
605! Description:
606! ------------
607!> Header output for synthetic turbulence generator
608!--------------------------------------------------------------------------------------------------!
609 SUBROUTINE stg_header ( io )
610
611    INTEGER(iwp), INTENT(IN) ::  io  !< Unit of the output file
612
613!
614!-- Write synthetic turbulence generator Header
615    WRITE( io, 1 )
616    IF ( use_syn_turb_gen )  THEN
617       WRITE( io, 2 )
618    ELSE
619       WRITE( io, 3 )
620    ENDIF
621
622    IF ( parametrize_inflow_turbulence )  THEN
623       WRITE( io, 4 ) dt_stg_adjust
624    ELSE
625       WRITE( io, 5 )
626    ENDIF
627
6281   FORMAT (//' Synthetic turbulence generator information:'/                                      &
629              ' ------------------------------------------'/)
6302   FORMAT ('    synthetic turbulence generator is switched on')
6313   FORMAT ('    synthetic turbulence generator is switched off')
6324   FORMAT ('    imposed turbulence statistics are parametrized and ajdusted to boundary-layer development each ', F8.2, ' s' )
6335   FORMAT ('    imposed turbulence is read from file' )
634
635 END SUBROUTINE stg_header
636
637
638!--------------------------------------------------------------------------------------------------!
639! Description:
640! ------------
641!> Initialization of the synthetic turbulence generator
642!--------------------------------------------------------------------------------------------------!
643 SUBROUTINE stg_init
644
645#if defined( __parallel )
646    INTEGER(KIND=MPI_ADDRESS_KIND) :: extent   !< extent of new MPI type
647    INTEGER(KIND=MPI_ADDRESS_KIND) :: tob = 0  !< dummy variable
648#endif
649
650    INTEGER(iwp) :: i         !> grid index in x-direction
651    INTEGER(iwp) :: j         !> loop index
652    INTEGER(iwp) :: k         !< index
653#if defined( __parallel )
654    INTEGER(iwp) :: newtype   !< dummy MPI type
655    INTEGER(iwp) :: realsize  !< size of REAL variables
656#endif
657
658    INTEGER(iwp), DIMENSION(3) ::  nr_non_topo_xz_l = 0 !< number of non-topography grid points at xz-cross-section on subdomain
659    INTEGER(iwp), DIMENSION(3) ::  nr_non_topo_yz_l = 0 !< number of non-topography grid points at yz-cross-section on subdomain
660
661
662    LOGICAL ::  file_stg_exist = .FALSE.  !< flag indicating whether parameter file for Reynolds stress and length scales exist
663
664!
665!-- Dummy variables used for reading profiles
666    REAL(wp) :: d1   !< u profile
667    REAL(wp) :: d2   !< v profile
668    REAL(wp) :: d3   !< w profile
669    REAL(wp) :: d5   !< e profile
670    REAL(wp) :: luy  !< length scale for u in y direction
671    REAL(wp) :: luz  !< length scale for u in z direction
672    REAL(wp) :: lvy  !< length scale for v in y direction
673    REAL(wp) :: lvz  !< length scale for v in z direction
674    REAL(wp) :: lwy  !< length scale for w in y direction
675    REAL(wp) :: lwz  !< length scale for w in z direction
676#if defined( __parallel )
677    REAL(wp) :: nnz  !< increment used to determine processor decomposition of z-axis along x and y direction
678#endif
679    REAL(wp) :: zz   !< height
680
681
682#if defined( __parallel )
683    CALL MPI_BARRIER( comm2d, ierr )
684#endif
685!
686!-- Create mpi-datatypes for exchange in case of non-local but distributed computation of the
687!-- velocity seeds. This option is useful in case large turbulent length scales are present, where
688!-- the computational effort becomes large and needs to be parallelized. For parameterized turbulence
689!-- the length scales are small and computing the velocity seeds locally is faster (no overhead by
690!-- communication).
691    IF ( .NOT. compute_velocity_seeds_local )  THEN
692#if defined( __parallel )
693!
694!--    Determine processor decomposition of z-axis along x- and y-direction
695       nnz = nz / pdims(1)
696       nzb_x_stg = 1 + myidx * INT( nnz )
697       nzt_x_stg = ( myidx + 1 ) * INT( nnz )
698
699       IF ( MOD( nz , pdims(1) ) /= 0  .AND.  myidx == id_stg_right )                              &
700          nzt_x_stg = nzt_x_stg + myidx * ( nnz - INT( nnz ) )
701
702       IF ( nesting_offline   .OR.  ( child_domain  .AND.  rans_mode_parent                        &
703                               .AND.  .NOT.  rans_mode ) )  THEN
704          nnz = nz / pdims(2)
705          nzb_y_stg = 1 + myidy * INT( nnz )
706          nzt_y_stg = ( myidy + 1 ) * INT( nnz )
707
708          IF ( MOD( nz , pdims(2) ) /= 0  .AND.  myidy == id_stg_north )                           &
709             nzt_y_stg = nzt_y_stg + myidy * ( nnz - INT( nnz ) )
710       ENDIF
711
712!
713!--    Define MPI type used in stg_generate_seed_yz to gather vertical splitted velocity seeds
714       CALL MPI_TYPE_SIZE( MPI_REAL, realsize, ierr )
715       extent = 1 * realsize
716!
717!--    Set-up MPI datatyp to involve all cores for turbulence generation at yz layer
718!--    stg_type_yz: yz-slice with vertical bounds nzb:nzt+1
719       CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt-nzb+2,nyn-nys+1],                                    &
720               [1,nyn-nys+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr )
721       CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_yz, ierr )
722       CALL MPI_TYPE_COMMIT( stg_type_yz, ierr )
723       CALL MPI_TYPE_FREE( newtype, ierr )
724
725       ! stg_type_yz_small: yz-slice with vertical bounds nzb_x_stg:nzt_x_stg+1
726       CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt_x_stg-nzb_x_stg+2,nyn-nys+1],                        &
727               [1,nyn-nys+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr )
728       CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_yz_small, ierr )
729       CALL MPI_TYPE_COMMIT( stg_type_yz_small, ierr )
730       CALL MPI_TYPE_FREE( newtype, ierr )
731
732       ! Receive count and displacement for MPI_GATHERV in stg_generate_seed_yz
733       ALLOCATE( recv_count_yz(pdims(1)), displs_yz(pdims(1)) )
734
735       recv_count_yz           = nzt_x_stg-nzb_x_stg + 1
736       recv_count_yz(pdims(1)) = recv_count_yz(pdims(1)) + 1
737
738       DO  j = 1, pdims(1)
739          displs_yz(j) = 0 + (nzt_x_stg-nzb_x_stg+1) * (j-1)
740       ENDDO
741!
742!--    Set-up MPI datatyp to involve all cores for turbulence generation at xz layer
743!--    stg_type_xz: xz-slice with vertical bounds nzb:nzt+1
744       IF ( nesting_offline  .OR.  ( child_domain .AND.  rans_mode_parent      &
745                              .AND.  .NOT.  rans_mode ) )  THEN
746          CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt-nzb+2,nxr-nxl+1],             &
747                  [1,nxr-nxl+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr )
748          CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_xz, ierr )
749          CALL MPI_TYPE_COMMIT( stg_type_xz, ierr )
750          CALL MPI_TYPE_FREE( newtype, ierr )
751
752          ! stg_type_yz_small: xz-slice with vertical bounds nzb_x_stg:nzt_x_stg+1
753          CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt_y_stg-nzb_y_stg+2,nxr-nxl+1], &
754                  [1,nxr-nxl+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr )
755          CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_xz_small, ierr )
756          CALL MPI_TYPE_COMMIT( stg_type_xz_small, ierr )
757          CALL MPI_TYPE_FREE( newtype, ierr )
758
759          ! Receive count and displacement for MPI_GATHERV in stg_generate_seed_yz
760          ALLOCATE( recv_count_xz(pdims(2)), displs_xz(pdims(2)) )
761
762          recv_count_xz           = nzt_y_stg-nzb_y_stg + 1
763          recv_count_xz(pdims(2)) = recv_count_xz(pdims(2)) + 1
764
765          DO  j = 1, pdims(2)
766             displs_xz(j) = 0 + (nzt_y_stg-nzb_y_stg+1) * (j-1)
767          ENDDO
768
769       ENDIF
770#endif
771    ENDIF
772!
773!-- Allocate required arrays.
774!-- In case no offline nesting or self nesting is used, the arrary mean_inflow profiles is required.
775!-- Check if it is already allocated, else allocate and initialize it appropriately. Note, in case
776!-- turbulence and inflow information is read from file, mean_inflow_profiles is already allocated
777!-- and initialized appropriately.
778    IF ( .NOT. nesting_offline  .AND.  .NOT.  child_domain )  THEN
779       IF ( .NOT. ALLOCATED( mean_inflow_profiles ) )  THEN
780          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,1:num_mean_inflow_profiles) )
781          mean_inflow_profiles = 0.0_wp
782          mean_inflow_profiles(:,1) = u_init
783          mean_inflow_profiles(:,2) = v_init
784!
785!--       Even though potential temperature and humidity are not perturbed, they need to be
786!--       initialized appropriately.
787          IF ( .NOT. neutral )                                                                     &
788             mean_inflow_profiles(:,4) = pt_init
789          IF ( humidity )                                                                          &
790             mean_inflow_profiles(:,6) = q_init
791       ENDIF
792    ENDIF
793!
794!-- Assign initial profiles. Note, this is only required if turbulent inflow from the left is
795!-- desired, not in case of any of the nesting (offline or self nesting) approaches.
796    IF ( .NOT. nesting_offline  .AND.  .NOT.  child_domain )  THEN
797       u_init = mean_inflow_profiles(:,1)
798       v_init = mean_inflow_profiles(:,2)
799       e_init = MAXVAL( mean_inflow_profiles(:,5) )
800    ENDIF
801
802    ALLOCATE ( a11(nzb:nzt+1), a21(nzb:nzt+1), a22(nzb:nzt+1),                                     &
803               a31(nzb:nzt+1), a32(nzb:nzt+1), a33(nzb:nzt+1),                                     &
804               nux(nzb:nzt+1), nuy(nzb:nzt+1), nuz(nzb:nzt+1),                                     &
805               nvx(nzb:nzt+1), nvy(nzb:nzt+1), nvz(nzb:nzt+1),                                     &
806               nwx(nzb:nzt+1), nwy(nzb:nzt+1), nwz(nzb:nzt+1),                                     &
807               r11(nzb:nzt+1), r21(nzb:nzt+1), r22(nzb:nzt+1),                                     &
808               r31(nzb:nzt+1), r32(nzb:nzt+1), r33(nzb:nzt+1),                                     &
809               tu(nzb:nzt+1),  tv(nzb:nzt+1),  tw(nzb:nzt+1) )
810
811    r11 = 0.0_wp
812    r21 = 0.0_wp
813    r22 = 0.0_wp
814    r31 = 0.0_wp
815    r32 = 0.0_wp
816    r33 = 0.0_wp
817    tu  = 0.0_wp
818    tv  = 0.0_wp
819    tw  = 0.0_wp
820
821    ALLOCATE ( dist_xz(nzb:nzt+1,nxl:nxr,3) )
822    ALLOCATE ( dist_yz(nzb:nzt+1,nys:nyn,3) )
823    dist_xz = 0.0_wp
824    dist_yz = 0.0_wp
825!
826!-- Read inflow profile
827!-- Height levels of profiles in input profile is as follows:
828!-- zu: luy, luz, tu, lvy, lvz, tv, r11, r21, r22, d1, d2, d5 zw: lwy, lwz, tw, r31, r32, r33, d3
829!-- WARNING: zz is not used at the moment
830    INQUIRE( FILE = 'STG_PROFILES' // TRIM( coupling_char ), EXIST = file_stg_exist  )
831
832    IF ( file_stg_exist )  THEN
833
834       OPEN( 90, FILE = 'STG_PROFILES' // TRIM( coupling_char ), STATUS = 'OLD', FORM = 'FORMATTED' )
835!
836!--    Skip header
837       READ( 90, * )
838
839       DO  k = nzb, nzt
840          READ( 90, * ) zz, luy, luz, tu(k), lvy, lvz, tv(k), lwy, lwz, tw(k), r11(k), r21(k),     &
841                        r22(k), r31(k), r32(k), r33(k), d1, d2, d3, d5
842
843!
844!--       Convert length scales from meter to number of grid points.
845          IF ( k /= nzb )  THEN
846             nuz(k) = INT( luz * ddzw(k) )
847             nvz(k) = INT( lvz * ddzw(k) )
848             nwz(k) = INT( lwz * ddzw(k) )
849          ELSE
850             nuz(k) = INT( luz * ddzw(k+1) )
851             nvz(k) = INT( lvz * ddzw(k+1) )
852             nwz(k) = INT( lwz * ddzw(k+1) )
853          ENDIF
854
855          nuy(k) = INT( luy * ddy )
856          nvy(k) = INT( lvy * ddy )
857          nwy(k) = INT( lwy * ddy )
858!
859!--       Set length scales in x-direction. As a workaround assume isotropic turbulence in x- and
860!--       y-direction.
861          nwx(k) = nwy(k)
862          nvx(k) = nvy(k)
863          nux(k) = nuy(k)
864!
865!--       Save Mean inflow profiles
866          IF ( TRIM( initializing_actions ) /= 'read_restart_data' ) THEN
867             mean_inflow_profiles(k,1) = d1
868             mean_inflow_profiles(k,2) = d2
869            !  mean_inflow_profiles(k,4) = d4
870             mean_inflow_profiles(k,5) = d5
871          ENDIF
872       ENDDO
873!
874!--    Set length scales at the surface and top boundary. At the surface the lengths scales are
875!--    simply overwritten.
876       nwx(nzb) = nwy(nzb+1)
877       nvx(nzb) = nvy(nzb+1)
878       nux(nzb) = nuy(nzb+1)
879       nuy(nzb) = nuy(nzb+1)
880       nuz(nzb) = nuz(nzb+1)
881       nvy(nzb) = nvy(nzb+1)
882       nvz(nzb) = nvz(nzb+1)
883       nwy(nzb) = nwy(nzb+1)
884       nwz(nzb) = nwz(nzb+1)
885
886       nwx(nzt+1) = nwy(nzt)
887       nvx(nzt+1) = nvy(nzt)
888       nux(nzt+1) = nuy(nzt)
889       nuy(nzt+1) = nuy(nzt)
890       nuz(nzt+1) = nuz(nzt)
891       nvy(nzt+1) = nvy(nzt)
892       nvz(nzt+1) = nvz(nzt)
893       nwy(nzt+1) = nwy(nzt)
894       nwz(nzt+1) = nwz(nzt)
895
896       CLOSE( 90 )
897!
898!--    Calculate coefficient matrix from Reynolds stress tensor (Lund rotation)
899       CALL calc_coeff_matrix
900!
901!-- No information about turbulence and its length scales are available. Instead, parametrize
902!-- turbulence which is imposed at the boundaries. Set flag which indicates that turbulence is
903!-- parametrized, which is done later when energy-balance models are already initialized. This is
904!-- because the STG needs information about surface properties such as roughness to build
905!-- 'realistic' turbulence profiles.
906    ELSE
907!
908!--    Precalulate the inversion number of xy-grid points - later used for mass conservation
909       d_nxy = 1.0_wp / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
910!
911!--    Set flag indicating that turbulence is parametrized
912       parametrize_inflow_turbulence = .TRUE.
913!
914!--    In case of dirichlet inflow boundary conditions only at one lateral boundary, i.e. in the
915!--    case no offline or self nesting is applied but synthetic turbulence shall be parametrized
916!--    nevertheless, the boundary-layer depth needs to be determined first.
917       IF ( .NOT. nesting_offline  .AND.  .NOT.  child_domain )                                    &
918          CALL nesting_offl_calc_zi
919!
920!--    Determine boundary-layer depth, which is used to initialize lenght scales
921       CALL calc_scaling_variables
922!
923!--    Parametrize Reynolds-stress tensor, diagonal elements as well as r21 (v'u'), r31 (w'u'),
924!--    r32 (w'v'). Parametrization follows Rotach et al. (1996) and is based on boundary-layer depth,
925!--    friction velocity and velocity scale.
926       CALL parametrize_turbulence
927!
928!--    Calculate coefficient matrix from Reynolds stress tensor (Lund rotation)
929       CALL calc_coeff_matrix
930
931    ENDIF
932!
933!-- Initial filter functions
934    CALL stg_setup_filter_function
935!
936!-- Allocate velocity seeds for turbulence at xz-layer
937    ALLOCATE ( fu_xz(nzb:nzt+1,nxl:nxr), fuo_xz(nzb:nzt+1,nxl:nxr),                                &
938               fv_xz(nzb:nzt+1,nxl:nxr), fvo_xz(nzb:nzt+1,nxl:nxr),                                &
939               fw_xz(nzb:nzt+1,nxl:nxr), fwo_xz(nzb:nzt+1,nxl:nxr) )
940
941!
942!-- Allocate velocity seeds for turbulence at yz-layer
943    ALLOCATE ( fu_yz(nzb:nzt+1,nys:nyn), fuo_yz(nzb:nzt+1,nys:nyn),                                &
944               fv_yz(nzb:nzt+1,nys:nyn), fvo_yz(nzb:nzt+1,nys:nyn),                                &
945               fw_yz(nzb:nzt+1,nys:nyn), fwo_yz(nzb:nzt+1,nys:nyn) )
946
947    fu_xz  = 0.0_wp
948    fuo_xz = 0.0_wp
949    fv_xz  = 0.0_wp
950    fvo_xz = 0.0_wp
951    fw_xz  = 0.0_wp
952    fwo_xz = 0.0_wp
953
954    fu_yz  = 0.0_wp
955    fuo_yz = 0.0_wp
956    fv_yz  = 0.0_wp
957    fvo_yz = 0.0_wp
958    fw_yz  = 0.0_wp
959    fwo_yz = 0.0_wp
960
961#if defined( __parallel )
962    CALL MPI_BARRIER( comm2d, ierr )
963#endif
964
965!
966!-- In case of restart, calculate velocity seeds fu, fv, fw from former time step.
967!   Bug: fu, fv, fw are different in those heights where a11, a22, a33 are 0 compared to the prerun.
968!-- This is mostly for k=nzt+1.
969    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
970       IF ( myidx == id_stg_left  .OR.  myidx == id_stg_right )  THEN
971
972          IF ( myidx == id_stg_left  )  i = -1
973          IF ( myidx == id_stg_right )  i = nxr+1
974
975          DO  j = nys, nyn
976             DO  k = nzb, nzt+1
977                IF  ( a11(k) > 10E-8_wp )  THEN
978                   fu_yz(k,j) = ( u(k,j,i) - u_init(k) ) / a11(k)
979                ELSE
980                   fu_yz(k,j) = 10E-8_wp
981                ENDIF
982
983                IF  ( a22(k) > 10E-8_wp )  THEN
984                   fv_yz(k,j) = ( v(k,j,i) - a21(k) * fu_yz(k,j) - v_init(k) ) / a22(k)
985                ELSE
986                   fv_yz(k,j) = 10E-8_wp
987                ENDIF
988
989                IF  ( a33(k) > 10E-8_wp )  THEN
990                   fw_yz(k,j) = ( w(k,j,i) - a31(k) * fu_yz(k,j) - a32(k) * fv_yz(k,j) ) / a33(k)
991                ELSE
992                   fw_yz(k,j) = 10E-8_wp
993                ENDIF
994             ENDDO
995          ENDDO
996       ENDIF
997
998       IF ( myidy == id_stg_south  .OR.  myidy == id_stg_north )  THEN
999
1000          IF ( myidy == id_stg_south )  j = -1
1001          IF ( myidy == id_stg_north )  j = nyn+1
1002
1003          DO  i = nxl, nxr
1004             DO  k = nzb, nzt+1
1005!
1006!--             In case the correlation coefficients are very small, the velocity seeds may become
1007!--             very large finally creating numerical instabilities in the synthetic turbulence
1008!--             generator. Empirically, a value of 10E-8 seems to be sufficient.
1009                IF  ( a11(k) > 10E-8_wp )  THEN
1010                   fu_xz(k,i) = ( u(k,j,i) - u_init(k) ) / a11(k)
1011                ELSE
1012                   fu_xz(k,i) = 10E-8_wp
1013                ENDIF
1014
1015                IF  ( a22(k) > 10E-8_wp )  THEN
1016                   fv_xz(k,i) = ( v(k,j,i) - a21(k) * fu_xz(k,i) - v_init(k) ) / a22(k)
1017                ELSE
1018                   fv_xz(k,i) = 10E-8_wp
1019                ENDIF
1020
1021                IF  ( a33(k) > 10E-8_wp )  THEN
1022                   fw_xz(k,i) = ( w(k,j,i) - a31(k) * fu_xz(k,i) - a32(k) * fv_xz(k,i) ) / a33(k)
1023                ELSE
1024                   fw_xz(k,i) = 10E-8_wp
1025                ENDIF
1026
1027             ENDDO
1028          ENDDO
1029       ENDIF
1030    ENDIF
1031!
1032!-- Count the number of non-topography grid points at the boundaries where perturbations are imposed.
1033!-- This number is later used for bias corrections of the perturbations, i.e. to force their mean to
1034!-- be zero. Please note, due to the asymetry of u and v along x and y direction, respectively,
1035!-- different cases must be distinguished.
1036    IF ( myidx == id_stg_left  .OR.  myidx == id_stg_right )  THEN
1037!
1038!--    Number of grid points where perturbations are imposed on u
1039       IF ( myidx == id_stg_left  )  i = nxl
1040       IF ( myidx == id_stg_right )  i = nxr+1
1041
1042       nr_non_topo_yz_l(1) = SUM( MERGE( 1, 0, BTEST( wall_flags_total_0(nzb:nzt,nys:nyn,i), 1 ) ) )
1043!
1044!--    Number of grid points where perturbations are imposed on v and w
1045       IF ( myidx == id_stg_left  )  i = nxl-1
1046       IF ( myidx == id_stg_right )  i = nxr+1
1047
1048       nr_non_topo_yz_l(2) = SUM( MERGE( 1, 0, BTEST( wall_flags_total_0(nzb:nzt,nysv:nyn,i), 2 ) ) )
1049       nr_non_topo_yz_l(3) = SUM( MERGE( 1, 0, BTEST( wall_flags_total_0(nzb:nzt,nys:nyn,i), 3 ) ) )
1050
1051#if defined( __parallel )
1052       CALL MPI_ALLREDUCE( nr_non_topo_yz_l, nr_non_topo_yz, 3, MPI_INTEGER, MPI_SUM, comm1dy, ierr )
1053#else
1054       nr_non_topo_yz = nr_non_topo_yz_l
1055#endif
1056    ENDIF
1057
1058    IF ( myidy == id_stg_south  .OR.  myidy == id_stg_north )  THEN
1059!
1060!--    Number of grid points where perturbations are imposed on v
1061       IF ( myidy == id_stg_south )  j = nys
1062       IF ( myidy == id_stg_north )  j = nyn+1
1063
1064       nr_non_topo_xz_l(2) = SUM( MERGE( 1, 0, BTEST( wall_flags_total_0(nzb:nzt,j,nxl:nxr), 2 ) ) )
1065!
1066!--    Number of grid points where perturbations are imposed on u and w
1067       IF ( myidy == id_stg_south )  j = nys-1
1068       IF ( myidy == id_stg_north )  j = nyn+1
1069
1070       nr_non_topo_xz_l(1) = SUM( MERGE( 1, 0, BTEST( wall_flags_total_0(nzb:nzt,j,nxlu:nxr), 1 ) ) )
1071       nr_non_topo_xz_l(3) = SUM( MERGE( 1, 0, BTEST( wall_flags_total_0(nzb:nzt,j,nxl:nxr), 3 ) ) )
1072
1073#if defined( __parallel )
1074       CALL MPI_ALLREDUCE( nr_non_topo_xz_l, nr_non_topo_xz, 3, MPI_INTEGER, MPI_SUM, comm1dx, ierr )
1075#else
1076       nr_non_topo_xz = nr_non_topo_xz_l
1077#endif
1078    ENDIF
1079!
1080!-- Initialize random number generator at xz- and yz-layers. Random numbers are initialized at each
1081!-- core. In case there is only inflow from the left, it is sufficient to generate only random
1082!-- numbers for the yz-layer, else random numbers for the xz-layer are also required.
1083    ALLOCATE ( id_rand_yz(-mergp_limit+nys:nyn+mergp_limit) )
1084    ALLOCATE ( seq_rand_yz(5,-mergp_limit+nys:nyn+mergp_limit) )
1085    id_rand_yz  = 0
1086    seq_rand_yz = 0
1087
1088    CALL init_parallel_random_generator( ny, -mergp_limit+nys, nyn+mergp_limit, id_rand_yz, seq_rand_yz )
1089
1090    IF ( nesting_offline  .OR.  ( child_domain .AND.  rans_mode_parent                             &
1091         .AND.  .NOT.  rans_mode ) )  THEN
1092       ALLOCATE ( id_rand_xz(-mergp_limit+nxl:nxr+mergp_limit) )
1093       ALLOCATE ( seq_rand_xz(5,-mergp_limit+nxl:nxr+mergp_limit) )
1094       id_rand_xz  = 0
1095       seq_rand_xz = 0
1096
1097       CALL init_parallel_random_generator( nx, -mergp_limit+nxl, nxr+mergp_limit, id_rand_xz, seq_rand_xz )
1098    ENDIF
1099
1100 END SUBROUTINE stg_init
1101
1102
1103!--------------------------------------------------------------------------------------------------!
1104! Description:
1105! ------------
1106!> Calculate filter function bxx from length scale nxx following Eg.9 and 10 (Xie and Castro, 2008)
1107!--------------------------------------------------------------------------------------------------!
1108 SUBROUTINE stg_filter_func( nxx, bxx, mergp )
1109
1110    INTEGER(iwp) ::  k     !< loop index
1111    INTEGER(iwp) ::  mergp !< passed length scale in grid points
1112    INTEGER(iwp) ::  n_k   !< length scale nXX in height k
1113    INTEGER(iwp) ::  nf    !< index for length scales
1114
1115    INTEGER(iwp), DIMENSION(nzb:nzt+1) ::  nxx  !< length scale (in gp)
1116
1117    REAL(wp) :: bdenom        !< denominator for filter functions bXX
1118    REAL(wp) :: qsi = 1.0_wp  !< minimization factor
1119
1120    REAL(wp), DIMENSION(-mergp:mergp,nzb:nzt+1) ::  bxx  !< filter function
1121
1122
1123    bxx = 0.0_wp
1124
1125    DO  k = nzb, nzt+1
1126       bdenom = 0.0_wp
1127       n_k    = nxx(k)
1128       IF ( n_k /= 0 )  THEN
1129
1130!
1131!--       ( Eq.10 )^2
1132          DO  nf = -n_k, n_k
1133             bdenom = bdenom + EXP( -qsi * pi * ABS( nf ) / n_k )**2
1134          ENDDO
1135
1136!
1137!--       ( Eq.9 )
1138          bdenom = SQRT( bdenom )
1139          DO  nf = -n_k, n_k
1140             bxx(nf,k) = EXP( -qsi * pi * ABS( nf ) / n_k ) / bdenom
1141          ENDDO
1142       ENDIF
1143    ENDDO
1144
1145 END SUBROUTINE stg_filter_func
1146
1147
1148!------------------------------------------------------------------------------!
1149! Description:
1150! ------------
1151!> Calculate filter function bxx from length scale nxx following Eg.9 and 10
1152!> (Xie and Castro, 2008)
1153!------------------------------------------------------------------------------!
1154 SUBROUTINE stg_setup_filter_function
1155
1156    INTEGER(iwp) :: j  !< dummy value to calculate maximum length scale index
1157    INTEGER(iwp) :: k  !< loop index along vertical direction
1158!
1159!-- Define the size of the filter functions and allocate them.
1160    mergp_x = 0
1161    mergp_y = 0
1162    mergp_z = 0
1163
1164    DO  k = nzb, nzt+1
1165       j = MAX( ABS(nux(k)), ABS(nvx(k)), ABS(nwx(k)) )
1166       IF ( j > mergp_x )  mergp_x = j
1167    ENDDO
1168    DO  k = nzb, nzt+1
1169       j = MAX( ABS(nuy(k)), ABS(nvy(k)), ABS(nwy(k)) )
1170       IF ( j > mergp_y )  mergp_y = j
1171    ENDDO
1172    DO  k = nzb, nzt+1
1173       j = MAX( ABS(nuz(k)), ABS(nvz(k)), ABS(nwz(k)) )
1174       IF ( j > mergp_z )  mergp_z = j
1175    ENDDO
1176
1177    mergp_xy = MAX( mergp_x, mergp_y )
1178
1179    IF ( ALLOCATED( bux ) )  DEALLOCATE( bux )
1180    IF ( ALLOCATED( bvx ) )  DEALLOCATE( bvx )
1181    IF ( ALLOCATED( bwx ) )  DEALLOCATE( bwx )
1182    IF ( ALLOCATED( buy ) )  DEALLOCATE( buy )
1183    IF ( ALLOCATED( bvy ) )  DEALLOCATE( bvy )
1184    IF ( ALLOCATED( bwy ) )  DEALLOCATE( bwy )
1185    IF ( ALLOCATED( buz ) )  DEALLOCATE( buz )
1186    IF ( ALLOCATED( bvz ) )  DEALLOCATE( bvz )
1187    IF ( ALLOCATED( bwz ) )  DEALLOCATE( bwz )
1188
1189    ALLOCATE ( bux(-mergp_x:mergp_x,nzb:nzt+1),                                                    &
1190               buy(-mergp_y:mergp_y,nzb:nzt+1),                                                    &
1191               buz(-mergp_z:mergp_z,nzb:nzt+1),                                                    &
1192               bvx(-mergp_x:mergp_x,nzb:nzt+1),                                                    &
1193               bvy(-mergp_y:mergp_y,nzb:nzt+1),                                                    &
1194               bvz(-mergp_z:mergp_z,nzb:nzt+1),                                                    &
1195               bwx(-mergp_x:mergp_x,nzb:nzt+1),                                                    &
1196               bwy(-mergp_y:mergp_y,nzb:nzt+1),                                                    &
1197               bwz(-mergp_z:mergp_z,nzb:nzt+1) )
1198!
1199!-- Create filter functions
1200    CALL stg_filter_func( nux, bux, mergp_x ) !filter ux
1201    CALL stg_filter_func( nuy, buy, mergp_y ) !filter uy
1202    CALL stg_filter_func( nuz, buz, mergp_z ) !filter uz
1203    CALL stg_filter_func( nvx, bvx, mergp_x ) !filter vx
1204    CALL stg_filter_func( nvy, bvy, mergp_y ) !filter vy
1205    CALL stg_filter_func( nvz, bvz, mergp_z ) !filter vz
1206    CALL stg_filter_func( nwx, bwx, mergp_x ) !filter wx
1207    CALL stg_filter_func( nwy, bwy, mergp_y ) !filter wy
1208    CALL stg_filter_func( nwz, bwz, mergp_z ) !filter wz
1209
1210 END SUBROUTINE stg_setup_filter_function
1211
1212!--------------------------------------------------------------------------------------------------!
1213! Description:
1214! ------------
1215!> Parin for &stg_par for synthetic turbulence generator
1216!--------------------------------------------------------------------------------------------------!
1217 SUBROUTINE stg_parin
1218
1219    CHARACTER(LEN=100) ::  line  !< dummy string that contains the current line of the parameter file
1220
1221    INTEGER(iwp) ::  io_status   !< status after reading the namelist file
1222
1223    LOGICAL ::  switch_off_module = .FALSE.  !< local namelist parameter to switch off the module
1224                                             !< although the respective module namelist appears in
1225                                             !< the namelist file
1226
1227    NAMELIST /stg_par/  compute_velocity_seeds_local,                                              &
1228                        dt_stg_adjust,                                                             &
1229                        dt_stg_call,                                                               &
1230                        switch_off_module,                                                         &
1231                        use_syn_turb_gen
1232
1233
1234!
1235!-- Move to the beginning of the namelist file and try to find and read the namelist.
1236    REWIND( 11 )
1237    READ( 11, stg_par, IOSTAT=io_status )
1238
1239!
1240!-- Action depending on the READ status
1241    IF ( io_status == 0 )  THEN
1242!
1243!--    stg_par namelist was found and read correctly. Set flag that indicates that the synthetic
1244!--    turbulence generator is switched on.
1245       IF ( .NOT. switch_off_module )  syn_turb_gen = .TRUE.
1246
1247    ELSEIF ( io_status > 0 )  THEN
1248!
1249!--    stg_par namelist was found but contained errors. Print an error message including the line
1250!--    that caused the problem.
1251       BACKSPACE( 11 )
1252       READ( 11 , '(A)') line
1253       CALL parin_fail_message( 'stg_par', line )
1254
1255    ENDIF
1256
1257 END SUBROUTINE stg_parin
1258
1259
1260!--------------------------------------------------------------------------------------------------!
1261! Description:
1262! ------------
1263!> Read module-specific global restart data (Fortran binary format).
1264!--------------------------------------------------------------------------------------------------!
1265 SUBROUTINE stg_rrd_global_ftn( found )
1266
1267    LOGICAL, INTENT(OUT) ::  found !< flag indicating if variable was found
1268
1269    found = .TRUE.
1270
1271
1272    SELECT CASE ( restart_string(1:length) )
1273
1274       CASE ( 'time_stg_adjust' )
1275          READ ( 13 )  time_stg_adjust
1276
1277       CASE ( 'time_stg_call' )
1278          READ ( 13 )  time_stg_call
1279
1280       CASE ( 'use_syn_turb_gen' )
1281          READ ( 13 )  use_syn_turb_gen
1282
1283       CASE DEFAULT
1284
1285          found = .FALSE.
1286
1287    END SELECT
1288
1289
1290 END SUBROUTINE stg_rrd_global_ftn
1291
1292
1293!--------------------------------------------------------------------------------------------------!
1294! Description:
1295! ------------
1296!> Read module-specific global restart data (MPI-IO).
1297!--------------------------------------------------------------------------------------------------!
1298 SUBROUTINE stg_rrd_global_mpi
1299
1300    CALL rrd_mpi_io( 'time_stg_adjust', time_stg_adjust )
1301    CALL rrd_mpi_io( 'time_stg_call', time_stg_call )
1302    CALL rrd_mpi_io( 'use_syn_turb_gen', use_syn_turb_gen )
1303
1304 END SUBROUTINE stg_rrd_global_mpi
1305
1306
1307!--------------------------------------------------------------------------------------------------!
1308! Description:
1309! ------------
1310!> This routine writes the respective restart data.
1311!--------------------------------------------------------------------------------------------------!
1312 SUBROUTINE stg_wrd_global
1313
1314    IF ( TRIM( restart_data_format_output ) == 'fortran_binary' )  THEN
1315
1316      CALL wrd_write_string( 'time_stg_adjust' )
1317       WRITE ( 14 )  time_stg_adjust
1318
1319       CALL wrd_write_string( 'time_stg_call' )
1320       WRITE ( 14 )  time_stg_call
1321
1322       CALL wrd_write_string( 'use_syn_turb_gen' )
1323       WRITE ( 14 )  use_syn_turb_gen
1324
1325    ELSEIF ( restart_data_format_output(1:3) == 'mpi' )  THEN
1326
1327       CALL wrd_mpi_io( 'time_stg_adjust', time_stg_adjust )
1328       CALL wrd_mpi_io( 'time_stg_call', time_stg_call )
1329       CALL wrd_mpi_io( 'use_syn_turb_gen', use_syn_turb_gen )
1330
1331    ENDIF
1332
1333 END SUBROUTINE stg_wrd_global
1334
1335
1336!--------------------------------------------------------------------------------------------------!
1337! Description:
1338! ------------
1339!> Create turbulent inflow fields for u, v, w with prescribed length scales and Reynolds stress
1340!> tensor after a method of Xie and Castro (2008), modified following suggestions of Kim et al.
1341!> (2013), and using a Lund rotation (Lund, 1998).
1342!--------------------------------------------------------------------------------------------------!
1343 SUBROUTINE stg_main
1344
1345    USE exchange_horiz_mod,                                                                        &
1346        ONLY:  exchange_horiz
1347
1348    INTEGER(iwp) ::  i  !< grid index in x-direction
1349    INTEGER(iwp) ::  j  !< loop index in y-direction
1350    INTEGER(iwp) ::  k  !< loop index in z-direction
1351
1352    LOGICAL  ::  stg_call = .FALSE.  !< control flag indicating whether turbulence was updated or only restored from last call
1353
1354    REAL(wp) ::  dt_stg  !< time interval the STG is called
1355
1356    REAL(wp), DIMENSION(3) ::  mc_factor_l  !< local mass flux correction factor
1357
1358    IF ( debug_output_timestep )  CALL debug_message( 'stg_main', 'start' )
1359!
1360!-- Calculate time step which is needed for filter functions
1361    dt_stg = MAX( dt_3d, dt_stg_call )
1362!
1363!-- Check if synthetic turbulence generator needs to be called and new perturbations need to be
1364!-- created or if old disturbances can be imposed again.
1365    IF ( time_stg_call >= dt_stg_call  .AND.                                                       &
1366         intermediate_timestep_count == intermediate_timestep_count_max  )  THEN
1367       stg_call = .TRUE.
1368    ELSE
1369       stg_call = .FALSE.
1370    ENDIF
1371!
1372!-- Initial value of fu, fv, fw
1373    IF ( time_since_reference_point == 0.0_wp .AND. .NOT. velocity_seed_initialized )  THEN
1374!
1375!--    Generate turbulence at the left and right boundary. Random numbers for the yz-planes at the
1376!--    left/right boundary are generated by the left-sided mpi ranks only. After random numbers are
1377!--    calculated, they are distributed to all other mpi ranks in the model, so that the velocity
1378!--    seed matrices are available on all mpi ranks (i.e. also on the right-sided boundary mpi ranks).
1379!--    In case of offline nesting, this implies, that the left- and the right-sided lateral boundary
1380!--    have the same initial seeds.
1381!--    Note, in case of inflow from the right only, only turbulence at the left boundary is required.
1382       IF ( .NOT. ( nesting_offline  .OR.  ( child_domain .AND.  rans_mode_parent                  &
1383            .AND.  .NOT.  rans_mode ) ) )  THEN
1384          CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fu_yz, id_stg_left )
1385          CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fv_yz, id_stg_left )
1386          CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fw_yz, id_stg_left )
1387       ELSE
1388          CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fu_yz, id_stg_left, id_stg_right )
1389          CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fv_yz, id_stg_left, id_stg_right )
1390          CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fw_yz, id_stg_left, id_stg_right )
1391!
1392!--       Generate turbulence at the south and north boundary. Random numbers for the xz-planes at
1393!--       the south/north boundary are generated by the south-sided mpi ranks only. Please see also
1394!--       comment above.
1395          CALL stg_generate_seed_xz( nux, nuz, bux, buz, fu_xz, id_stg_south, id_stg_north )
1396          CALL stg_generate_seed_xz( nvx, nvz, bvx, bvz, fv_xz, id_stg_south, id_stg_north )
1397          CALL stg_generate_seed_xz( nwx, nwz, bwx, bwz, fw_xz, id_stg_south, id_stg_north )
1398       ENDIF
1399       velocity_seed_initialized = .TRUE.
1400    ENDIF
1401!
1402!-- New set of fu, fv, fw. Note, for inflow from the left side only, velocity seeds are only
1403!-- required at the left boundary, while in case of offline nesting or RANS-LES nesting velocity
1404!-- seeds are required also at the right, south and north boundaries.
1405    IF ( stg_call )  THEN
1406       IF ( .NOT. ( nesting_offline  .OR.                                                          &
1407                  ( child_domain .AND.  rans_mode_parent                                           &
1408                                 .AND.  .NOT.  rans_mode ) ) )  THEN
1409          CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fuo_yz, id_stg_left )
1410          CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fvo_yz, id_stg_left )
1411          CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fwo_yz, id_stg_left )
1412
1413       ELSE
1414          CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fuo_yz, id_stg_left, id_stg_right )
1415          CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fvo_yz, id_stg_left, id_stg_right )
1416          CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fwo_yz, id_stg_left, id_stg_right )
1417
1418          CALL stg_generate_seed_xz( nux, nuz, bux, buz, fuo_xz, id_stg_south, id_stg_north )
1419          CALL stg_generate_seed_xz( nvx, nvz, bvx, bvz, fvo_xz, id_stg_south, id_stg_north )
1420          CALL stg_generate_seed_xz( nwx, nwz, bwx, bwz, fwo_xz, id_stg_south, id_stg_north )
1421       ENDIF
1422    ENDIF
1423
1424!
1425!-- Turbulence generation at left and/or right boundary
1426    IF ( myidx == id_stg_left  .OR.  myidx == id_stg_right )  THEN
1427!
1428!--    Calculate new set of perturbations. Do this only at last RK3-substep and when dt_stg_call is
1429!--    exceeded. Else the old set of perturbations is imposed
1430       IF ( stg_call )  THEN
1431
1432          DO  j = nys, nyn
1433             DO  k = nzb, nzt + 1
1434!
1435!--             Update fu, fv, fw following Eq. 14 of Xie and Castro (2008)
1436                IF ( tu(k) == 0.0_wp  .OR.  adjustment_step )  THEN
1437                   fu_yz(k,j) = fuo_yz(k,j)
1438                ELSE
1439                   fu_yz(k,j) = fu_yz(k,j) * EXP( -pi * dt_stg * 0.5_wp / tu(k) ) +                &
1440                                fuo_yz(k,j) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tu(k) ) )
1441                ENDIF
1442
1443                IF ( tv(k) == 0.0_wp  .OR.  adjustment_step )  THEN
1444                   fv_yz(k,j) = fvo_yz(k,j)
1445                ELSE
1446                   fv_yz(k,j) = fv_yz(k,j) * EXP( -pi * dt_stg * 0.5_wp / tv(k) ) +                &
1447                                fvo_yz(k,j) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tv(k) ) )
1448                ENDIF
1449
1450                IF ( tw(k) == 0.0_wp  .OR.  adjustment_step )  THEN
1451                   fw_yz(k,j) = fwo_yz(k,j)
1452                ELSE
1453                   fw_yz(k,j) = fw_yz(k,j) * EXP( -pi * dt_stg * 0.5_wp / tw(k) ) +                &
1454                                fwo_yz(k,j) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tw(k) ) )
1455                ENDIF
1456             ENDDO
1457          ENDDO
1458
1459          dist_yz(nzb,:,1) = 0.0_wp
1460          dist_yz(nzb,:,2) = 0.0_wp
1461          dist_yz(nzb,:,3) = 0.0_wp
1462
1463          IF ( myidx == id_stg_left  )  i = nxl
1464          IF ( myidx == id_stg_right )  i = nxr+1
1465          DO  j = nys, nyn
1466             DO  k = nzb+1, nzt + 1
1467!
1468!--             Lund rotation following Eq. 17 in Xie and Castro (2008).
1469                dist_yz(k,j,1) = a11(k) * fu_yz(k,j) *                                             &
1470                                 MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 1 ) )
1471
1472             ENDDO
1473          ENDDO
1474
1475          IF ( myidx == id_stg_left  )  i = nxl-1
1476          IF ( myidx == id_stg_right )  i = nxr+1
1477          DO  j = nys, nyn
1478             DO  k = nzb+1, nzt + 1
1479!
1480!--             Lund rotation following Eq. 17 in Xie and Castro (2008).
1481!--             Additional factors are added to improve the variance of v and w experimental test
1482!--             of 1.2
1483                dist_yz(k,j,2) = ( a21(k) * fu_yz(k,j) + a22(k) * fv_yz(k,j) ) *                   &
1484                                 MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 2 ) )
1485                dist_yz(k,j,3) = ( a31(k) * fu_yz(k,j) + a32(k) * fv_yz(k,j) +                     &
1486                                   a33(k) * fw_yz(k,j) ) *                                         &
1487                                 MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 3 ) )
1488
1489             ENDDO
1490          ENDDO
1491       ENDIF
1492!
1493!--    Mass flux correction following Kim et al. (2013)
1494!--    This correction factor ensures that the mass flux is preserved at the inflow boundary. First,
1495!--    calculate mean value of the imposed perturbations at yz boundary. Note, this needs to be done
1496!--    only after the last RK3-substep, else the boundary values won't be accessed.
1497       IF ( intermediate_timestep_count == intermediate_timestep_count_max )  THEN
1498          mc_factor_l = 0.0_wp
1499          mc_factor   = 0.0_wp
1500!
1501!--       Sum up the original volume flows (with and without perturbations).
1502!--       Note, for non-normal components (here v and w) it is actually no volume flow.
1503          IF ( myidx == id_stg_left  )  i = nxl
1504          IF ( myidx == id_stg_right )  i = nxr+1
1505
1506          mc_factor_l(1) = SUM( dist_yz(nzb:nzt,nys:nyn,1) *                                       &
1507                       MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(nzb:nzt,nys:nyn,i), 1 ) ) )
1508
1509          IF ( myidx == id_stg_left  )  i = nxl-1
1510          IF ( myidx == id_stg_right )  i = nxr+1
1511
1512          mc_factor_l(2) = SUM( dist_yz(nzb:nzt,nysv:nyn,2) *                                      &
1513                       MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(nzb:nzt,nysv:nyn,i), 2 ) ) )
1514          mc_factor_l(3) = SUM( dist_yz(nzb:nzt,nys:nyn,3) *                                       &
1515                       MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(nzb:nzt,nys:nyn,i), 3 ) ) )
1516
1517#if defined( __parallel )
1518          CALL MPI_ALLREDUCE( mc_factor_l, mc_factor, 3, MPI_REAL, MPI_SUM, comm1dy, ierr )
1519#else
1520          mc_factor   = mc_factor_l
1521#endif
1522!
1523!--       Calculate correction factor and force zero mean perturbations.
1524          mc_factor = mc_factor / REAL( nr_non_topo_yz, KIND = wp )
1525
1526          IF ( myidx == id_stg_left  )  i = nxl
1527          IF ( myidx == id_stg_right )  i = nxr+1
1528
1529          dist_yz(:,nys:nyn,1) = ( dist_yz(:,nys:nyn,1) - mc_factor(1) ) *                         &
1530                                MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(:,nys:nyn,i), 1 ) )
1531
1532
1533          IF ( myidx == id_stg_left  )  i = nxl-1
1534          IF ( myidx == id_stg_right )  i = nxr+1
1535
1536          dist_yz(:,nys:nyn,2) = ( dist_yz(:,nys:nyn,2) - mc_factor(2) ) *                         &
1537                                MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(:,nys:nyn,i), 2 ) )
1538
1539          dist_yz(:,nys:nyn,3) = ( dist_yz(:,nys:nyn,3) - mc_factor(3) ) *                         &
1540                                MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(:,nys:nyn,i), 3 ) )
1541!
1542!--       Add disturbances
1543          IF ( myidx == id_stg_left  )  THEN
1544!
1545!--          For the left boundary distinguish between mesoscale offline / self nesting and
1546!--          turbulent inflow at the left boundary only. In the latter case turbulence is imposed on
1547!--          the mean inflow profiles.
1548             IF ( .NOT. nesting_offline  .AND.  .NOT. child_domain )  THEN
1549!
1550!--             Add disturbance at the inflow
1551                DO  j = nys, nyn
1552                   DO  k = nzb, nzt+1
1553                      u(k,j,-nbgp+1:0) = ( mean_inflow_profiles(k,1) + dist_yz(k,j,1) ) *          &
1554                                     MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,0), 1 ) )
1555                      v(k,j,-nbgp:-1)  = ( mean_inflow_profiles(k,2) + dist_yz(k,j,2) ) *          &
1556                                     MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,-1), 2 ) )
1557                      w(k,j,-nbgp:-1)  = dist_yz(k,j,3) *                                          &
1558                                     MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,-1), 3 ) )
1559                   ENDDO
1560                ENDDO
1561             ELSE
1562
1563                DO  j = nys, nyn
1564                   DO  k = nzb+1, nzt
1565                      u(k,j,0)   = ( u(k,j,0)  + dist_yz(k,j,1) ) *                                &
1566                                   MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,0), 1 ) )
1567                      u(k,j,-1)  = u(k,j,0)
1568                      v(k,j,-1)  = ( v(k,j,-1) + dist_yz(k,j,2) ) *                                &
1569                                   MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,-1), 2 ) )
1570                      w(k,j,-1)  = ( w(k,j,-1) + dist_yz(k,j,3) ) *                                &
1571                                   MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,-1), 3 ) )
1572                   ENDDO
1573                ENDDO
1574             ENDIF
1575          ENDIF
1576          IF ( myidx == id_stg_right  )  THEN
1577             DO  j = nys, nyn
1578                DO  k = nzb+1, nzt
1579                   u(k,j,nxr+1) = ( u(k,j,nxr+1) + dist_yz(k,j,1) ) *                              &
1580                                  MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,nxr+1), 1 ) )
1581                   v(k,j,nxr+1) = ( v(k,j,nxr+1) + dist_yz(k,j,2) ) *                              &
1582                                  MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,nxr+1), 2 ) )
1583                   w(k,j,nxr+1) = ( w(k,j,nxr+1) + dist_yz(k,j,3) ) *                              &
1584                                  MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,nxr+1), 3 ) )
1585                ENDDO
1586             ENDDO
1587          ENDIF
1588       ENDIF
1589    ENDIF
1590!
1591!-- Turbulence generation at north and south boundary
1592    IF ( myidy == id_stg_north  .OR.  myidy == id_stg_south )  THEN
1593!
1594!--    Calculate new set of perturbations. Do this only at last RK3-substep and when dt_stg_call is
1595!--    exceeded. Else the old set of perturbations is imposed
1596       IF ( stg_call )  THEN
1597          DO  i = nxl, nxr
1598             DO  k = nzb, nzt + 1
1599!
1600!--             Update fu, fv, fw following Eq. 14 of Xie and Castro (2008)
1601                IF ( tu(k) == 0.0_wp  .OR.  adjustment_step )  THEN
1602                   fu_xz(k,i) = fuo_xz(k,i)
1603                ELSE
1604                   fu_xz(k,i) = fu_xz(k,i) * EXP( -pi * dt_stg * 0.5_wp / tu(k) ) +                &
1605                                fuo_xz(k,i) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tu(k) ) )
1606                ENDIF
1607
1608                IF ( tv(k) == 0.0_wp  .OR.  adjustment_step )  THEN
1609                   fv_xz(k,i) = fvo_xz(k,i)
1610                ELSE
1611                   fv_xz(k,i) = fv_xz(k,i) * EXP( -pi * dt_stg * 0.5_wp / tv(k) ) +                &
1612                                fvo_xz(k,i) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tv(k) ) )
1613                ENDIF
1614
1615                IF ( tw(k) == 0.0_wp  .OR.  adjustment_step )  THEN
1616                   fw_xz(k,i) = fwo_xz(k,i)
1617                ELSE
1618                   fw_xz(k,i) = fw_xz(k,i) * EXP( -pi * dt_stg * 0.5_wp / tw(k) ) +                &
1619                                fwo_xz(k,i) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tw(k) ) )
1620                ENDIF
1621             ENDDO
1622          ENDDO
1623
1624
1625          dist_xz(nzb,:,1) = 0.0_wp
1626          dist_xz(nzb,:,2) = 0.0_wp
1627          dist_xz(nzb,:,3) = 0.0_wp
1628
1629          IF ( myidy == id_stg_south  ) j = nys
1630          IF ( myidy == id_stg_north )  j = nyn+1
1631          DO  i = nxl, nxr
1632             DO  k = nzb+1, nzt + 1
1633!
1634!--             Lund rotation following Eq. 17 in Xie and Castro (2008).
1635!--             Additional factors are added to improve the variance of v and w
1636                !experimental test of 1.2
1637                dist_xz(k,i,2) = ( a21(k) * fu_xz(k,i) + a22(k) * fv_xz(k,i) ) *                   &
1638                                 MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 2 ) )
1639             ENDDO
1640          ENDDO
1641
1642          IF ( myidy == id_stg_south  ) j = nys-1
1643          IF ( myidy == id_stg_north )  j = nyn+1
1644          DO  i = nxl, nxr
1645             DO  k = nzb+1, nzt + 1
1646!
1647!--             Lund rotation following Eq. 17 in Xie and Castro (2008).
1648!--             Additional factors are added to improve the variance of v and w
1649                dist_xz(k,i,1) = a11(k) * fu_xz(k,i) *                                             &
1650                                 MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 1 ) )
1651
1652                dist_xz(k,i,3) = ( a31(k) * fu_xz(k,i) + a32(k) * fv_xz(k,i) +                     &
1653                                   a33(k) * fw_xz(k,i) ) *                                         &
1654                                 MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 3 ) )
1655             ENDDO
1656          ENDDO
1657       ENDIF
1658
1659!
1660!--    Mass flux correction following Kim et al. (2013). This correction factor ensures that the
1661!--    mass flux is preserved at the inflow boundary. First, calculate mean value of the imposed
1662!--    perturbations at yz boundary. Note, this needs to be done only after the last RK3-substep,
1663!--    else the boundary values won't be accessed.
1664       IF ( intermediate_timestep_count == intermediate_timestep_count_max )  THEN
1665          mc_factor_l = 0.0_wp
1666          mc_factor   = 0.0_wp
1667
1668          IF ( myidy == id_stg_south  ) j = nys
1669          IF ( myidy == id_stg_north )  j = nyn+1
1670
1671          mc_factor_l(2) = SUM( dist_xz(nzb:nzt,nxl:nxr,2) *                                       &
1672                       MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(nzb:nzt,j,nxl:nxr), 2 ) ) )
1673
1674          IF ( myidy == id_stg_south ) j = nys-1
1675          IF ( myidy == id_stg_north ) j = nyn+1
1676
1677          mc_factor_l(1) = SUM( dist_xz(nzb:nzt,nxlu:nxr,1) *                                      &
1678                       MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(nzb:nzt,j,nxlu:nxr), 1 ) ) )
1679          mc_factor_l(3) = SUM( dist_xz(nzb:nzt,nxl:nxr,3) *                                       &
1680                       MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(nzb:nzt,j,nxl:nxr), 3 ) ) )
1681
1682#if defined( __parallel )
1683          CALL MPI_ALLREDUCE( mc_factor_l, mc_factor, 3, MPI_REAL, MPI_SUM, comm1dx, ierr )
1684#else
1685          mc_factor   = mc_factor_l
1686#endif
1687
1688          mc_factor = mc_factor / REAL( nr_non_topo_xz, KIND = wp )
1689
1690          IF ( myidy == id_stg_south ) j = nys
1691          IF ( myidy == id_stg_north ) j = nyn+1
1692
1693          dist_xz(:,nxl:nxr,2) = ( dist_xz(:,nxl:nxr,2) - mc_factor(2) ) *                         &
1694                                MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(:,j,nxl:nxr), 2 ) )
1695
1696
1697          IF ( myidy == id_stg_south ) j = nys-1
1698          IF ( myidy == id_stg_north ) j = nyn+1
1699
1700          dist_xz(:,nxl:nxr,1) = ( dist_xz(:,nxl:nxr,1) - mc_factor(1) ) *                         &
1701                                MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(:,j,nxl:nxr), 1 ) )
1702
1703          dist_xz(:,nxl:nxr,3) = ( dist_xz(:,nxl:nxr,3) - mc_factor(3) ) *                         &
1704                                MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(:,j,nxl:nxr), 3 ) )
1705!
1706!--       Add disturbances
1707          IF ( myidy == id_stg_south )  THEN
1708             DO  i = nxl, nxr
1709                DO  k = nzb+1, nzt
1710                   u(k,-1,i) = ( u(k,-1,i) + dist_xz(k,i,1) )                                      &
1711                               * MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,-1,i), 1 ) )
1712                   v(k,0,i)  = ( v(k,0,i)  + dist_xz(k,i,2) )                                      &
1713                               * MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,0,i), 2 ) )
1714                   v(k,-1,i) = v(k,0,i)
1715                   w(k,-1,i) = ( w(k,-1,i) + dist_xz(k,i,3) )                                      &
1716                               * MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,-1,i), 3 ) )
1717                ENDDO
1718             ENDDO
1719          ENDIF
1720          IF ( myidy == id_stg_north  )  THEN
1721
1722             DO  i = nxl, nxr
1723                DO  k = nzb+1, nzt
1724                   u(k,nyn+1,i) = ( u(k,nyn+1,i) + dist_xz(k,i,1) ) *                              &
1725                                  MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,nyn+1,i), 1 ) )
1726                   v(k,nyn+1,i) = ( v(k,nyn+1,i) + dist_xz(k,i,2) ) *                              &
1727                                  MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,nyn+1,i), 2 ) )
1728                   w(k,nyn+1,i) = ( w(k,nyn+1,i) + dist_xz(k,i,3) ) *                              &
1729                                  MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,nyn+1,i), 3 ) )
1730                ENDDO
1731             ENDDO
1732          ENDIF
1733       ENDIF
1734    ENDIF
1735!
1736!-- Exchange ghost points.
1737    CALL exchange_horiz( u, nbgp )
1738    CALL exchange_horiz( v, nbgp )
1739    CALL exchange_horiz( w, nbgp )
1740!
1741!-- Finally, set time counter for calling STG to zero
1742    IF ( stg_call )  time_stg_call = 0.0_wp
1743!
1744!-- Set adjustment step to False to indicate that time correlation can be
1745!-- switched-on again.
1746    adjustment_step = .FALSE.
1747
1748    IF ( debug_output_timestep )  CALL debug_message( 'stg_main', 'end' )
1749
1750 END SUBROUTINE stg_main
1751
1752!--------------------------------------------------------------------------------------------------!
1753! Description:
1754! ------------
1755!> Generate a set of random number rand_it wich is equal on each PE and calculate the velocity seed
1756!> f_n. f_n is splitted in vertical direction by the number of PEs in x-direction and each PE
1757!> calculates a vertical subsection of f_n. At the the end, all parts are collected to form the full
1758!> array.
1759!--------------------------------------------------------------------------------------------------!
1760 SUBROUTINE stg_generate_seed_yz( n_y, n_z, b_y, b_z, f_n, id_left, id_right )
1761
1762    INTEGER(iwp)           ::  id_left     !< core ids at respective boundaries
1763    INTEGER(iwp), OPTIONAL ::  id_right    !< core ids at respective boundaries
1764    INTEGER(iwp)           ::  j           !< loop index in y-direction
1765    INTEGER(iwp)           ::  jj          !< loop index in y-direction
1766    INTEGER(iwp)           ::  k           !< loop index in z-direction
1767    INTEGER(iwp)           ::  kk          !< loop index in z-direction
1768    INTEGER(iwp)           ::  send_count  !< send count for MPI_GATHERV
1769
1770    INTEGER(iwp), DIMENSION(nzb:nzt+1) ::  n_y  !< length scale in y-direction
1771    INTEGER(iwp), DIMENSION(nzb:nzt+1) ::  n_z  !< length scale in z-direction
1772
1773    REAL(wp) ::  nyz_inv         !< inverse of number of grid points in yz-slice
1774    REAL(wp) ::  rand_av         !< average of random number
1775    REAL(wp) ::  rand_sigma_inv  !< inverse of stdev of random number
1776
1777    REAL(wp), DIMENSION(-mergp_y:mergp_y,nzb:nzt+1) :: b_y     !< filter function in y-direction
1778    REAL(wp), DIMENSION(-mergp_z:mergp_z,nzb:nzt+1) :: b_z     !< filter function in z-direction
1779
1780    REAL(wp), DIMENSION(nzb_x_stg:nzt_x_stg+1,nys:nyn) ::  f_n_l  !<  local velocity seed
1781    REAL(wp), DIMENSION(nzb:nzt+1,nys:nyn)             ::  f_n    !<  velocity seed
1782
1783    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  rand_it  !< global array of random numbers
1784!
1785!-- Generate random numbers using the parallel random generator. The set of random numbers are
1786!-- modified to have an average of 0 and unit variance. Note, at the end the random number array
1787!-- must be defined globally in order to compute the correlation matrices. However, the calculation
1788!-- and normalization of random numbers is done locally, while the result is later distributed to
1789!-- all processes. Further, please note, a set of random numbers is only calculated for the west
1790!-- boundary, while the east boundary uses the same random numbers and thus also computes the same
1791!-- correlation matrix.
1792    ALLOCATE( rand_it(nzb-mergp_z:nzt+1+mergp_z,-mergp_y+nys:nyn+mergp_y) )
1793    rand_it = 0.0_wp
1794
1795    rand_av        = 0.0_wp
1796    rand_sigma_inv = 0.0_wp
1797    nyz_inv        = 1.0_wp / REAL( ( nzt + 1 + mergp_z - ( nzb - mergp_z ) + 1 )                  &
1798                     * ( ny + mergp_y - ( 0 - mergp_y ) + 1 ), KIND = wp )
1799!
1800!-- Compute and normalize random numbers.
1801    DO  j = nys - mergp_y, nyn + mergp_y
1802!
1803!--    Put the random seeds at grid point j
1804       CALL random_seed_parallel( put=seq_rand_yz(:,j) )
1805       DO  k = nzb - mergp_z, nzt + 1 + mergp_z
1806          CALL random_number_parallel( random_dummy )
1807          rand_it(k,j) = random_dummy
1808       ENDDO
1809!
1810!--    Get the new random seeds from last call at grid point j
1811       CALL random_seed_parallel( get=seq_rand_yz(:,j) )
1812    ENDDO
1813!
1814!-- For normalization to zero mean, sum-up the global random numers. To normalize the global set of
1815!-- random numbers, the inner ghost layers mergp must not be summed-up, else the random numbers on
1816!-- the ghost layers will be stronger weighted as they also occur on the inner subdomains.
1817    DO  j = MERGE( nys, nys - mergp_y, nys /= 0 ), MERGE( nyn, nyn + mergp_y, nyn /= ny )
1818       DO  k = nzb - mergp_z, nzt + 1 + mergp_z
1819          rand_av = rand_av + rand_it(k,j)
1820       ENDDO
1821    ENDDO
1822
1823#if defined( __parallel )
1824!
1825!-- Sum-up the local averages of the random numbers
1826    CALL MPI_ALLREDUCE( MPI_IN_PLACE, rand_av, 1, MPI_REAL, MPI_SUM, comm1dy, ierr )
1827#endif
1828    rand_av = rand_av * nyz_inv
1829!
1830!-- Obtain zero mean
1831    rand_it= rand_it - rand_av
1832!
1833!-- Now, compute the variance
1834    DO  j = MERGE( nys, nys - mergp_y, nys /= 0 ), MERGE( nyn, nyn + mergp_y, nyn /= ny )
1835       DO  k = nzb - mergp_z, nzt + 1 + mergp_z
1836          rand_sigma_inv = rand_sigma_inv + rand_it(k,j)**2
1837       ENDDO
1838    ENDDO
1839
1840#if defined( __parallel )
1841!
1842!-- Sum-up the local quadratic averages of the random numbers
1843    CALL MPI_ALLREDUCE( MPI_IN_PLACE, rand_sigma_inv, 1, MPI_REAL, MPI_SUM, comm1dy, ierr )
1844#endif
1845!
1846!-- Compute standard deviation
1847    IF ( rand_sigma_inv /= 0.0_wp )  THEN
1848       rand_sigma_inv = 1.0_wp / SQRT( rand_sigma_inv * nyz_inv )
1849    ELSE
1850       rand_sigma_inv = 1.0_wp
1851    ENDIF
1852!
1853!-- Normalize with standard deviation to obtain unit variance
1854    rand_it = rand_it * rand_sigma_inv
1855
1856    CALL cpu_log( log_point_s(31), 'STG f_n factors', 'start' )
1857!
1858!-- Generate velocity seed following Eq.6 of Xie and Castro (2008). There are two options. In the
1859!-- first one, the computation of the seeds is distributed to all processes along the communicator
1860!-- comm1dy and gathered on the leftmost and, if necessary, on the rightmost process. For huge
1861!-- length scales the computational effort can become quite huge (it scales with the turbulent
1862!-- length scales), so that gain by parallelization exceeds the costs by the subsequent
1863!-- communication. In the second option, which performs better when the turbulent length scales
1864!-- are parametrized and thus the loops are smaller, the seeds are computed locally and no
1865!-- communication is necessary.
1866    IF ( compute_velocity_seeds_local )  THEN
1867
1868       f_n  = 0.0_wp
1869       DO  j = nys, nyn
1870          DO  k = nzb, nzt+1
1871             DO  jj = -n_y(k), n_y(k)
1872                DO  kk = -n_z(k), n_z(k)
1873                   f_n(k,j) = f_n(k,j) + b_y(jj,k) * b_z(kk,k) * rand_it(k+kk,j+jj)
1874                ENDDO
1875             ENDDO
1876          ENDDO
1877       ENDDO
1878
1879    ELSE
1880
1881       f_n_l  = 0.0_wp
1882       DO  j = nys, nyn
1883          DO  k = nzb_x_stg, nzt_x_stg+1
1884             DO  jj = -n_y(k), n_y(k)
1885                DO  kk = -n_z(k), n_z(k)
1886                   f_n_l(k,j) = f_n_l(k,j) + b_y(jj,k) * b_z(kk,k) * rand_it(k+kk,j+jj)
1887                ENDDO
1888             ENDDO
1889          ENDDO
1890       ENDDO
1891!
1892!--    Gather velocity seeds of full subdomain
1893       send_count = nzt_x_stg - nzb_x_stg + 1
1894       IF ( nzt_x_stg == nzt )  send_count = send_count + 1
1895
1896#if defined( __parallel )
1897!
1898!--    Gather the velocity seed matrix on left boundary mpi ranks.
1899       CALL MPI_GATHERV( f_n_l(nzb_x_stg,nys), send_count, stg_type_yz_small, f_n(nzb+1,nys),      &
1900                         recv_count_yz, displs_yz, stg_type_yz, id_left, comm1dx, ierr )
1901!
1902!--    If required, gather the same velocity seed matrix on right boundary mpi ranks (in offline
1903!--    nesting for example).
1904       IF ( PRESENT( id_right ) )  THEN
1905          CALL MPI_GATHERV( f_n_l(nzb_x_stg,nys), send_count, stg_type_yz_small, f_n(nzb+1,nys),   &
1906                            recv_count_yz, displs_yz, stg_type_yz, id_right, comm1dx, ierr )
1907       ENDIF
1908#else
1909       f_n(nzb+1:nzt+1,nys:nyn) = f_n_l(nzb_x_stg:nzt_x_stg+1,nys:nyn)
1910!
1911!--    Next line required to avoid compile errors because of unused dummy arguments
1912       IF ( id_left == 0 )  id_right = 0
1913#endif
1914
1915    ENDIF
1916
1917    DEALLOCATE( rand_it )
1918
1919    CALL cpu_log( log_point_s(31), 'STG f_n factors', 'stop' )
1920
1921 END SUBROUTINE stg_generate_seed_yz
1922
1923
1924!--------------------------------------------------------------------------------------------------!
1925! Description:
1926! ------------
1927!> Generate a set of random number rand_it wich is equal on each PE and calculate the velocity seed
1928!> f_n.
1929!> f_n is splitted in vertical direction by the number of PEs in y-direction and and each PE
1930!> calculates a vertical subsection of f_n. At the the end, all parts are collected to form the
1931!> full array.
1932!--------------------------------------------------------------------------------------------------!
1933 SUBROUTINE stg_generate_seed_xz( n_x, n_z, b_x, b_z, f_n, id_south, id_north )
1934
1935    INTEGER(iwp) :: i           !< loop index in x-direction
1936    INTEGER(iwp) :: id_north    !< core ids at respective boundaries
1937    INTEGER(iwp) :: id_south    !< core ids at respective boundaries
1938    INTEGER(iwp) :: ii          !< loop index in x-direction
1939    INTEGER(iwp) :: k           !< loop index in z-direction
1940    INTEGER(iwp) :: kk          !< loop index in z-direction
1941    INTEGER(iwp) :: send_count  !< send count for MPI_GATHERV
1942
1943    INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_x    !< length scale in x-direction
1944    INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_z    !< length scale in z-direction
1945
1946    REAL(wp) :: nxz_inv         !< inverse of number of grid points in xz-slice
1947    REAL(wp) :: rand_av         !< average of random number
1948    REAL(wp) :: rand_sigma_inv  !< inverse of stdev of random number
1949
1950    REAL(wp), DIMENSION(-mergp_x:mergp_x,nzb:nzt+1) :: b_x     !< filter function in x-direction
1951    REAL(wp), DIMENSION(-mergp_z:mergp_z,nzb:nzt+1) :: b_z     !< filter function in z-direction
1952
1953    REAL(wp), DIMENSION(nzb_y_stg:nzt_y_stg+1,nxl:nxr) :: f_n_l   !<  local velocity seed
1954    REAL(wp), DIMENSION(nzb:nzt+1,nxl:nxr)             :: f_n     !<  velocity seed
1955
1956    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  rand_it   !< global array of random numbers
1957
1958!
1959!-- Generate random numbers using the parallel random generator. The set of random numbers are
1960!-- modified to have an average of 0 and unit variance. Note, at the end the random number array
1961!-- must be defined globally in order to compute the correlation matrices. However, the calculation
1962!-- and normalization of random numbers is done locally, while the result is later distributed to
1963!-- all processes. Further, please note, a set of random numbers is only calculated for the south
1964!-- boundary, while the north boundary uses the same random numbers and thus also computes the same
1965!-- correlation matrix.
1966    ALLOCATE( rand_it(nzb-mergp_z:nzt+1+mergp_z,-mergp_x+nxl:nxr+mergp_x) )
1967    rand_it = 0.0_wp
1968
1969    rand_av        = 0.0_wp
1970    rand_sigma_inv = 0.0_wp
1971    nxz_inv        = 1.0_wp / REAL( ( nzt + 1 + mergp_z - ( nzb - mergp_z ) + 1 )                  &
1972                     * ( nx + mergp_x - ( 0 - mergp_x ) +1 ), KIND = wp )
1973!
1974!-- Compute and normalize random numbers.
1975    DO  i = nxl - mergp_x, nxr + mergp_x
1976!
1977!--    Put the random seeds at grid point ii
1978       CALL random_seed_parallel( put=seq_rand_xz(:,i) )
1979       DO  k = nzb - mergp_z, nzt + 1 + mergp_z
1980          CALL random_number_parallel( random_dummy )
1981          rand_it(k,i) = random_dummy
1982       ENDDO
1983!
1984!--    Get the new random seeds from last call at grid point ii
1985       CALL random_seed_parallel( get=seq_rand_xz(:,i) )
1986    ENDDO
1987!
1988!-- For normalization to zero mean, sum-up the global random numers.
1989!-- To normalize the global set of random numbers, the inner ghost layers mergp must not be
1990!-- summed-up, else the random numbers on the ghost layers will be stronger weighted as they
1991!-- also occur on the inner subdomains.
1992    DO  i = MERGE( nxl, nxl - mergp_x, nxl /= 0 ), MERGE( nxr, nxr + mergp_x, nxr /= nx )
1993       DO  k = nzb - mergp_z, nzt + 1 + mergp_z
1994          rand_av = rand_av + rand_it(k,i)
1995       ENDDO
1996    ENDDO
1997
1998#if defined( __parallel )
1999!
2000!-- Sum-up the local averages of the random numbers
2001    CALL MPI_ALLREDUCE( MPI_IN_PLACE, rand_av, 1, MPI_REAL, MPI_SUM, comm1dx, ierr )
2002#endif
2003    rand_av = rand_av * nxz_inv
2004!
2005!-- Obtain zero mean
2006    rand_it= rand_it - rand_av
2007!
2008!-- Now, compute the variance
2009    DO  i = MERGE( nxl, nxl - mergp_x, nxl /= 0 ), MERGE( nxr, nxr + mergp_x, nxr /= nx )
2010       DO  k = nzb - mergp_z, nzt + 1 + mergp_z
2011          rand_sigma_inv = rand_sigma_inv + rand_it(k,i)**2
2012       ENDDO
2013    ENDDO
2014
2015#if defined( __parallel )
2016!
2017!-- Sum-up the local quadratic averages of the random numbers
2018    CALL MPI_ALLREDUCE( MPI_IN_PLACE, rand_sigma_inv, 1, MPI_REAL, MPI_SUM, comm1dx, ierr )
2019#endif
2020!
2021!-- Compute standard deviation
2022    IF ( rand_sigma_inv /= 0.0_wp )  THEN
2023       rand_sigma_inv = 1.0_wp / SQRT( rand_sigma_inv * nxz_inv )
2024    ELSE
2025       rand_sigma_inv = 1.0_wp
2026    ENDIF
2027!
2028!-- Normalize with standard deviation to obtain unit variance
2029    rand_it = rand_it * rand_sigma_inv
2030
2031    CALL cpu_log( log_point_s(31), 'STG f_n factors', 'start' )
2032!
2033!-- Generate velocity seed following Eq.6 of Xie and Castro (2008). There are two options. In the
2034!-- first one, the computation of the seeds is distributed to all processes along the communicator
2035!-- comm1dx and gathered on the southmost and, if necessary, on the northmost process. For huge
2036!-- length scales the computational effort can become quite huge (it scales with the turbulent
2037!-- length scales), so that gain by parallelization exceeds the costs by the subsequent communication.
2038!-- In the second option, which performs better when the turbulent length scales are parametrized
2039!-- and thus the loops are smaller, the seeds are computed locally and no communication is necessary.
2040    IF ( compute_velocity_seeds_local )  THEN
2041
2042       f_n  = 0.0_wp
2043       DO  i = nxl, nxr
2044          DO  k = nzb, nzt+1
2045             DO  ii = -n_x(k), n_x(k)
2046                DO  kk = -n_z(k), n_z(k)
2047                   f_n(k,i) = f_n(k,i) + b_x(ii,k) * b_z(kk,k) * rand_it(k+kk,i+ii)
2048                ENDDO
2049             ENDDO
2050          ENDDO
2051       ENDDO
2052
2053    ELSE
2054
2055       f_n_l  = 0.0_wp
2056       DO  i = nxl, nxr
2057          DO  k = nzb_y_stg, nzt_y_stg+1
2058             DO  ii = -n_x(k), n_x(k)
2059                DO  kk = -n_z(k), n_z(k)
2060                   f_n_l(k,i) = f_n_l(k,i) + b_x(ii,k) * b_z(kk,k) * rand_it(k+kk,i+ii)
2061                ENDDO
2062             ENDDO
2063          ENDDO
2064       ENDDO
2065!
2066!--    Gather velocity seeds of full subdomain
2067       send_count = nzt_y_stg - nzb_y_stg + 1
2068       IF ( nzt_y_stg == nzt )  send_count = send_count + 1
2069
2070#if defined( __parallel )
2071!
2072!--    Gather the processed velocity seed on south boundary mpi ranks.
2073       CALL MPI_GATHERV( f_n_l(nzb_y_stg,nxl), send_count, stg_type_xz_small, f_n(nzb+1,nxl),      &
2074                         recv_count_xz, displs_xz, stg_type_xz, id_south, comm1dy, ierr )
2075!
2076!--    Gather the processed velocity seed on north boundary mpi ranks.
2077       CALL MPI_GATHERV( f_n_l(nzb_y_stg,nxl), send_count, stg_type_xz_small, f_n(nzb+1,nxl),      &
2078                         recv_count_xz, displs_xz, stg_type_xz, id_north, comm1dy, ierr )
2079#else
2080       f_n(nzb+1:nzt+1,nxl:nxr) = f_n_l(nzb_y_stg:nzt_y_stg+1,nxl:nxr)
2081!
2082!--    Next line required to avoid compile errors because of unused dummy arguments
2083       IF ( id_north == 0 )  id_south = 0
2084#endif
2085
2086    ENDIF
2087
2088    DEALLOCATE( rand_it )
2089
2090    CALL cpu_log( log_point_s(31), 'STG f_n factors', 'stop' )
2091
2092 END SUBROUTINE stg_generate_seed_xz
2093
2094!--------------------------------------------------------------------------------------------------!
2095! Description:
2096! ------------
2097!> Parametrization of the Reynolds-stress componentes, turbulent length- and time scales. The
2098!> parametrization follows Brost et al. (1982) with modifications described in Rotach et al. (1996),
2099!> which is applied in state-of-the-art dispserion modelling.
2100!--------------------------------------------------------------------------------------------------!
2101  SUBROUTINE parametrize_turbulence
2102
2103    INTEGER(iwp) :: k                 !< loop index in z-direction
2104
2105    REAL(wp) ::  corr_term_uh         !< correction term in parametrization of horizontal variances for unstable stratification
2106    REAL(wp) ::  d_zi                 !< inverse boundary-layer depth
2107    REAL(wp) ::  length_scale_lon     !< length scale in flow direction
2108    REAL(wp) ::  length_scale_lon_zi  !< length scale in flow direction at boundary-layer top
2109    REAL(wp) ::  length_scale_lat     !< length scale in crosswind direction
2110    REAL(wp) ::  length_scale_lat_zi  !< length scale in crosswind direction at boundary-layer top
2111    REAL(wp) ::  length_scale_vert    !< length scale in vertical direction
2112    REAL(wp) ::  length_scale_vert_zi !< length scale in vertical direction at boundary-layer top
2113    REAL(wp) ::  time_scale_zi        !< time scale at boundary-layer top
2114    REAL(wp) ::  r11_zi               !< longitudinal variance at boundary-layer top
2115    REAL(wp) ::  r22_zi               !< crosswind variance at boundary-layer top
2116    REAL(wp) ::  r33_zi               !< vertical variance at boundary-layer top
2117    REAL(wp) ::  r31_zi               !< r31 at boundary-layer top
2118    REAL(wp) ::  r32_zi               !< r32 at boundary-layer top
2119    REAL(wp) ::  rlat1                !< first dummy argument for crosswind compoment of reynolds stress
2120    REAL(wp) ::  rlat2                !< second dummy argument for crosswind compoment of reynolds stress
2121    REAL(wp) ::  rlon1                !< first dummy argument for longitudinal compoment of reynolds stress
2122    REAL(wp) ::  rlon2                !< second dummy argument for longitudinal compoment of reynolds stress
2123    REAL(wp) ::  zzi                  !< ratio of z/zi
2124
2125    REAL(wp), DIMENSION(nzb+1:nzt+1) ::  cos_phi    !< cosine of angle between mean flow direction and x-axis
2126    REAL(wp), DIMENSION(nzb+1:nzt+1) ::  phi        !< angle between mean flow direction and x-axis
2127    REAL(wp), DIMENSION(nzb+1:nzt+1) ::  sin_phi    !< sine of angle between mean flow direction and x-axis
2128
2129!
2130!-- Calculate the boundary-layer top index. Boundary-layer top is calculated by Richardson-bulk
2131!-- criterion according to Heinze et al. (2017).
2132    k_zi = MAX( 1, MINLOC( ABS( zu - zi_ribulk ), DIM = 1 ) )
2133    IF ( zu(k_zi) > zi_ribulk )  k_zi = k_zi - 1
2134    k_zi = MIN( k_zi, nzt )
2135!
2136!-- Calculate angle between flow direction and x-axis.
2137    DO  k = nzb+1, nzt + 1
2138       IF ( u_init(k) /= 0.0_wp )  THEN
2139          phi(k) = ATAN( v_init(k) / u_init(k) )
2140       ELSE
2141          phi(k) = 0.5 * pi
2142       ENDIF
2143!
2144!--    Pre-calculate sine and cosine.
2145       cos_phi(k) = COS( phi(k) )
2146       sin_phi(k) = SIN( phi(k) )
2147    ENDDO
2148!
2149!-- Parametrize Reynolds-stress components. Please note, parametrization is formulated for the
2150!-- stream- and spanwise components, while stream- and spanwise direction do not necessarily
2151!-- coincide with the grid axis. Hence, components are rotated after computation.
2152    d_zi = 1.0_wp / zi_ribulk
2153    DO  k = nzb+1, k_zi
2154
2155       corr_term_uh = MERGE( 0.35_wp * ABS(                                                        &
2156                        - zi_ribulk / ( kappa * scale_l - 10E-4_wp )                               &
2157                                       )**( 2.0_wp / 3.0_wp ),                                     &
2158                          0.0_wp,                                                                  &
2159                          scale_l < -5.0_wp )
2160!
2161!--    Determine normalized height coordinate
2162       zzi = zu(k) * d_zi
2163!
2164!--    Compute longitudinal and crosswind component of reynolds stress. Rotate these components by
2165!--    the wind direction to map onto the u- and v-component. Note, since reynolds stress for
2166!--    variances cannot become negative, take the absolute values.
2167       rlon1 = scale_us**2 * ( corr_term_uh + 5.0_wp - 4.0_wp * zzi )
2168       rlat1 = scale_us**2 * ( corr_term_uh + 2.0_wp - zzi )
2169
2170       r11(k) = ABS(   cos_phi(k) * rlon1 + sin_phi(k) * rlat1 )
2171       r22(k) = ABS( - sin_phi(k) * rlon1 + cos_phi(k) * rlat1 )
2172!
2173!--    w'w'
2174       r33(k) = scale_wm**2 * (                                                                    &
2175                   1.5_wp * zzi**( 2.0_wp / 3.0_wp ) * EXP( -2.0_wp * zzi )                        &
2176                 + ( 1.7_wp - zzi ) * ( scale_us / scale_wm )**2                                   &
2177                              )
2178!
2179!--    u'w' and v'w'. After calculation of the longitudinal and crosswind component
2180!--    these are projected along the x- and y-direction. Note, it is assumed that
2181!--    the flux within the boundary points opposite to the vertical gradient.
2182       rlon2 = scale_us**2 * ( zzi - 1.0_wp )
2183       rlat2 = scale_us**2 * ( 0.4 * zzi * ( 1.0_wp - zzi ) )
2184
2185       r31(k) = SIGN( ABS(   cos_phi(k) * rlon2 + sin_phi(k) * rlat2 ),                            &
2186                   - ( u_init(k+1) - u_init(k-1) ) )
2187       r32(k) = SIGN( ABS( - sin_phi(k) * rlon2 + cos_phi(k) * rlat2 ),                            &
2188                   - ( v_init(k+1) - v_init(k-1) ) )
2189!
2190!--    For u'v' no parametrization exist so far. For simplicity assume a similar profile as
2191!--    for the vertical transport.
2192       r21(k) = 0.5_wp * ( r31(k) + r32(k) )
2193!
2194!--    Compute turbulent time scales according to Brost et al. (1982). Note, time scales are
2195!--    limited to the adjustment time scales.
2196       tu(k)  = MIN( dt_stg_adjust,                                                                &
2197                     3.33_wp * zzi * ( 1.0 - 0.67_wp * zzi ) / scale_wm * zi_ribulk )
2198
2199       tv(k)  = tu(k)
2200       tw(k)  = tu(k)
2201
2202       length_scale_lon  = MIN( SQRT( r11(k) ) * tu(k), zi_ribulk )
2203       length_scale_lat  = MIN( SQRT( r22(k) ) * tv(k), zi_ribulk )
2204       length_scale_vert = MIN( SQRT( r33(k) ) * tw(k), zi_ribulk )
2205!
2206!--    Assume isotropic turbulence length scales
2207       nux(k) = MAX( INT( length_scale_lon  * ddx     ), 1 )
2208       nuy(k) = MAX( INT( length_scale_lat  * ddy     ), 1 )
2209       nvx(k) = MAX( INT( length_scale_lon  * ddx     ), 1 )
2210       nvy(k) = MAX( INT( length_scale_lat  * ddy     ), 1 )
2211       nwx(k) = MAX( INT( length_scale_lon  * ddx     ), 1 )
2212       nwy(k) = MAX( INT( length_scale_lat  * ddy     ), 1 )
2213       nuz(k) = MAX( INT( length_scale_vert * ddzw(k) ), 1 )
2214       nvz(k) = MAX( INT( length_scale_vert * ddzw(k) ), 1 )
2215       nwz(k) = MAX( INT( length_scale_vert * ddzw(k) ), 1 )
2216
2217    ENDDO
2218!
2219!-- Above boundary-layer top length- and timescales as well as reynolds-stress components are
2220!-- reduced to zero by a blending function.
2221    length_scale_lon_zi  = SQRT( r11(k_zi) ) * tu(k_zi)
2222    length_scale_lat_zi  = SQRT( r22(k_zi) ) * tv(k_zi)
2223    length_scale_vert_zi = SQRT( r33(k_zi) ) * tu(k_zi)
2224
2225    time_scale_zi        = tu(k_zi)
2226    r11_zi               = r11(k_zi)
2227    r22_zi               = r22(k_zi)
2228    r33_zi               = r33(k_zi)
2229    r31_zi               = r31(k_zi)
2230    r32_zi               = r32(k_zi)
2231
2232    d_l = blend_coeff / MAX( length_scale_vert_zi, dx, dy, MINVAL( dzw ) )
2233
2234    DO  k = k_zi+1, nzt+1
2235!
2236!--    Calculate function to gradually decrease Reynolds stress above ABL top.
2237       blend = MIN( 1.0_wp, EXP( d_l * zu(k) - d_l * zi_ribulk ) )
2238!
2239!--    u'u' and v'v'. Assume isotropy. Note, add a small negative number to the denominator, else
2240!--    the mergpe-function can crash if scale_l is zero.
2241       r11(k) = r11_zi * blend
2242       r22(k) = r22_zi * blend
2243       r33(k) = r33_zi * blend
2244       r31(k) = r31_zi * blend
2245       r32(k) = r32_zi * blend
2246       r21(k) = 0.5_wp * ( r31(k) + r32(k) )
2247
2248!
2249!--    Compute turbulent time scales according to Brost et al. (1982).
2250!--    Note, time scales are limited to the adjustment time scales.
2251       tu(k)  = time_scale_zi * blend
2252       tv(k)  = time_scale_zi * blend
2253       tw(k)  = time_scale_zi * blend
2254
2255       length_scale_lon  = length_scale_lon_zi  * blend
2256       length_scale_lat  = length_scale_lat_zi  * blend
2257       length_scale_vert = length_scale_vert_zi * blend
2258!
2259!--    Assume isotropic turbulence length scales
2260       nux(k) = MAX( INT( length_scale_lon  * ddx ),     1 )
2261       nuy(k) = MAX( INT( length_scale_lat  * ddy ),     1 )
2262       nvx(k) = MAX( INT( length_scale_lon  * ddx ),     1 )
2263       nvy(k) = MAX( INT( length_scale_lat  * ddy ),     1 )
2264       nwx(k) = MAX( INT( length_scale_lon  * ddx ),     1 )
2265       nwy(k) = MAX( INT( length_scale_lat  * ddy ),     1 )
2266       nuz(k) = MAX( INT( length_scale_vert * ddzw(k) ), 1 )
2267       nvz(k) = MAX( INT( length_scale_vert * ddzw(k) ), 1 )
2268       nwz(k) = MAX( INT( length_scale_vert * ddzw(k) ), 1 )
2269    ENDDO
2270!
2271!-- Set bottom boundary condition for reynolds stresses
2272    r11(nzb) = r11(nzb+1)
2273    r22(nzb) = r22(nzb+1)
2274    r33(nzb) = r33(nzb+1)
2275
2276    r21(nzb) = r11(nzb+1)
2277    r31(nzb) = r31(nzb+1)
2278    r32(nzb) = r32(nzb+1)
2279
2280!
2281!-- Set bottom boundary condition for the length and time scales
2282    nux(nzb) = nux(nzb+1)
2283    nuy(nzb) = nuy(nzb+1)
2284    nuz(nzb) = nuz(nzb+1)
2285    nvx(nzb) = nvx(nzb+1)
2286    nvy(nzb) = nvy(nzb+1)
2287    nvz(nzb) = nvz(nzb+1)
2288    nwx(nzb) = nwx(nzb+1)
2289    nwy(nzb) = nwy(nzb+1)
2290    nwz(nzb) = nwz(nzb+1)
2291
2292    tu(nzb)  = tu(nzb+1)
2293    tv(nzb)  = tv(nzb+1)
2294    tw(nzb)  = tw(nzb+1)
2295
2296    adjustment_step = .TRUE.
2297
2298 END SUBROUTINE parametrize_turbulence
2299
2300!--------------------------------------------------------------------------------------------------!
2301! Description:
2302! ------------
2303!> Calculate the coefficient matrix from the Lund rotation.
2304!--------------------------------------------------------------------------------------------------!
2305 SUBROUTINE calc_coeff_matrix
2306
2307    INTEGER(iwp) ::  k  !< loop index in z-direction
2308
2309!
2310!-- Calculate coefficient matrix. Split loops to allow for loop vectorization.
2311    DO  k = nzb+1, nzt+1
2312       IF ( r11(k) > 10E-6_wp )  THEN
2313          a11(k) = SQRT( r11(k) )
2314          a21(k) = r21(k) / a11(k)
2315          a31(k) = r31(k) / a11(k)
2316       ELSE
2317          a11(k) = 10E-8_wp
2318          a21(k) = 10E-8_wp
2319          a31(k) = 10E-8_wp
2320       ENDIF
2321    ENDDO
2322    DO  k = nzb+1, nzt+1
2323       a22(k) = r22(k) - a21(k)**2
2324       IF ( a22(k) > 10E-6_wp )  THEN
2325          a22(k) = SQRT( a22(k) )
2326          a32(k) = ( r32(k) - a21(k) * a31(k) ) / a22(k)
2327       ELSE
2328          a22(k) = 10E-8_wp
2329          a32(k) = 10E-8_wp
2330       ENDIF
2331    ENDDO
2332    DO  k = nzb+1, nzt+1
2333       a33(k) = r33(k) - a31(k)**2 - a32(k)**2
2334       IF ( a33(k) > 10E-6_wp )  THEN
2335          a33(k) =  SQRT( a33(k) )
2336       ELSE
2337          a33(k) = 10E-8_wp
2338       ENDIF
2339    ENDDO
2340!
2341!-- Set bottom boundary condition
2342    a11(nzb) = a11(nzb+1)
2343    a22(nzb) = a22(nzb+1)
2344    a21(nzb) = a21(nzb+1)
2345    a33(nzb) = a33(nzb+1)
2346    a31(nzb) = a31(nzb+1)
2347    a32(nzb) = a32(nzb+1)
2348
2349 END SUBROUTINE calc_coeff_matrix
2350
2351!--------------------------------------------------------------------------------------------------!
2352! Description:
2353! ------------
2354!> This routine controls the re-adjustment of the turbulence statistics used for generating
2355!> turbulence at the lateral boundaries.
2356!--------------------------------------------------------------------------------------------------!
2357 SUBROUTINE stg_adjust
2358
2359    IF ( debug_output_timestep )  CALL debug_message( 'stg_adjust', 'start' )
2360!
2361!-- In case of dirichlet inflow boundary conditions only at one lateral boundary, i.e. in the case
2362!-- no offline or self nesting is applied but synthetic turbulence shall be parametrized
2363!-- nevertheless, the boundary-layer depth need to determined first.
2364    IF ( .NOT. nesting_offline  .AND.  .NOT.  child_domain )                                       &
2365       CALL nesting_offl_calc_zi
2366!
2367!-- Compute scaling parameters (domain-averaged), such as friction velocity are calculated.
2368    CALL calc_scaling_variables
2369!
2370!-- Parametrize Reynolds-stress tensor. Parametrization follows Brost et al. (1982) with
2371!-- modifications described in Rotach et al. (1996) and is based on boundary-layer depth, friction
2372!-- velocity and velocity scale.
2373    CALL parametrize_turbulence
2374!
2375!-- Calculate coefficient matrix from Reynolds stress tensor
2376!-- (Lund rotation)
2377    CALL calc_coeff_matrix
2378!
2379!-- Setup filter functions according to updated length scales.
2380    CALL stg_setup_filter_function
2381!
2382!-- Reset time counter for controlling next adjustment to zero
2383    time_stg_adjust = 0.0_wp
2384
2385    IF ( debug_output_timestep )  CALL debug_message( 'stg_adjust', 'end' )
2386
2387 END SUBROUTINE stg_adjust
2388
2389
2390!--------------------------------------------------------------------------------------------------!
2391! Description:
2392! ------------
2393!> Calculate scaling variables which are used for turbulence parametrization according to Rotach
2394!> et al. (1996). Scaling variables are: friction velocity, boundary-layer depth, momentum velocity
2395!> scale, and Obukhov length.
2396!--------------------------------------------------------------------------------------------------!
2397 SUBROUTINE calc_scaling_variables
2398
2399    INTEGER(iwp) ::  i  !< loop index in x-direction
2400    INTEGER(iwp) ::  j  !< loop index in y-direction
2401    INTEGER(iwp) ::  k  !< loop index in z-direction
2402    INTEGER(iwp) ::  m  !< surface element index
2403
2404    REAL(wp) ::  friction_vel_l  !< mean friction veloctiy on subdomain
2405    REAL(wp) ::  pt_surf_mean    !< mean near surface temperature (at 1st grid point)
2406    REAL(wp) ::  pt_surf_mean_l  !< mean near surface temperature (at 1st grid point) on subdomain
2407    REAL(wp) ::  scale_l_l       !< mean Obukhov lenght on subdomain
2408    REAL(wp) ::  shf_mean        !< mean surface sensible heat flux
2409    REAL(wp) ::  shf_mean_l      !< mean surface sensible heat flux on subdomain
2410    REAL(wp) ::  w_convective    !< convective velocity scale
2411
2412!
2413!-- Calculate mean friction velocity, velocity scale, heat flux and near-surface temperature in the
2414!-- model domain.
2415    pt_surf_mean_l = 0.0_wp
2416    shf_mean_l     = 0.0_wp
2417    scale_l_l      = 0.0_wp
2418    friction_vel_l = 0.0_wp
2419    DO  m = 1, surf_def_h(0)%ns
2420       i = surf_def_h(0)%i(m)
2421       j = surf_def_h(0)%j(m)
2422       k = surf_def_h(0)%k(m)
2423       friction_vel_l = friction_vel_l  + surf_def_h(0)%us(m)
2424       shf_mean_l     = shf_mean_l      + surf_def_h(0)%shf(m) * drho_air(k)
2425       scale_l_l      = scale_l_l       + surf_def_h(0)%ol(m)
2426       pt_surf_mean_l = pt_surf_mean_l  + pt(k,j,i)
2427    ENDDO
2428    DO  m = 1, surf_lsm_h(0)%ns
2429       i = surf_lsm_h(0)%i(m)
2430       j = surf_lsm_h(0)%j(m)
2431       k = surf_lsm_h(0)%k(m)
2432       friction_vel_l = friction_vel_l  + surf_lsm_h(0)%us(m)
2433       shf_mean_l     = shf_mean_l      + surf_lsm_h(0)%shf(m) * drho_air(k)
2434       scale_l_l      = scale_l_l       + surf_lsm_h(0)%ol(m)
2435       pt_surf_mean_l = pt_surf_mean_l  + pt(k,j,i)
2436    ENDDO
2437    DO  m = 1, surf_usm_h(0)%ns
2438       i = surf_usm_h(0)%i(m)
2439       j = surf_usm_h(0)%j(m)
2440       k = surf_usm_h(0)%k(m)
2441       friction_vel_l = friction_vel_l  + surf_usm_h(0)%us(m)
2442       shf_mean_l     = shf_mean_l      + surf_usm_h(0)%shf(m) * drho_air(k)
2443       scale_l_l      = scale_l_l       + surf_usm_h(0)%ol(m)
2444       pt_surf_mean_l = pt_surf_mean_l  + pt(k,j,i)
2445    ENDDO
2446
2447#if defined( __parallel )
2448    CALL MPI_ALLREDUCE( friction_vel_l, scale_us,     1, MPI_REAL, MPI_SUM, comm2d, ierr )
2449    CALL MPI_ALLREDUCE( shf_mean_l, shf_mean,         1, MPI_REAL, MPI_SUM, comm2d, ierr )
2450    CALL MPI_ALLREDUCE( scale_l_l, scale_l,           1, MPI_REAL, MPI_SUM, comm2d, ierr )
2451    CALL MPI_ALLREDUCE( pt_surf_mean_l, pt_surf_mean, 1, MPI_REAL, MPI_SUM, comm2d, ierr )
2452#else
2453    scale_us     = friction_vel_l
2454    shf_mean     = shf_mean_l
2455    scale_l      = scale_l_l
2456    pt_surf_mean = pt_surf_mean_l
2457#endif
2458
2459    scale_us     = scale_us     * d_nxy
2460    shf_mean     = shf_mean     * d_nxy
2461    scale_l      = scale_l      * d_nxy
2462    pt_surf_mean = pt_surf_mean * d_nxy
2463!
2464!-- Compute mean convective velocity scale. Note, in case the mean heat flux is negative, set
2465!-- convective velocity scale to zero.
2466    IF ( shf_mean > 0.0_wp )  THEN
2467       w_convective = ( g * shf_mean * zi_ribulk / pt_surf_mean )**( 1.0_wp / 3.0_wp )
2468    ELSE
2469       w_convective = 0.0_wp
2470    ENDIF
2471!
2472!-- At the initial run the friction velocity is initialized close to zero, leading to almost zero
2473!-- disturbances at the boundaries. In order to obtain disturbances nevertheless, set a minium
2474!-- value of friction velocity for the reynolds-stress parametrization.
2475    IF ( time_since_reference_point <= 0.0_wp )  scale_us = MAX( scale_us, 0.05_wp )
2476!
2477!-- Finally, in order to consider also neutral or stable stratification, compute momentum velocity
2478!-- scale from u* and convective velocity scale, according to Rotach et al. (1996).
2479    scale_wm = ( scale_us**3 + 0.6_wp * w_convective**3 )**( 1.0_wp / 3.0_wp )
2480
2481 END SUBROUTINE calc_scaling_variables
2482
2483 END MODULE synthetic_turbulence_generator_mod
Note: See TracBrowser for help on using the repository browser.