!> @synthetic_turbulence_generator_mod.f90 !------------------------------------------------------------------------------! ! This file is part of the PALM model system. ! ! PALM is free software: you can redistribute it and/or modify it under the ! terms of the GNU General Public License as published by the Free Software ! Foundation, either version 3 of the License, or (at your option) any later ! version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 2017 Leibniz Universitaet Hannover !------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: synthetic_turbulence_generator_mod.f90 3775 2019-03-04 12:40:20Z eckhard $ ! removed unused variables ! ! 3719 2019-02-06 13:10:18Z kanani ! Removed log_point measurement from stg_init, since this part is counted to ! log_point(2) 'initialisation' already. Moved other log_points to calls of ! the subroutines in time_integration for better overview. ! ! 3646 2018-12-28 17:58:49Z kanani ! Bugfix: use time_since_reference_point instead of simulated_time (relevant ! when using wall/soil spinup) ! ! 3579 2018-11-29 15:32:39Z suehring ! - Bugfix in calculation of turbulence scaling parameters for turbulence ! generator in case of topography ! - Additional checks implemented - no STG in RANS-RANS nesting or LES-LES ! nesting ! ! 3376 2018-10-19 10:15:32Z suehring ! Error messages and numbers reivsed. ! ! 3349 2018-10-15 16:39:41Z suehring ! Fix for format descriptor ! ! 3348 2018-10-15 14:30:51Z suehring ! - Revise structure of init routine ! - introduce new parameters to skip STG for some timesteps ! - introduce time-dependent parametrization of Reynolds-stress tensor ! - Bugfix in allocation of mean_inflow_profiles ! ! 3341 2018-10-15 10:31:27Z suehring ! Introduce module parameter for number of inflow profiles ! ! 3288 2018-09-28 10:23:08Z suehring ! Modularization of all bulk cloud physics code components ! ! 3248 2018-09-14 09:42:06Z sward ! Minor formating changes ! ! 3246 2018-09-13 15:14:50Z sward ! Added error handling for input namelist via parin_fail_message ! ! 3241 2018-09-12 15:02:00Z raasch ! unused variables removed ! ! 3186 2018-07-30 17:07:14Z suehring ! Mask topography while imposing inflow perturbations at the boundaries; do not ! impose perturbations at top boundary as well as ghost points ! ! 3183 2018-07-27 14:25:55Z suehring ! Rename variables and extend error message ! Enable geneartor also for stretched grids ! ! 3182 2018-07-27 13:36:03Z suehring ! Error message related to vertical stretching has been added, dz was replaced ! by dz(1) ! ! 3051 2018-05-30 17:43:55Z suehring ! Bugfix in calculation of initial Reynolds-stress tensor. ! ! 3045 2018-05-28 07:55:41Z Giersch ! Error messages revised ! ! 3044 2018-05-25 10:59:41Z gronemeier ! Add missing variable descriptions ! ! 3038 2018-05-24 10:54:00Z gronemeier ! updated variable descriptions ! ! 2967 2018-04-13 11:22:08Z raasch ! bugfix: missing parallel cpp-directives added ! ! 2946 2018-04-04 17:01:23Z suehring ! Remove unused module load ! ! 2945 2018-04-04 16:27:14Z suehring ! - Bugfix in parallelization of synthetic turbulence generator in case the ! z-dimension is not an integral divisor of the number of processors along ! the x- and y-dimension ! - Revision in control mimic in case of RAN-LES nesting ! ! 2938 2018-03-27 15:52:42Z suehring ! Apply turbulence generator at all non-cyclic lateral boundaries in case of ! realistic Inifor large-scale forcing or RANS-LES nesting ! ! 2936 2018-03-27 14:49:27Z suehring ! variable named found has been introduced for checking if restart data was found, ! reading of restart strings has been moved completely to read_restart_data_mod, ! redundant skipping function has been removed, stg_read/write_restart_data ! have been renamed to stg_r/wrd_global, stg_rrd_global is called in ! read_restart_data_mod now, flag syn_turb_gen_prerun and marker *** end stg ! *** have been removed (Giersch), strings and their respective lengths are ! written out and read now in case of restart runs to get rid of prescribed ! character lengths (Giersch), CASE DEFAULT was added if restart data is read ! ! 2841 2018-02-27 15:02:57Z suehring ! Bugfix: wrong placement of include 'mpif.h' corrected ! ! 2836 2018-02-26 13:40:05Z Giersch ! The variables synthetic_turbulence_generator and ! use_synthetic_turbulence_generator have been abbreviated + syn_turb_gen_prerun ! flag is used to define if module related parameters were outputted as restart ! data ! ! 2716 2017-12-29 16:35:59Z kanani ! Corrected "Former revisions" section ! ! 2696 2017-12-14 17:12:51Z kanani ! Change in file header (GPL part) ! ! 2669 2017-12-06 16:03:27Z raasch ! unit number for file containing turbulence generator data changed to 90 ! bugfix: preprocessor directives added for MPI specific code ! ! 2576 2017-10-24 13:49:46Z Giersch ! Definition of a new function called stg_skip_global to skip module ! parameters during reading restart data ! ! 2563 2017-10-19 15:36:10Z Giersch ! stg_read_restart_data is called in stg_parin in the case of a restart run ! ! 2259 2017-06-08 09:09:11Z gronemeier ! Initial revision ! ! ! ! Authors: ! -------- ! @author Tobias Gronemeier, Matthias Suehring, Atsushi Inagaki, Micha Gryschka, Christoph Knigge ! ! ! ! Description: ! ------------ !> The module generates turbulence at the inflow boundary based on a method by !> Xie and Castro (2008) utilizing a Lund rotation (Lund, 1998) and a mass-flux !> correction by Kim et al. (2013). !> The turbulence is correlated based on length scales in y- and z-direction and !> a time scale for each velocity component. The profiles of length and time !> scales, mean u, v, w, e and pt, and all components of the Reynolds stress !> tensor can be either read from file STG_PROFILES, or will be parametrized !> within the boundary layer. !> !> @todo test restart !> enable cyclic_fill !> implement turbulence generation for e and pt !> @todo Input of height-constant length scales via namelist !> @note !> @bug Height information from input file is not used. Profiles from input !> must match with current PALM grid. !> In case of restart, velocity seeds differ from precursor run if a11, !> a22, or a33 are zero. !------------------------------------------------------------------------------! MODULE synthetic_turbulence_generator_mod USE arrays_3d, & ONLY: mean_inflow_profiles, q, pt, u, v, w, zu, zw USE basic_constants_and_equations_mod, & ONLY: g, kappa, pi USE control_parameters, & ONLY: initializing_actions, num_mean_inflow_profiles, message_string, & syn_turb_gen USE indices, & ONLY: nbgp, nzb, nzt, nxl, nxlg, nxr, nxrg, nys, nyn, nyng, nysg USE kinds #if defined( __parallel ) && !defined( __mpifh ) USE MPI #endif USE nesting_offl_mod, & ONLY: zi_ribulk USE pegrid, & ONLY: comm1dx, comm1dy, comm2d, ierr, myidx, myidy, pdims USE transpose_indices, & ONLY: nzb_x, nzt_x IMPLICIT NONE #if defined( __parallel ) && defined( __mpifh ) INCLUDE "mpif.h" #endif LOGICAL :: velocity_seed_initialized = .FALSE. !< true after first call of stg_main LOGICAL :: parametrize_inflow_turbulence = .FALSE. !< flag indicating that inflow turbulence is either read from file (.FALSE.) or if it parametrized LOGICAL :: use_syn_turb_gen = .FALSE. !< switch to use synthetic turbulence generator INTEGER(iwp) :: id_stg_left !< left lateral boundary core id in case of turbulence generator INTEGER(iwp) :: id_stg_north !< north lateral boundary core id in case of turbulence generator INTEGER(iwp) :: id_stg_right !< right lateral boundary core id in case of turbulence generator INTEGER(iwp) :: id_stg_south !< south lateral boundary core id in case of turbulence generator INTEGER(iwp) :: stg_type_xz !< MPI type for full z range INTEGER(iwp) :: stg_type_xz_small !< MPI type for small z range INTEGER(iwp) :: stg_type_yz !< MPI type for full z range INTEGER(iwp) :: stg_type_yz_small !< MPI type for small z range INTEGER(iwp) :: merg !< maximum length scale (in gp) INTEGER(iwp) :: mergp !< merg + nbgp INTEGER(iwp) :: nzb_x_stg !< lower bound of z coordinate (required for transposing z on PEs along x) INTEGER(iwp) :: nzt_x_stg !< upper bound of z coordinate (required for transposing z on PEs along x) INTEGER(iwp) :: nzb_y_stg !< lower bound of z coordinate (required for transposing z on PEs along y) INTEGER(iwp) :: nzt_y_stg !< upper bound of z coordinate (required for transposing z on PEs along y) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: displs_xz !< displacement for MPI_GATHERV INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: recv_count_xz !< receive count for MPI_GATHERV INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: displs_yz !< displacement for MPI_GATHERV INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: recv_count_yz !< receive count for MPI_GATHERV INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nux !< length scale of u in x direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nuy !< length scale of u in y direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nuz !< length scale of u in z direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nvx !< length scale of v in x direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nvy !< length scale of v in y direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nvz !< length scale of v in z direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nwx !< length scale of w in x direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nwy !< length scale of w in y direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: nwz !< length scale of w in z direction (in gp) INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: seed !< seed of random number for rn-generator REAL(wp) :: cosmo_ref = 10.0_wp !< height of first vertical grid level in mesoscale model, used for calculation of scaling parameters REAL(wp) :: dt_stg_adjust = 300.0_wp !< time interval for adjusting turbulence statistics REAL(wp) :: dt_stg_call = 5.0_wp !< time interval for calling synthetic turbulence generator REAL(wp) :: mc_factor !< mass flux correction factor REAL(wp) :: scale_l !< scaling parameter used for turbulence parametrization - Obukhov length REAL(wp) :: scale_us !< scaling parameter used for turbulence parametrization - friction velocity REAL(wp) :: scale_wm !< scaling parameter used for turbulence parametrization - momentum scale REAL(wp) :: time_stg_adjust = 0.0_wp !< time counter for adjusting turbulence information REAL(wp) :: time_stg_call = 0.0_wp !< time counter for calling generator REAL(wp),DIMENSION(:), ALLOCATABLE :: r11 !< Reynolds parameter REAL(wp),DIMENSION(:), ALLOCATABLE :: r21 !< Reynolds parameter REAL(wp),DIMENSION(:), ALLOCATABLE :: r22 !< Reynolds parameter REAL(wp),DIMENSION(:), ALLOCATABLE :: r31 !< Reynolds parameter REAL(wp),DIMENSION(:), ALLOCATABLE :: r32 !< Reynolds parameter REAL(wp),DIMENSION(:), ALLOCATABLE :: r33 !< Reynolds parameter REAL(wp), DIMENSION(:), ALLOCATABLE :: a11 !< coefficient for Lund rotation REAL(wp), DIMENSION(:), ALLOCATABLE :: a21 !< coefficient for Lund rotation REAL(wp), DIMENSION(:), ALLOCATABLE :: a22 !< coefficient for Lund rotation REAL(wp), DIMENSION(:), ALLOCATABLE :: a31 !< coefficient for Lund rotation REAL(wp), DIMENSION(:), ALLOCATABLE :: a32 !< coefficient for Lund rotation REAL(wp), DIMENSION(:), ALLOCATABLE :: a33 !< coefficient for Lund rotation REAL(wp), DIMENSION(:), ALLOCATABLE :: tu !< Lagrangian time scale of u REAL(wp), DIMENSION(:), ALLOCATABLE :: tv !< Lagrangian time scale of v REAL(wp), DIMENSION(:), ALLOCATABLE :: tw !< Lagrangian time scale of w REAL(wp), DIMENSION(:,:), ALLOCATABLE :: bux !< filter function for u in x direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: buy !< filter function for u in y direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: buz !< filter function for u in z direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: bvx !< filter function for v in x direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: bvy !< filter function for v in y direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: bvz !< filter function for v in z direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: bwx !< filter function for w in y direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: bwy !< filter function for w in y direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: bwz !< filter function for w in z direction REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fu_xz !< velocity seed for u at xz plane REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fuo_xz !< velocity seed for u at xz plane with new random number REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fu_yz !< velocity seed for u at yz plane REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fuo_yz !< velocity seed for u at yz plane with new random number REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fv_xz !< velocity seed for v at xz plane REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fvo_xz !< velocity seed for v at xz plane with new random number REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fv_yz !< velocity seed for v at yz plane REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fvo_yz !< velocity seed for v at yz plane with new random number REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fw_xz !< velocity seed for w at xz plane REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fwo_xz !< velocity seed for w at xz plane with new random number REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fw_yz !< velocity seed for w at yz plane REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fwo_yz !< velocity seed for w at yz plane with new random number REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: dist_xz !< imposed disturbances at north/south boundary REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: dist_yz !< imposed disturbances at north/south boundary ! !-- PALM interfaces: !-- Adjust time and lenght scales, Reynolds stress, and filter functions INTERFACE stg_adjust MODULE PROCEDURE stg_adjust END INTERFACE stg_adjust ! !-- Input parameter checks to be done in check_parameters INTERFACE stg_check_parameters MODULE PROCEDURE stg_check_parameters END INTERFACE stg_check_parameters ! !-- Calculate filter functions INTERFACE stg_filter_func MODULE PROCEDURE stg_filter_func END INTERFACE stg_filter_func ! !-- Generate velocity seeds at south and north domain boundary INTERFACE stg_generate_seed_xz MODULE PROCEDURE stg_generate_seed_xz END INTERFACE stg_generate_seed_xz ! !-- Generate velocity seeds at left and/or right domain boundary INTERFACE stg_generate_seed_yz MODULE PROCEDURE stg_generate_seed_yz END INTERFACE stg_generate_seed_yz ! !-- Output of information to the header file INTERFACE stg_header MODULE PROCEDURE stg_header END INTERFACE stg_header ! !-- Initialization actions INTERFACE stg_init MODULE PROCEDURE stg_init END INTERFACE stg_init ! !-- Main procedure of synth. turb. gen. INTERFACE stg_main MODULE PROCEDURE stg_main END INTERFACE stg_main ! !-- Reading of NAMELIST parameters INTERFACE stg_parin MODULE PROCEDURE stg_parin END INTERFACE stg_parin ! !-- Reading of parameters for restart runs INTERFACE stg_rrd_global MODULE PROCEDURE stg_rrd_global END INTERFACE stg_rrd_global ! !-- Writing of binary output for restart runs INTERFACE stg_wrd_global MODULE PROCEDURE stg_wrd_global END INTERFACE stg_wrd_global SAVE PRIVATE ! !-- Public interfaces PUBLIC stg_adjust, stg_check_parameters, stg_header, stg_init, stg_main, & stg_parin, stg_rrd_global, stg_wrd_global ! !-- Public variables PUBLIC dt_stg_call, dt_stg_adjust, id_stg_left, id_stg_north, & id_stg_right, id_stg_south, parametrize_inflow_turbulence, & time_stg_adjust, time_stg_call, use_syn_turb_gen CONTAINS !------------------------------------------------------------------------------! ! Description: ! ------------ !> Check parameters routine for synthetic turbulence generator !------------------------------------------------------------------------------! SUBROUTINE stg_check_parameters USE control_parameters, & ONLY: bc_lr, bc_ns, child_domain, nesting_offline, rans_mode, & turbulent_inflow USE pmc_interface, & ONLY : rans_mode_parent IMPLICIT NONE IF ( .NOT. use_syn_turb_gen .AND. .NOT. rans_mode .AND. & nesting_offline ) THEN message_string = 'Synthetic turbulence generator is required ' // & 'if offline nesting is applied and PALM operates ' // & 'in LES mode.' CALL message( 'stg_check_parameters', 'PA0520', 0, 0, 0, 6, 0 ) ENDIF IF ( .NOT. use_syn_turb_gen .AND. child_domain & .AND. rans_mode_parent .AND. .NOT. rans_mode ) THEN message_string = 'Synthetic turbulence generator is required ' // & 'when nesting is applied and parent operates in ' // & 'RANS-mode but current child in LES mode.' CALL message( 'stg_check_parameters', 'PA0524', 1, 2, 0, 6, 0 ) ENDIF IF ( use_syn_turb_gen ) THEN IF ( child_domain .AND. .NOT. rans_mode .AND. & .NOT. rans_mode_parent ) THEN message_string = 'Using synthetic turbulence generator ' // & 'is not allowed in LES-LES nesting.' CALL message( 'stg_check_parameters', 'PA0620', 1, 2, 0, 6, 0 ) ENDIF IF ( child_domain .AND. rans_mode .AND. & rans_mode_parent ) THEN message_string = 'Using synthetic turbulence generator ' // & 'is not allowed in RANS-RANS nesting.' CALL message( 'stg_check_parameters', 'PA0621', 1, 2, 0, 6, 0 ) ENDIF IF ( .NOT. nesting_offline .AND. .NOT. child_domain ) THEN IF ( INDEX( initializing_actions, 'set_constant_profiles' ) == 0 & .AND. INDEX( initializing_actions, 'read_restart_data' ) == 0 ) THEN message_string = 'Using synthetic turbulence generator ' // & 'requires %initializing_actions = ' // & '"set_constant_profiles" or "read_restart_data"' //& ', if not offline nesting is applied.' CALL message( 'stg_check_parameters', 'PA0015', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_lr /= 'dirichlet/radiation' ) THEN message_string = 'Using synthetic turbulence generator ' // & 'requires &bc_lr = "dirichlet/radiation", ' // & 'if not offline nesting is applied.' CALL message( 'stg_check_parameters', 'PA0035', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_ns /= 'cyclic' ) THEN message_string = 'Using synthetic turbulence generator ' // & 'requires &bc_ns = "cyclic", ' // & 'if not offline nesting is applied.' CALL message( 'stg_check_parameters', 'PA0037', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( turbulent_inflow ) THEN message_string = 'Using synthetic turbulence generator ' // & 'in combination &with turbulent_inflow = .T. '// & 'is not allowed' CALL message( 'stg_check_parameters', 'PA0039', 1, 2, 0, 6, 0 ) ENDIF ENDIF END SUBROUTINE stg_check_parameters !------------------------------------------------------------------------------! ! Description: ! ------------ !> Header output for synthetic turbulence generator !------------------------------------------------------------------------------! SUBROUTINE stg_header ( io ) IMPLICIT NONE INTEGER(iwp), INTENT(IN) :: io !< Unit of the output file ! !-- Write synthetic turbulence generator Header WRITE( io, 1 ) IF ( use_syn_turb_gen ) THEN WRITE( io, 2 ) ELSE WRITE( io, 3 ) ENDIF IF ( parametrize_inflow_turbulence ) THEN WRITE( io, 4 ) dt_stg_adjust ELSE WRITE( io, 5 ) ENDIF 1 FORMAT (//' Synthetic turbulence generator information:'/ & ' ------------------------------------------'/) 2 FORMAT (' synthetic turbulence generator is switched on') 3 FORMAT (' synthetic turbulence generator is switched off') 4 FORMAT (' imposed turbulence statistics are parametrized and ajdusted to boundary-layer development each ', F8.2, ' s' ) 5 FORMAT (' imposed turbulence is read from file' ) END SUBROUTINE stg_header !------------------------------------------------------------------------------! ! Description: ! ------------ !> Initialization of the synthetic turbulence generator !------------------------------------------------------------------------------! SUBROUTINE stg_init USE arrays_3d, & ONLY: ddzw, u_init, v_init USE control_parameters, & ONLY: child_domain, coupling_char, e_init, nesting_offline, rans_mode USE grid_variables, & ONLY: ddy USE indices, & ONLY: nz USE pmc_interface, & ONLY : rans_mode_parent IMPLICIT NONE LOGICAL :: file_stg_exist = .FALSE. !< flag indicating whether parameter file for Reynolds stress and length scales exist #if defined( __parallel ) INTEGER(KIND=MPI_ADDRESS_KIND) :: extent !< extent of new MPI type INTEGER(KIND=MPI_ADDRESS_KIND) :: tob=0 !< dummy variable #endif INTEGER(iwp) :: i !> grid index in x-direction INTEGER(iwp) :: j !> loop index INTEGER(iwp) :: k !< index INTEGER(iwp) :: newtype !< dummy MPI type INTEGER(iwp) :: realsize !< size of REAL variables INTEGER(iwp) :: nseed !< dimension of random number seed INTEGER(iwp) :: startseed = 1234567890 !< start seed for random number generator ! !-- Dummy variables used for reading profiles REAL(wp) :: d1 !< u profile REAL(wp) :: d2 !< v profile REAL(wp) :: d3 !< w profile REAL(wp) :: d5 !< e profile REAL(wp) :: luy !< length scale for u in y direction REAL(wp) :: luz !< length scale for u in z direction REAL(wp) :: lvy !< length scale for v in y direction REAL(wp) :: lvz !< length scale for v in z direction REAL(wp) :: lwy !< length scale for w in y direction REAL(wp) :: lwz !< length scale for w in z direction REAL(wp) :: nnz !< increment used to determine processor decomposition of z-axis along x and y direction REAL(wp) :: zz !< height #if defined( __parallel ) CALL MPI_BARRIER( comm2d, ierr ) #endif #if defined( __parallel ) ! !-- Determine processor decomposition of z-axis along x- and y-direction nnz = nz / pdims(1) nzb_x_stg = 1 + myidx * INT( nnz ) nzt_x_stg = ( myidx + 1 ) * INT( nnz ) IF ( MOD( nz , pdims(1) ) /= 0 .AND. myidx == id_stg_right ) & nzt_x_stg = nzt_x_stg + myidx * ( nnz - INT( nnz ) ) IF ( nesting_offline .OR. ( child_domain .AND. rans_mode_parent & .AND. .NOT. rans_mode ) ) THEN nnz = nz / pdims(2) nzb_y_stg = 1 + myidy * INT( nnz ) nzt_y_stg = ( myidy + 1 ) * INT( nnz ) IF ( MOD( nz , pdims(2) ) /= 0 .AND. myidy == id_stg_north ) & nzt_y_stg = nzt_y_stg + myidy * ( nnz - INT( nnz ) ) ENDIF ! !-- Define MPI type used in stg_generate_seed_yz to gather vertical splitted !-- velocity seeds CALL MPI_TYPE_SIZE( MPI_REAL, realsize, ierr ) extent = 1 * realsize ! !-- Set-up MPI datatyp to involve all cores for turbulence generation at yz !-- layer !-- stg_type_yz: yz-slice with vertical bounds nzb:nzt+1 CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt-nzb+2,nyng-nysg+1], & [1,nyng-nysg+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr ) CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_yz, ierr ) CALL MPI_TYPE_COMMIT( stg_type_yz, ierr ) CALL MPI_TYPE_FREE( newtype, ierr ) ! stg_type_yz_small: yz-slice with vertical bounds nzb_x_stg:nzt_x_stg+1 CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt_x_stg-nzb_x_stg+2,nyng-nysg+1], & [1,nyng-nysg+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr ) CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_yz_small, ierr ) CALL MPI_TYPE_COMMIT( stg_type_yz_small, ierr ) CALL MPI_TYPE_FREE( newtype, ierr ) ! receive count and displacement for MPI_GATHERV in stg_generate_seed_yz ALLOCATE( recv_count_yz(pdims(1)), displs_yz(pdims(1)) ) recv_count_yz = nzt_x_stg-nzb_x_stg + 1 recv_count_yz(pdims(1)) = recv_count_yz(pdims(1)) + 1 DO j = 1, pdims(1) displs_yz(j) = 0 + (nzt_x_stg-nzb_x_stg+1) * (j-1) ENDDO ! !-- Set-up MPI datatyp to involve all cores for turbulence generation at xz !-- layer !-- stg_type_xz: xz-slice with vertical bounds nzb:nzt+1 IF ( nesting_offline .OR. ( child_domain .AND. rans_mode_parent & .AND. .NOT. rans_mode ) ) THEN CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt-nzb+2,nxrg-nxlg+1], & [1,nxrg-nxlg+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr ) CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_xz, ierr ) CALL MPI_TYPE_COMMIT( stg_type_xz, ierr ) CALL MPI_TYPE_FREE( newtype, ierr ) ! stg_type_yz_small: xz-slice with vertical bounds nzb_x_stg:nzt_x_stg+1 CALL MPI_TYPE_CREATE_SUBARRAY( 2, [nzt_y_stg-nzb_y_stg+2,nxrg-nxlg+1], & [1,nxrg-nxlg+1], [0,0], MPI_ORDER_FORTRAN, MPI_REAL, newtype, ierr ) CALL MPI_TYPE_CREATE_RESIZED( newtype, tob, extent, stg_type_xz_small, ierr ) CALL MPI_TYPE_COMMIT( stg_type_xz_small, ierr ) CALL MPI_TYPE_FREE( newtype, ierr ) ! receive count and displacement for MPI_GATHERV in stg_generate_seed_yz ALLOCATE( recv_count_xz(pdims(2)), displs_xz(pdims(2)) ) recv_count_xz = nzt_y_stg-nzb_y_stg + 1 recv_count_xz(pdims(2)) = recv_count_xz(pdims(2)) + 1 DO j = 1, pdims(2) displs_xz(j) = 0 + (nzt_y_stg-nzb_y_stg+1) * (j-1) ENDDO ENDIF #endif ! !-- Define seed of random number CALL RANDOM_SEED() CALL RANDOM_SEED( size=nseed ) ALLOCATE( seed(1:nseed) ) DO j = 1, nseed seed(j) = startseed + j ENDDO CALL RANDOM_SEED(put=seed) ! !-- Allocate required arrays !-- mean_inflow profiles must not be allocated in offline nesting IF ( .NOT. nesting_offline .AND. .NOT. child_domain ) THEN IF ( .NOT. ALLOCATED( mean_inflow_profiles ) ) & ALLOCATE( mean_inflow_profiles(nzb:nzt+1,1:num_mean_inflow_profiles) ) ENDIF ALLOCATE ( a11(nzb:nzt+1), a21(nzb:nzt+1), a22(nzb:nzt+1), & a31(nzb:nzt+1), a32(nzb:nzt+1), a33(nzb:nzt+1), & nux(nzb:nzt+1), nuy(nzb:nzt+1), nuz(nzb:nzt+1), & nvx(nzb:nzt+1), nvy(nzb:nzt+1), nvz(nzb:nzt+1), & nwx(nzb:nzt+1), nwy(nzb:nzt+1), nwz(nzb:nzt+1), & r11(nzb:nzt+1), r21(nzb:nzt+1), r22(nzb:nzt+1), & r31(nzb:nzt+1), r32(nzb:nzt+1), r33(nzb:nzt+1), & tu(nzb:nzt+1), tv(nzb:nzt+1), tw(nzb:nzt+1) ) ALLOCATE ( dist_xz(nzb:nzt+1,nxlg:nxrg,3) ) ALLOCATE ( dist_yz(nzb:nzt+1,nysg:nyng,3) ) dist_xz = 0.0_wp dist_yz = 0.0_wp ! !-- Read inflow profile !-- Height levels of profiles in input profile is as follows: !-- zu: luy, luz, tu, lvy, lvz, tv, r11, r21, r22, d1, d2, d5 !-- zw: lwy, lwz, tw, r31, r32, r33, d3 !-- WARNING: zz is not used at the moment INQUIRE( FILE = 'STG_PROFILES' // TRIM( coupling_char ), & EXIST = file_stg_exist ) IF ( file_stg_exist ) THEN OPEN( 90, FILE='STG_PROFILES'//TRIM( coupling_char ), STATUS='OLD', & FORM='FORMATTED') ! !-- Skip header READ( 90, * ) DO k = nzb+1, nzt+1 READ( 90, * ) zz, luy, luz, tu(k), lvy, lvz, tv(k), lwy, lwz, tw(k), & r11(k), r21(k), r22(k), r31(k), r32(k), r33(k), & d1, d2, d3, d5 ! !-- Convert length scales from meter to number of grid points. nuy(k) = INT( luy * ddy ) nuz(k) = INT( luz * ddzw(k) ) nvy(k) = INT( lvy * ddy ) nvz(k) = INT( lvz * ddzw(k) ) nwy(k) = INT( lwy * ddy ) nwz(k) = INT( lwz * ddzw(k) ) ! !-- Workaround, assume isotropic turbulence nwx(k) = nwy(k) nvx(k) = nvy(k) nux(k) = nuy(k) ! !-- Save Mean inflow profiles IF ( TRIM( initializing_actions ) /= 'read_restart_data' ) THEN mean_inflow_profiles(k,1) = d1 mean_inflow_profiles(k,2) = d2 ! mean_inflow_profiles(k,4) = d4 mean_inflow_profiles(k,5) = d5 ENDIF ENDDO ! !-- Set lenght scales at surface grid point nuy(nzb) = nuy(nzb+1) nuz(nzb) = nuz(nzb+1) nvy(nzb) = nvy(nzb+1) nvz(nzb) = nvz(nzb+1) nwy(nzb) = nwy(nzb+1) nwz(nzb) = nwz(nzb+1) CLOSE( 90 ) ! !-- Calculate coefficient matrix from Reynolds stress tensor !-- (Lund rotation) CALL calc_coeff_matrix ! !-- No information about turbulence and its length scales are available. !-- Instead, parametrize turbulence which is imposed at the boundaries. !-- Set flag which indicates that turbulence is parametrized, which is done !-- later when energy-balance models are already initialized. This is !-- because the STG needs information about surface properties such as !-- roughness to build 'realistic' turbulence profiles. ELSE ! !-- Set flag indicating that turbulence is parametrized parametrize_inflow_turbulence = .TRUE. ! !-- Determine boundary-layer depth, which is used to initialize lenght !-- scales CALL calc_scaling_variables ! !-- Initialize lenght and time scales, which in turn are used !-- to initialize the filter functions. CALL calc_length_and_time_scale ENDIF ! !-- Assign initial profiles. Note, this is only required if turbulent !-- inflow from the left is desired, not in case of any of the !-- nesting (offline or self nesting) approaches. IF ( .NOT. nesting_offline .AND. .NOT. child_domain ) THEN u_init = mean_inflow_profiles(:,1) v_init = mean_inflow_profiles(:,2) !pt_init = mean_inflow_profiles(:,4) e_init = MAXVAL( mean_inflow_profiles(:,5) ) ENDIF ! !-- Define the size of the filter functions and allocate them. merg = 0 ! arrays must be large enough to cover the largest length scale DO k = nzb, nzt+1 j = MAX( ABS(nux(k)), ABS(nuy(k)), ABS(nuz(k)), & ABS(nvx(k)), ABS(nvy(k)), ABS(nvz(k)), & ABS(nwx(k)), ABS(nwy(k)), ABS(nwz(k)) ) IF ( j > merg ) merg = j ENDDO merg = 2 * merg mergp = merg + nbgp ALLOCATE ( bux(-merg:merg,nzb:nzt+1), & buy(-merg:merg,nzb:nzt+1), & buz(-merg:merg,nzb:nzt+1), & bvx(-merg:merg,nzb:nzt+1), & bvy(-merg:merg,nzb:nzt+1), & bvz(-merg:merg,nzb:nzt+1), & bwx(-merg:merg,nzb:nzt+1), & bwy(-merg:merg,nzb:nzt+1), & bwz(-merg:merg,nzb:nzt+1) ) ! !-- Allocate velocity seeds for turbulence at xz-layer ALLOCATE ( fu_xz( nzb:nzt+1,nxlg:nxrg), fuo_xz(nzb:nzt+1,nxlg:nxrg), & fv_xz( nzb:nzt+1,nxlg:nxrg), fvo_xz(nzb:nzt+1,nxlg:nxrg), & fw_xz( nzb:nzt+1,nxlg:nxrg), fwo_xz(nzb:nzt+1,nxlg:nxrg) ) ! !-- Allocate velocity seeds for turbulence at yz-layer ALLOCATE ( fu_yz( nzb:nzt+1,nysg:nyng), fuo_yz(nzb:nzt+1,nysg:nyng), & fv_yz( nzb:nzt+1,nysg:nyng), fvo_yz(nzb:nzt+1,nysg:nyng), & fw_yz( nzb:nzt+1,nysg:nyng), fwo_yz(nzb:nzt+1,nysg:nyng) ) fu_xz = 0.0_wp fuo_xz = 0.0_wp fv_xz = 0.0_wp fvo_xz = 0.0_wp fw_xz = 0.0_wp fwo_xz = 0.0_wp fu_yz = 0.0_wp fuo_yz = 0.0_wp fv_yz = 0.0_wp fvo_yz = 0.0_wp fw_yz = 0.0_wp fwo_yz = 0.0_wp ! !-- Create filter functions CALL stg_filter_func( nux, bux ) !filter ux CALL stg_filter_func( nuy, buy ) !filter uy CALL stg_filter_func( nuz, buz ) !filter uz CALL stg_filter_func( nvx, bvx ) !filter vx CALL stg_filter_func( nvy, bvy ) !filter vy CALL stg_filter_func( nvz, bvz ) !filter vz CALL stg_filter_func( nwx, bwx ) !filter wx CALL stg_filter_func( nwy, bwy ) !filter wy CALL stg_filter_func( nwz, bwz ) !filter wz #if defined( __parallel ) CALL MPI_BARRIER( comm2d, ierr ) #endif ! !-- In case of restart, calculate velocity seeds fu, fv, fw from former ! time step. ! Bug: fu, fv, fw are different in those heights where a11, a22, a33 ! are 0 compared to the prerun. This is mostly for k=nzt+1. IF ( TRIM( initializing_actions ) == 'read_restart_data' ) THEN IF ( myidx == id_stg_left .OR. myidx == id_stg_right ) THEN IF ( myidx == id_stg_left ) i = -1 IF ( myidx == id_stg_right ) i = nxr+1 DO j = nysg, nyng DO k = nzb, nzt+1 IF ( a11(k) .NE. 0._wp ) THEN fu_yz(k,j) = ( u(k,j,i) / mc_factor - u_init(k) ) / a11(k) ELSE fu_yz(k,j) = 0._wp ENDIF IF ( a22(k) .NE. 0._wp ) THEN fv_yz(k,j) = ( v(k,j,i) / mc_factor - a21(k) * fu_yz(k,j) - & v_init(k) ) / a22(k) ELSE fv_yz(k,j) = 0._wp ENDIF IF ( a33(k) .NE. 0._wp ) THEN fw_yz(k,j) = ( w(k,j,i) / mc_factor - a31(k) * fu_yz(k,j) - & a32(k) * fv_yz(k,j) ) / a33(k) ELSE fw_yz = 0._wp ENDIF ENDDO ENDDO ENDIF ENDIF END SUBROUTINE stg_init !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate filter function bxx from length scale nxx following Eg.9 and 10 !> (Xie and Castro, 2008) !------------------------------------------------------------------------------! SUBROUTINE stg_filter_func( nxx, bxx ) IMPLICIT NONE INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: n_k !< length scale nXX in height k INTEGER(iwp) :: n_k2 !< n_k * 2 INTEGER(iwp) :: nf !< index for length scales REAL(wp) :: bdenom !< denominator for filter functions bXX REAL(wp) :: qsi = 1.0_wp !< minimization factor INTEGER(iwp), DIMENSION(:) :: nxx(nzb:nzt+1) !< length scale (in gp) REAL(wp), DIMENSION(:,:) :: bxx(-merg:merg,nzb:nzt+1) !< filter function bxx = 0.0_wp DO k = nzb, nzt+1 bdenom = 0.0_wp n_k = nxx(k) IF ( n_k /= 0 ) THEN n_k2 = n_k * 2 ! !-- ( Eq.10 )^2 DO nf = -n_k2, n_k2 bdenom = bdenom + EXP( -qsi * pi * ABS(nf) / n_k )**2 ENDDO ! !-- ( Eq.9 ) bdenom = SQRT( bdenom ) DO nf = -n_k2, n_k2 bxx(nf,k) = EXP( -qsi * pi * ABS(nf) / n_k ) / bdenom ENDDO ENDIF ENDDO END SUBROUTINE stg_filter_func !------------------------------------------------------------------------------! ! Description: ! ------------ !> Parin for &stg_par for synthetic turbulence generator !------------------------------------------------------------------------------! SUBROUTINE stg_parin IMPLICIT NONE CHARACTER (LEN=80) :: line !< dummy string that contains the current line of the parameter file NAMELIST /stg_par/ dt_stg_adjust, dt_stg_call, use_syn_turb_gen line = ' ' ! !-- Try to find stg package REWIND ( 11 ) line = ' ' DO WHILE ( INDEX( line, '&stg_par' ) == 0 ) READ ( 11, '(A)', END=20 ) line ENDDO BACKSPACE ( 11 ) ! !-- Read namelist READ ( 11, stg_par, ERR = 10, END = 20 ) ! !-- Set flag that indicates that the synthetic turbulence generator is switched !-- on syn_turb_gen = .TRUE. GOTO 20 10 BACKSPACE( 11 ) READ( 11 , '(A)') line CALL parin_fail_message( 'stg_par', line ) 20 CONTINUE END SUBROUTINE stg_parin !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine reads the respective restart data. !------------------------------------------------------------------------------! SUBROUTINE stg_rrd_global( found ) USE control_parameters, & ONLY: length, restart_string IMPLICIT NONE LOGICAL, INTENT(OUT) :: found !< flag indicating if variable was found found = .TRUE. SELECT CASE ( restart_string(1:length) ) CASE ( 'mc_factor' ) READ ( 13 ) mc_factor CASE ( 'use_syn_turb_gen' ) READ ( 13 ) use_syn_turb_gen CASE DEFAULT found = .FALSE. END SELECT END SUBROUTINE stg_rrd_global !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine writes the respective restart data. !------------------------------------------------------------------------------! SUBROUTINE stg_wrd_global IMPLICIT NONE CALL wrd_write_string( 'mc_factor' ) WRITE ( 14 ) mc_factor CALL wrd_write_string( 'use_syn_turb_gen' ) WRITE ( 14 ) use_syn_turb_gen END SUBROUTINE stg_wrd_global !------------------------------------------------------------------------------! ! Description: ! ------------ !> Create turbulent inflow fields for u, v, w with prescribed length scales and !> Reynolds stress tensor after a method of Xie and Castro (2008), modified !> following suggestions of Kim et al. (2013), and using a Lund rotation !> (Lund, 1998). !------------------------------------------------------------------------------! SUBROUTINE stg_main USE arrays_3d, & ONLY: dzw USE control_parameters, & ONLY: child_domain, dt_3d, & nesting_offline, rans_mode, time_since_reference_point, & volume_flow_initial USE grid_variables, & ONLY: dx, dy USE indices, & ONLY: wall_flags_0 USE pmc_interface, & ONLY : rans_mode_parent IMPLICIT NONE INTEGER(iwp) :: i !< grid index in x-direction INTEGER(iwp) :: j !< loop index in y-direction INTEGER(iwp) :: k !< loop index in z-direction REAL(wp) :: dt_stg !< wheighted subtimestep REAL(wp) :: mc_factor_l !< local mass flux correction factor REAL(wp) :: volume_flow !< mass flux through lateral boundary REAL(wp) :: volume_flow_l !< local mass flux through lateral boundary ! !-- Calculate time step which is needed for filter functions dt_stg = MAX( dt_3d, dt_stg_call ) ! !-- Initial value of fu, fv, fw IF ( time_since_reference_point == 0.0_wp .AND. .NOT. velocity_seed_initialized ) THEN CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fu_yz, id_stg_left ) CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fv_yz, id_stg_left ) CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fw_yz, id_stg_left ) IF ( nesting_offline .OR. ( child_domain .AND. rans_mode_parent & .AND. .NOT. rans_mode ) ) THEN ! !-- Generate turbulence at right boundary CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fu_yz, id_stg_right ) CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fv_yz, id_stg_right ) CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fw_yz, id_stg_right ) ! !-- Generate turbulence at north boundary CALL stg_generate_seed_xz( nux, nuz, bux, buz, fu_xz, id_stg_north ) CALL stg_generate_seed_xz( nvx, nvz, bvx, bvz, fv_xz, id_stg_north ) CALL stg_generate_seed_xz( nwx, nwz, bwx, bwz, fw_xz, id_stg_north ) ! !-- Generate turbulence at south boundary CALL stg_generate_seed_xz( nux, nuz, bux, buz, fu_xz, id_stg_south ) CALL stg_generate_seed_xz( nvx, nvz, bvx, bvz, fv_xz, id_stg_south ) CALL stg_generate_seed_xz( nwx, nwz, bwx, bwz, fw_xz, id_stg_south ) ENDIF velocity_seed_initialized = .TRUE. ENDIF ! !-- New set of fu, fv, fw CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fuo_yz, id_stg_left ) CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fvo_yz, id_stg_left ) CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fwo_yz, id_stg_left ) IF ( nesting_offline .OR. ( child_domain .AND. rans_mode_parent & .AND. .NOT. rans_mode ) ) THEN ! !-- Generate turbulence at right boundary CALL stg_generate_seed_yz( nuy, nuz, buy, buz, fuo_yz, id_stg_right ) CALL stg_generate_seed_yz( nvy, nvz, bvy, bvz, fvo_yz, id_stg_right ) CALL stg_generate_seed_yz( nwy, nwz, bwy, bwz, fwo_yz, id_stg_right ) ! !-- Generate turbulence at north boundary CALL stg_generate_seed_xz( nux, nuz, bux, buz, fuo_xz, id_stg_north ) CALL stg_generate_seed_xz( nvx, nvz, bvx, bvz, fvo_xz, id_stg_north ) CALL stg_generate_seed_xz( nwx, nwz, bwx, bwz, fwo_xz, id_stg_north ) ! !-- Generate turbulence at south boundary CALL stg_generate_seed_xz( nux, nuz, bux, buz, fuo_xz, id_stg_south ) CALL stg_generate_seed_xz( nvx, nvz, bvx, bvz, fvo_xz, id_stg_south ) CALL stg_generate_seed_xz( nwx, nwz, bwx, bwz, fwo_xz, id_stg_south ) ENDIF ! !-- Turbulence generation at left and or right boundary IF ( myidx == id_stg_left .OR. myidx == id_stg_right ) THEN DO j = nysg, nyng DO k = nzb, nzt + 1 ! !-- Update fu, fv, fw following Eq. 14 of Xie and Castro (2008) IF ( tu(k) == 0.0_wp ) THEN fu_yz(k,j) = fuo_yz(k,j) ELSE fu_yz(k,j) = fu_yz(k,j) * EXP( -pi * dt_stg * 0.5_wp / tu(k) ) + & fuo_yz(k,j) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tu(k) ) ) ENDIF IF ( tv(k) == 0.0_wp ) THEN fv_yz(k,j) = fvo_yz(k,j) ELSE fv_yz(k,j) = fv_yz(k,j) * EXP( -pi * dt_stg * 0.5_wp / tv(k) ) + & fvo_yz(k,j) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tv(k) ) ) ENDIF IF ( tw(k) == 0.0_wp ) THEN fw_yz(k,j) = fwo_yz(k,j) ELSE fw_yz(k,j) = fw_yz(k,j) * EXP( -pi * dt_stg * 0.5_wp / tw(k) ) + & fwo_yz(k,j) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tw(k) ) ) ENDIF ! !-- Lund rotation following Eq. 17 in Xie and Castro (2008). !-- Additional factors are added to improve the variance of v and w IF( k == 0 ) THEN dist_yz(k,j,1) = 0.0_wp dist_yz(k,j,2) = 0.0_wp dist_yz(k,j,3) = 0.0_wp ELSE dist_yz(k,j,1) = MIN( a11(k) * fu_yz(k,j), 3.0_wp ) !experimental test of 1.2 dist_yz(k,j,2) = MIN( ( SQRT( a22(k) / MAXVAL(a22) ) & * 1.2_wp ) & * ( a21(k) * fu_yz(k,j) & + a22(k) * fv_yz(k,j) ), 3.0_wp ) dist_yz(k,j,3) = MIN( ( SQRT(a33(k) / MAXVAL(a33) ) & * 1.3_wp ) & * ( a31(k) * fu_yz(k,j) & + a32(k) * fv_yz(k,j) & + a33(k) * fw_yz(k,j) ), 3.0_wp ) ENDIF ENDDO ENDDO ! !-- Mass flux correction following Kim et al. (2013) !-- This correction factor insures that the mass flux is preserved at the !-- inflow boundary IF ( .NOT. nesting_offline .AND. .NOT. child_domain ) THEN mc_factor_l = 0.0_wp mc_factor = 0.0_wp DO j = nys, nyn DO k = nzb+1, nzt mc_factor_l = mc_factor_l + dzw(k) * & ( mean_inflow_profiles(k,1) + dist_yz(k,j,1) ) ENDDO ENDDO #if defined( __parallel ) CALL MPI_ALLREDUCE( mc_factor_l, mc_factor, 1, MPI_REAL, MPI_SUM, & comm1dy, ierr ) #else mc_factor = mc_factor_l #endif mc_factor = volume_flow_initial(1) / mc_factor ! !-- Add disturbance at the inflow DO j = nysg, nyng DO k = nzb, nzt+1 u(k,j,-nbgp+1:0) = ( mean_inflow_profiles(k,1) + & dist_yz(k,j,1) ) * mc_factor v(k,j,-nbgp:-1) = ( mean_inflow_profiles(k,2) + & dist_yz(k,j,2) ) * mc_factor w(k,j,-nbgp:-1) = dist_yz(k,j,3) * mc_factor ENDDO ENDDO ELSE ! !-- First, calculate volume flow at yz boundary IF ( myidx == id_stg_left ) i = nxl IF ( myidx == id_stg_right ) i = nxr+1 volume_flow_l = 0.0_wp volume_flow = 0.0_wp mc_factor_l = 0.0_wp mc_factor = 0.0_wp DO j = nys, nyn DO k = nzb+1, nzt volume_flow_l = volume_flow_l + u(k,j,i) * dzw(k) * dy & * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 1 ) ) mc_factor_l = mc_factor_l + ( u(k,j,i) + dist_yz(k,j,1) ) & * dzw(k) * dy & * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 1 ) ) ENDDO ENDDO #if defined( __parallel ) CALL MPI_ALLREDUCE( volume_flow_l, volume_flow, & 1, MPI_REAL, MPI_SUM, comm1dy, ierr ) CALL MPI_ALLREDUCE( mc_factor_l, mc_factor, & 1, MPI_REAL, MPI_SUM, comm1dy, ierr ) #else volume_flow = volume_flow_l mc_factor = mc_factor_l #endif mc_factor = volume_flow / mc_factor ! !-- Add disturbances IF ( myidx == id_stg_left ) THEN DO j = nys, nyn DO k = nzb+1, nzt u(k,j,0) = ( u(k,j,0) + dist_yz(k,j,1) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,0), 1 ) ) u(k,j,-1) = u(k,j,0) v(k,j,-1) = ( v(k,j,-1) + dist_yz(k,j,2) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,-1), 2 ) ) w(k,j,-1) = ( w(k,j,-1) + dist_yz(k,j,3) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,-1), 3 ) ) ENDDO ENDDO ENDIF IF ( myidx == id_stg_right ) THEN DO j = nys, nyn DO k = nzb+1, nzt u(k,j,nxr+1) = ( u(k,j,nxr+1) + dist_yz(k,j,1) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,nxr+1), 1 ) ) v(k,j,nxr+1) = ( v(k,j,nxr+1) + dist_yz(k,j,2) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,nxr+1), 2 ) ) w(k,j,nxr+1) = ( w(k,j,nxr+1) + dist_yz(k,j,3) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,nxr+1), 3 ) ) ENDDO ENDDO ENDIF ENDIF ENDIF ! !-- Turbulence generation at north and south boundary IF ( myidy == id_stg_north .OR. myidy == id_stg_south ) THEN DO i = nxlg, nxrg DO k = nzb, nzt + 1 ! !-- Update fu, fv, fw following Eq. 14 of Xie and Castro (2008) IF ( tu(k) == 0.0_wp ) THEN fu_xz(k,i) = fuo_xz(k,i) ELSE fu_xz(k,i) = fu_xz(k,i) * EXP( -pi * dt_stg * 0.5_wp / tu(k) ) + & fuo_xz(k,i) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tu(k) ) ) ENDIF IF ( tv(k) == 0.0_wp ) THEN fv_xz(k,i) = fvo_xz(k,i) ELSE fv_xz(k,i) = fv_xz(k,i) * EXP( -pi * dt_stg * 0.5_wp / tv(k) ) + & fvo_xz(k,i) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tv(k) ) ) ENDIF IF ( tw(k) == 0.0_wp ) THEN fw_xz(k,i) = fwo_xz(k,i) ELSE fw_xz(k,i) = fw_xz(k,i) * EXP( -pi * dt_stg * 0.5_wp / tw(k) ) + & fwo_xz(k,i) * SQRT( 1.0_wp - EXP( -pi * dt_stg / tw(k) ) ) ENDIF ! !-- Lund rotation following Eq. 17 in Xie and Castro (2008). !-- Additional factors are added to improve the variance of v and w IF( k == 0 ) THEN dist_xz(k,i,1) = 0.0_wp dist_xz(k,i,2) = 0.0_wp dist_xz(k,i,3) = 0.0_wp ELSE dist_xz(k,i,1) = MIN( a11(k) * fu_xz(k,i), 3.0_wp ) !experimental test of 1.2 dist_xz(k,i,2) = MIN( ( SQRT( a22(k) / MAXVAL(a22) ) & * 1.2_wp ) & * ( a21(k) * fu_xz(k,i) & + a22(k) * fv_xz(k,i) ), 3.0_wp ) dist_xz(k,i,3) = MIN( ( SQRT(a33(k) / MAXVAL(a33) ) & * 1.3_wp ) & * ( a31(k) * fu_xz(k,i) & + a32(k) * fv_xz(k,i) & + a33(k) * fw_xz(k,i) ), 3.0_wp ) ENDIF ENDDO ENDDO ! !-- Mass flux correction following Kim et al. (2013) !-- This correction factor insures that the mass flux is preserved at the !-- inflow boundary. !-- First, calculate volume flow at xz boundary IF ( myidy == id_stg_south ) j = nys IF ( myidy == id_stg_north ) j = nyn+1 volume_flow_l = 0.0_wp volume_flow = 0.0_wp mc_factor_l = 0.0_wp mc_factor = 0.0_wp DO i = nxl, nxr DO k = nzb+1, nzt volume_flow_l = volume_flow_l + v(k,j,i) * dzw(k) * dx & * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 2 ) ) mc_factor_l = mc_factor_l + ( v(k,j,i) + dist_xz(k,i,2) ) & * dzw(k) * dx & * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 2 ) ) ENDDO ENDDO #if defined( __parallel ) CALL MPI_ALLREDUCE( volume_flow_l, volume_flow, & 1, MPI_REAL, MPI_SUM, comm1dx, ierr ) CALL MPI_ALLREDUCE( mc_factor_l, mc_factor, & 1, MPI_REAL, MPI_SUM, comm1dx, ierr ) #else volume_flow = volume_flow_l mc_factor = mc_factor_l #endif mc_factor = volume_flow / mc_factor ! !-- Add disturbances IF ( myidy == id_stg_south ) THEN DO i = nxl, nxr DO k = nzb+1, nzt u(k,-1,i) = ( u(k,-1,i) + dist_xz(k,i,1) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,-1,i), 1 ) ) v(k,0,i) = ( v(k,0,i) + dist_xz(k,i,2) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,0,i), 2 ) ) v(k,-1,i) = v(k,0,i) w(k,-1,i) = ( w(k,-1,i) + dist_xz(k,i,3) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,-1,i), 3 ) ) ENDDO ENDDO ENDIF IF ( myidy == id_stg_north ) THEN DO i = nxl, nxr DO k = nzb+1, nzt u(k,nyn+1,i) = ( u(k,nyn+1,i) + dist_xz(k,i,1) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,nyn+1,i), 1 ) ) v(k,nyn+1,i) = ( v(k,nyn+1,i) + dist_xz(k,i,2) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,nyn+1,i), 2 ) ) w(k,nyn+1,i) = ( w(k,nyn+1,i) + dist_xz(k,i,3) ) & * mc_factor * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,nyn+1,i), 3 ) ) ENDDO ENDDO ENDIF ENDIF ! !-- Finally, set time counter for calling STG to zero time_stg_call = 0.0_wp END SUBROUTINE stg_main !------------------------------------------------------------------------------! ! Description: ! ------------ !> Generate a set of random number rand_it wich is equal on each PE !> and calculate the velocity seed f_n. !> f_n is splitted in vertical direction by the number of PEs in x-direction and !> and each PE calculates a vertical subsection of f_n. At the the end, all !> parts are collected to form the full array. !------------------------------------------------------------------------------! SUBROUTINE stg_generate_seed_yz( n_y, n_z, b_y, b_z, f_n, id ) USE indices, & ONLY: ny IMPLICIT NONE INTEGER(iwp) :: id !< core ids at respective boundaries INTEGER(iwp) :: j !< loop index in y-direction INTEGER(iwp) :: jj !< loop index in y-direction INTEGER(iwp) :: k !< loop index in z-direction INTEGER(iwp) :: kk !< loop index in z-direction INTEGER(iwp) :: send_count !< send count for MPI_GATHERV INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_y !< length scale in y-direction INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_z !< length scale in z-direction INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_y2 !< n_y*2 INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_z2 !< n_z*2 REAL(wp) :: nyz_inv !< inverse of number of grid points in yz-slice REAL(wp) :: rand_av !< average of random number REAL(wp) :: rand_sigma_inv !< inverse of stdev of random number REAL(wp), DIMENSION(-merg:merg,nzb:nzt+1) :: b_y !< filter function in y-direction REAL(wp), DIMENSION(-merg:merg,nzb:nzt+1) :: b_z !< filter function in z-direction REAL(wp), DIMENSION(nzb_x_stg:nzt_x_stg+1,nysg:nyng) :: f_n_l !< local velocity seed REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng) :: f_n !< velocity seed REAL(wp), DIMENSION(:,:), ALLOCATABLE :: rand_it !< random number ! !-- Generate random numbers using a seed generated in stg_init. !-- The set of random numbers are modified to have an average of 0 and !-- unit variance. ALLOCATE( rand_it(nzb-mergp:nzt+1+mergp,-mergp:ny+mergp) ) rand_av = 0.0_wp rand_sigma_inv = 0.0_wp nyz_inv = 1.0_wp / REAL( ( nzt+1 - nzb+1 ) * ( ny+1 ), KIND=wp ) DO j = 0, ny DO k = nzb, nzt+1 CALL RANDOM_NUMBER( rand_it(k,j) ) rand_av = rand_av + rand_it(k,j) ENDDO ENDDO rand_av = rand_av * nyz_inv DO j = 0, ny DO k = nzb, nzt+1 rand_it(k,j) = rand_it(k,j) - rand_av rand_sigma_inv = rand_sigma_inv + rand_it(k,j) ** 2 ENDDO ENDDO rand_sigma_inv = 1.0_wp / SQRT(rand_sigma_inv * nyz_inv) DO j = 0, ny DO k = nzb, nzt+1 rand_it(k,j) = rand_it(k,j) * rand_sigma_inv ENDDO ENDDO ! !-- Periodic fill of random number in space DO j = 0, ny DO k = 1, mergp rand_it(nzb -k,j) = rand_it(nzt+2-k,j) ! bottom margin rand_it(nzt+1+k,j) = rand_it(nzb+k-1,j) ! top margin ENDDO ENDDO DO j = 1, mergp DO k = nzb-mergp, nzt+1+mergp rand_it(k, -j) = rand_it(k,ny-j+1) ! south margin rand_it(k,ny+j) = rand_it(k, j-1) ! north margin ENDDO ENDDO ! !-- Generate velocity seed following Eq.6 of Xie and Castro (2008) n_y2 = n_y * 2 n_z2 = n_z * 2 f_n_l = 0.0_wp DO j = nysg, nyng DO k = nzb_x_stg, nzt_x_stg+1 DO jj = -n_y2(k), n_y2(k) DO kk = -n_z2(k), n_z2(k) f_n_l(k,j) = f_n_l(k,j) & + b_y(jj,k) * b_z(kk,k) * rand_it(k+kk,j+jj) ENDDO ENDDO ENDDO ENDDO DEALLOCATE( rand_it ) ! !-- Gather velocity seeds of full subdomain send_count = nzt_x_stg - nzb_x_stg + 1 IF ( nzt_x_stg == nzt ) send_count = send_count + 1 #if defined( __parallel ) CALL MPI_GATHERV( f_n_l(nzb_x_stg,nysg), send_count, stg_type_yz_small, & f_n(nzb+1,nysg), recv_count_yz, displs_yz, stg_type_yz, & id, comm1dx, ierr ) #else f_n(nzb+1:nzt+1,nysg:nyng) = f_n_l(nzb_x_stg:nzt_x_stg+1,nysg:nyng) #endif END SUBROUTINE stg_generate_seed_yz !------------------------------------------------------------------------------! ! Description: ! ------------ !> Generate a set of random number rand_it wich is equal on each PE !> and calculate the velocity seed f_n. !> f_n is splitted in vertical direction by the number of PEs in y-direction and !> and each PE calculates a vertical subsection of f_n. At the the end, all !> parts are collected to form the full array. !------------------------------------------------------------------------------! SUBROUTINE stg_generate_seed_xz( n_x, n_z, b_x, b_z, f_n, id ) USE indices, & ONLY: nx IMPLICIT NONE INTEGER(iwp) :: id !< core ids at respective boundaries INTEGER(iwp) :: i !< loop index in x-direction INTEGER(iwp) :: ii !< loop index in x-direction INTEGER(iwp) :: k !< loop index in z-direction INTEGER(iwp) :: kk !< loop index in z-direction INTEGER(iwp) :: send_count !< send count for MPI_GATHERV INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_x !< length scale in x-direction INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_z !< length scale in z-direction INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_x2 !< n_x*2 INTEGER(iwp), DIMENSION(nzb:nzt+1) :: n_z2 !< n_z*2 REAL(wp) :: nxz_inv !< inverse of number of grid points in xz-slice REAL(wp) :: rand_av !< average of random number REAL(wp) :: rand_sigma_inv !< inverse of stdev of random number REAL(wp), DIMENSION(-merg:merg,nzb:nzt+1) :: b_x !< filter function in x-direction REAL(wp), DIMENSION(-merg:merg,nzb:nzt+1) :: b_z !< filter function in z-direction REAL(wp), DIMENSION(nzb_y_stg:nzt_y_stg+1,nxlg:nxrg) :: f_n_l !< local velocity seed REAL(wp), DIMENSION(nzb:nzt+1,nxlg:nxrg) :: f_n !< velocity seed REAL(wp), DIMENSION(:,:), ALLOCATABLE :: rand_it !< random number ! !-- Generate random numbers using a seed generated in stg_init. !-- The set of random numbers are modified to have an average of 0 and !-- unit variance. ALLOCATE( rand_it(nzb-mergp:nzt+1+mergp,-mergp:nx+mergp) ) rand_av = 0.0_wp rand_sigma_inv = 0.0_wp nxz_inv = 1.0_wp / REAL( ( nzt+1 - nzb+1 ) * ( nx+1 ), KIND=wp ) DO i = 0, nx DO k = nzb, nzt+1 CALL RANDOM_NUMBER( rand_it(k,i) ) rand_av = rand_av + rand_it(k,i) ENDDO ENDDO rand_av = rand_av * nxz_inv DO i = 0, nx DO k = nzb, nzt+1 rand_it(k,i) = rand_it(k,i) - rand_av rand_sigma_inv = rand_sigma_inv + rand_it(k,i) ** 2 ENDDO ENDDO rand_sigma_inv = 1.0_wp / SQRT(rand_sigma_inv * nxz_inv) DO i = 0, nx DO k = nzb, nzt+1 rand_it(k,i) = rand_it(k,i) * rand_sigma_inv ENDDO ENDDO ! !-- Periodic fill of random number in space DO i = 0, nx DO k = 1, mergp rand_it(nzb-k,i) = rand_it(nzt+2-k,i) ! bottom margin rand_it(nzt+1+k,i) = rand_it(nzb+k-1,i) ! top margin ENDDO ENDDO DO i = 1, mergp DO k = nzb-mergp, nzt+1+mergp rand_it(k,-i) = rand_it(k,nx-i+1) ! left margin rand_it(k,nx+i) = rand_it(k,i-1) ! right margin ENDDO ENDDO ! !-- Generate velocity seed following Eq.6 of Xie and Castro (2008) n_x2 = n_x * 2 n_z2 = n_z * 2 f_n_l = 0.0_wp DO i = nxlg, nxrg DO k = nzb_y_stg, nzt_y_stg+1 DO ii = -n_x2(k), n_x2(k) DO kk = -n_z2(k), n_z2(k) f_n_l(k,i) = f_n_l(k,i) & + b_x(ii,k) * b_z(kk,k) * rand_it(k+kk,i+ii) ENDDO ENDDO ENDDO ENDDO DEALLOCATE( rand_it ) ! !-- Gather velocity seeds of full subdomain send_count = nzt_y_stg - nzb_y_stg + 1 IF ( nzt_y_stg == nzt ) send_count = send_count + 1 #if defined( __parallel ) CALL MPI_GATHERV( f_n_l(nzb_y_stg,nxlg), send_count, stg_type_xz_small, & f_n(nzb+1,nxlg), recv_count_xz, displs_xz, stg_type_xz, & id, comm1dy, ierr ) #else f_n(nzb+1:nzt+1,nxlg:nxrg) = f_n_l(nzb_y_stg:nzt_y_stg+1,nxlg:nxrg) #endif END SUBROUTINE stg_generate_seed_xz !------------------------------------------------------------------------------! ! Description: ! ------------ !> Parametrization of the Reynolds stress tensor, following the parametrization !> described in Rotach et al. (1996), which is applied in state-of-the-art !> dispserion modelling. Please note, the parametrization does not distinguish !> between along-wind and cross-wind turbulence. !------------------------------------------------------------------------------! SUBROUTINE parametrize_reynolds_stress USE arrays_3d, & ONLY: zu IMPLICIT NONE INTEGER(iwp) :: k !< loop index in z-direction REAL(wp) :: zzi !< ratio of z/zi ! !-- DO k = nzb+1, nzt+1 IF ( zu(k) <= zi_ribulk ) THEN ! !-- Determine normalized height coordinate zzi = zu(k) / zi_ribulk ! !-- u'u' and v'v'. Assume isotropy. Note, add a small negative number !-- to the denominator, else the merge-function can crash if scale_l is !-- zero. r11(k) = scale_us**2 * ( & MERGE( 0.35_wp * ( & - zi_ribulk / ( kappa * scale_l - 10E-4_wp ) & )**( 2.0_wp / 3.0_wp ), & 0.0_wp, & scale_l < 0.0_wp ) & + 5.0_wp - 4.0_wp * zzi & ) r22(k) = r11(k) ! !-- w'w' r33(k) = scale_wm**2 * ( & 1.5_wp * zzi**( 2.0_wp / 3.0_wp ) * EXP( -2.0_wp * zzi ) & + ( 1.7_wp - zzi ) * ( scale_us / scale_wm )**2 & ) ! !-- u'w' and v'w'. Assume isotropy. r31(k) = - scale_us**2 * ( & 1.0_wp - EXP( 3.0_wp * ( zzi - 1.0_wp ) ) & ) r32(k) = r31(k) ! !-- For u'v' no parametrization exist so far - ?. For simplicity assume !-- a similar profile as for u'w'. r21(k) = r31(k) ! !-- Above the boundary layer, assmume laminar flow conditions. ELSE r11(k) = 10E-8_wp r22(k) = 10E-8_wp r33(k) = 10E-8_wp r21(k) = 10E-8_wp r31(k) = 10E-8_wp r32(k) = 10E-8_wp ENDIF ENDDO ! !-- Set bottom boundary condition r11(nzb) = r11(nzb+1) r22(nzb) = r22(nzb+1) r33(nzb) = r33(nzb+1) r21(nzb) = r11(nzb+1) r31(nzb) = r31(nzb+1) r32(nzb) = r32(nzb+1) END SUBROUTINE parametrize_reynolds_stress !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the coefficient matrix from the Lund rotation. !------------------------------------------------------------------------------! SUBROUTINE calc_coeff_matrix IMPLICIT NONE INTEGER(iwp) :: k !< loop index in z-direction ! !-- Calculate coefficient matrix. Split loops to allow for loop vectorization. DO k = nzb+1, nzt+1 IF ( r11(k) > 0.0_wp ) THEN a11(k) = SQRT( r11(k) ) a21(k) = r21(k) / a11(k) a31(k) = r31(k) / a11(k) ELSE a11(k) = 10E-8_wp a21(k) = 10E-8_wp a31(k) = 10E-8_wp ENDIF ENDDO DO k = nzb+1, nzt+1 a22(k) = r22(k) - a21(k)**2 IF ( a22(k) > 0.0_wp ) THEN a22(k) = SQRT( a22(k) ) a32(k) = r32(k) - a21(k) * a31(k) / a22(k) ELSE a22(k) = 10E-8_wp a32(k) = 10E-8_wp ENDIF ENDDO DO k = nzb+1, nzt+1 a33(k) = r33(k) - a31(k)**2 - a32(k)**2 IF ( a33(k) > 0.0_wp ) THEN a33(k) = SQRT( a33(k) ) ELSE a33(k) = 10E-8_wp ENDIF ENDDO ! !-- Set bottom boundary condition a11(nzb) = a11(nzb+1) a22(nzb) = a22(nzb+1) a21(nzb) = a21(nzb+1) a33(nzb) = a33(nzb+1) a31(nzb) = a31(nzb+1) a32(nzb) = a32(nzb+1) END SUBROUTINE calc_coeff_matrix !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine controls the re-adjustment of the turbulence statistics used !> for generating turbulence at the lateral boundaries. !------------------------------------------------------------------------------! SUBROUTINE stg_adjust IMPLICIT NONE ! !-- Compute mean boundary layer height according to Richardson-Bulk !-- criterion using the inflow profiles. Further velocity scale as well as !-- mean friction velocity are calculated. CALL calc_scaling_variables ! !-- Set length and time scales depending on boundary-layer height CALL calc_length_and_time_scale ! !-- Parametrize Reynolds-stress tensor, diagonal elements as well !-- as r21 (v'u'), r31 (w'u'), r32 (w'v'). Parametrization follows !-- Rotach et al. (1996) and is based on boundary-layer depth, !-- friction velocity and velocity scale. CALL parametrize_reynolds_stress ! !-- Calculate coefficient matrix from Reynolds stress tensor !-- (Lund rotation) CALL calc_coeff_matrix ! !-- Determine filter functions on basis of updated length scales CALL stg_filter_func( nux, bux ) !filter ux CALL stg_filter_func( nuy, buy ) !filter uy CALL stg_filter_func( nuz, buz ) !filter uz CALL stg_filter_func( nvx, bvx ) !filter vx CALL stg_filter_func( nvy, bvy ) !filter vy CALL stg_filter_func( nvz, bvz ) !filter vz CALL stg_filter_func( nwx, bwx ) !filter wx CALL stg_filter_func( nwy, bwy ) !filter wy CALL stg_filter_func( nwz, bwz ) !filter wz ! !-- Reset time counter for controlling next adjustment to zero time_stg_adjust = 0.0_wp END SUBROUTINE stg_adjust !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates turbuluent length and time scales if these are not available !> from measurements. !------------------------------------------------------------------------------! SUBROUTINE calc_length_and_time_scale USE arrays_3d, & ONLY: dzw, ddzw, u_init, v_init USE grid_variables, & ONLY: ddx, ddy, dx, dy IMPLICIT NONE INTEGER(iwp) :: k !< loop index in z-direction REAL(wp) :: length_scale !< typical length scale ! !-- In initial call the boundary-layer depth can be zero. This case, set !-- minimum value for boundary-layer depth, to setup length scales correctly. zi_ribulk = MAX( zi_ribulk, zw(nzb+2) ) ! !-- Set-up default turbulent length scales. From the numerical point of !-- view the imposed perturbations should not be immediately dissipated !-- by the numerics. The numerical dissipation, however, acts on scales !-- up to 8 x the grid spacing. For this reason, set the turbulence !-- length scale to 8 time the grid spacing. Further, above the boundary !-- layer height, set turbulence lenght scales to zero (equivalent to not !-- imposing any perturbations) in order to save computational costs. !-- Typical time scales are derived by assuming Taylors's hypothesis, !-- using the length scales and the mean profiles of u- and v-component. DO k = nzb+1, nzt+1 length_scale = 8.0_wp * MIN( dx, dy, dzw(k) ) IF ( zu(k) <= zi_ribulk ) THEN ! !-- Assume isotropic turbulence length scales nux(k) = MAX( INT( length_scale * ddx ), 1 ) nuy(k) = MAX( INT( length_scale * ddy ), 1 ) nuz(k) = MAX( INT( length_scale * ddzw(k) ), 1 ) nvx(k) = MAX( INT( length_scale * ddx ), 1 ) nvy(k) = MAX( INT( length_scale * ddy ), 1 ) nvz(k) = MAX( INT( length_scale * ddzw(k) ), 1 ) nwx(k) = MAX( INT( length_scale * ddx ), 1 ) nwy(k) = MAX( INT( length_scale * ddy ), 1 ) nwz(k) = MAX( INT( length_scale * ddzw(k) ), 1 ) ! !-- Limit time scales, else they become very larger for low wind speed, !-- imposing long-living inflow perturbations which in turn propagate !-- further into the model domain. Use u_init and v_init to calculate !-- the time scales, which will be equal to the inflow profiles, both, !-- in offline nesting mode or in dirichlet/radiation mode. tu(k) = MIN( dt_stg_adjust, length_scale / & ( ABS( u_init(k) ) + 0.1_wp ) ) tv(k) = MIN( dt_stg_adjust, length_scale / & ( ABS( v_init(k) ) + 0.1_wp ) ) ! !-- Time scale of w-component is a mixture from u- and v-component. tw(k) = SQRT( tu(k)**2 + tv(k)**2 ) ! !-- Above the boundary layer length scales are zero, i.e. imposed turbulence !-- is not correlated in space and time, just white noise. This saves !-- computations power. ELSE nux(k) = 0.0_wp nuy(k) = 0.0_wp nuz(k) = 0.0_wp nvx(k) = 0.0_wp nvy(k) = 0.0_wp nvz(k) = 0.0_wp nwx(k) = 0.0_wp nwy(k) = 0.0_wp nwz(k) = 0.0_wp tu(k) = 0.0_wp tv(k) = 0.0_wp tw(k) = 0.0_wp ENDIF ENDDO ! !-- Set bottom boundary condition for the length and time scales nux(nzb) = nux(nzb+1) nuy(nzb) = nuy(nzb+1) nuz(nzb) = nuz(nzb+1) nvx(nzb) = nvx(nzb+1) nvy(nzb) = nvy(nzb+1) nvz(nzb) = nvz(nzb+1) nwx(nzb) = nwx(nzb+1) nwy(nzb) = nwy(nzb+1) nwz(nzb) = nwz(nzb+1) tu(nzb) = tu(nzb+1) tv(nzb) = tv(nzb+1) tw(nzb) = tw(nzb+1) END SUBROUTINE calc_length_and_time_scale !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate scaling variables which are used for turbulence parametrization !> according to Rotach et al. (1996). Scaling variables are: friction velocity, !> boundary-layer depth, momentum velocity scale, and Obukhov length. !------------------------------------------------------------------------------! SUBROUTINE calc_scaling_variables USE control_parameters, & ONLY: bc_dirichlet_l, bc_dirichlet_n, bc_dirichlet_r, bc_dirichlet_s, & pt_surface USE indices, & ONLY: nx, ny USE surface_mod, & ONLY: get_topography_top_index_ji, surf_def_h, surf_lsm_h, surf_usm_h IMPLICIT NONE INTEGER(iwp) :: i !< loop index in x-direction INTEGER(iwp) :: j !< loop index in y-direction INTEGER(iwp) :: k !< loop index in z-direction INTEGER(iwp) :: k_ref !< index in z-direction for COSMO reference height INTEGER(iwp) :: k_topo !< topography top index INTEGER(iwp) :: m !< surface element index REAL(wp) :: friction_vel_l !< mean friction veloctiy on subdomain REAL(wp) :: shf_mean !< mean surface sensible heat flux REAL(wp) :: shf_mean_l !< mean surface sensible heat flux on subdomain REAL(wp) :: u_int !< u-component REAL(wp) :: v_int !< v-component REAL(wp) :: w_convective !< convective velocity scale REAL(wp) :: z0_mean !< mean roughness length REAL(wp) :: z0_mean_l !< mean roughness length on subdomain ! !-- Mean friction velocity and velocity scale. Therefore, !-- pre-calculate mean roughness length and surface sensible heat flux !-- in the model domain, which are further used to estimate friction !-- velocity and velocity scale. Note, for z0 linear averaging is applied, !-- even though this is known to unestimate the effective roughness. !-- This need to be revised later. z0_mean_l = 0.0_wp shf_mean_l = 0.0_wp DO m = 1, surf_def_h(0)%ns z0_mean_l = z0_mean_l + surf_def_h(0)%z0(m) shf_mean_l = shf_mean_l + surf_def_h(0)%shf(m) ENDDO DO m = 1, surf_lsm_h%ns z0_mean_l = z0_mean_l + surf_lsm_h%z0(m) shf_mean_l = shf_mean_l + surf_lsm_h%shf(m) ENDDO DO m = 1, surf_usm_h%ns z0_mean_l = z0_mean_l + surf_usm_h%z0(m) shf_mean_l = shf_mean_l + surf_usm_h%shf(m) ENDDO #if defined( __parallel ) CALL MPI_ALLREDUCE( z0_mean_l, z0_mean, 1, MPI_REAL, MPI_SUM, & comm2d, ierr ) CALL MPI_ALLREDUCE( shf_mean_l, shf_mean, 1, MPI_REAL, MPI_SUM, & comm2d, ierr ) #else z0_mean = z0_mean_l shf_mean = shf_mean_l #endif z0_mean = z0_mean / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp ) shf_mean = shf_mean / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp ) ! !-- Note, Inifor does not use logarithmic interpolation of the !-- velocity components near the ground, so that near-surface !-- wind speed and thus the friction velocity is overestimated. !-- However, friction velocity is used for turbulence !-- parametrization, so that more physically meaningful values are important. !-- Hence, derive friction velocity from wind speed at a reference height, !-- which is 10 m, according to the height of the 1st vertical grid level !-- in the COSMO level. However, in case of topography that is higher than !-- the reference level, the k index is determined from the 1st vertical !-- PALM grid level instead. !-- For a first guess use 20 m, which is in the range of the first !-- COSMO vertical level. k_ref = MINLOC( ABS( zu - cosmo_ref ), DIM = 1 ) - 1 friction_vel_l = 0.0_wp IF ( bc_dirichlet_l .OR. bc_dirichlet_r ) THEN i = MERGE( -1, nxr + 1, bc_dirichlet_l ) DO j = nys, nyn ! !-- Determine the k index and topography top index k_topo = MAX( get_topography_top_index_ji( j, i, 'u' ), & get_topography_top_index_ji( j, i, 'v' ) ) k = MAX( k_ref, k_topo + 1 ) ! !-- Note, in u- and v- component the imposed perturbations !-- from the STG are already included. Check whether this !-- makes any difference compared to using the pure-mean !-- inflow profiles. u_int = MERGE( u(k,j,i+1), u(k,j,i), bc_dirichlet_l ) v_int = v(k,j,i) ! !-- Calculate friction velocity and sum-up. Therefore, assume !-- neutral condtions. friction_vel_l = friction_vel_l + kappa * & SQRT( u_int * u_int + v_int * v_int ) / & LOG( ( zu(k) - zu(k_topo) ) / z0_mean ) ENDDO ENDIF IF ( bc_dirichlet_s .OR. bc_dirichlet_n ) THEN j = MERGE( -1, nyn + 1, bc_dirichlet_s ) DO i = nxl, nxr k_topo = MAX( get_topography_top_index_ji( j, i, 'u' ), & get_topography_top_index_ji( j, i, 'v' ) ) k = MAX( k_ref, k_topo + 1 ) u_int = u(k,j,i) v_int = MERGE( v(k,j+1,i), v(k,j,i), bc_dirichlet_s ) friction_vel_l = friction_vel_l + kappa * & SQRT( u_int * u_int + v_int * v_int ) / & LOG( ( zu(k) - zu(k_topo) ) / z0_mean ) ENDDO ENDIF #if defined( __parallel ) CALL MPI_ALLREDUCE( friction_vel_l, scale_us, 1, MPI_REAL, MPI_SUM, & comm2d, ierr ) #else scale_us = friction_vel_l #endif scale_us = scale_us / REAL( 2 * nx + 2 * ny, KIND = wp ) ! !-- Compute mean Obukhov length IF ( shf_mean > 0.0_wp ) THEN scale_l = - pt_surface / ( kappa * g ) * scale_us**3 / shf_mean ELSE scale_l = 0.0_wp ENDIF ! !-- Compute mean convective velocity scale. Note, in case the mean heat flux !-- is negative, set convective velocity scale to zero. IF ( shf_mean > 0.0_wp ) THEN w_convective = ( g * shf_mean * zi_ribulk / pt_surface )**( 1.0_wp / 3.0_wp ) ELSE w_convective = 0.0_wp ENDIF ! !-- Finally, in order to consider also neutral or stable stratification, !-- compute momentum velocity scale from u* and convective velocity scale, !-- according to Rotach et al. (1996). scale_wm = ( scale_us**3 + 0.6_wp * w_convective**3 )**( 1.0_wp / 3.0_wp ) END SUBROUTINE calc_scaling_variables END MODULE synthetic_turbulence_generator_mod