source: palm/trunk/SOURCE/init_3d_model.f90 @ 4845

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

maximum phase velocities are alwasy used for radiation boundary conditions, parameter use_cmax removed

  • Property svn:keywords set to Id
  • Property svn:mergeinfo set to (toggle deleted branches)
    /palm/branches/chemistry/SOURCE/init_3d_model.f902047-3190,​3218-3297
    /palm/branches/rans/SOURCE/init_3d_model.f902078-3128
    /palm/branches/resler/SOURCE/init_3d_model.f902023-4492
    /palm/branches/salsa/SOURCE/init_3d_model.f902503-3581
    /palm/branches/forwind/SOURCE/init_3d_model.f901564-1913
    /palm/branches/mosaik_M2/init_3d_model.f902360-3471
    /palm/branches/palm4u/SOURCE/init_3d_model.f902540-2692
File size: 84.3 KB
RevLine 
[1682]1!> @file init_3d_model.f90
[4648]2!--------------------------------------------------------------------------------------------------!
[2696]3! This file is part of the PALM model system.
[1036]4!
[4648]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.
[1036]8!
[4648]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.
[1036]12!
[4648]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/>.
[1036]15!
[4828]16! Copyright 1997-2021 Leibniz Universitaet Hannover
[4648]17!--------------------------------------------------------------------------------------------------!
[1036]18!
[254]19! Current revisions:
[732]20! ------------------
[4648]21!
22!
[2233]23! Former revisions:
24! -----------------
25! $Id: init_3d_model.f90 4845 2021-01-18 11:15:37Z raasch $
[4845]26! removed allocation and initialization of arrays required for radiation boundary conditions
27!
28! 4828 2021-01-05 11:21:41Z Giersch
[4791]29! Add possibility to initialize surface flux of passive scalar via static driver
30!
31! 4783 2020-11-13 13:58:45Z raasch
[4783]32! bugfix for reading restart data with MPI-I/O (does not work with blockwise I/O)
33!
34! 4680 2020-09-16 10:20:34Z gronemeier
[4680]35! Add option to fix date or time of the simulation
36!
37! 4671 2020-09-09 20:27:58Z pavelkrc
[4671]38! Implementation of downward facing USM and LSM surfaces
39!
40! 4648 2020-08-25 07:52:08Z raasch
[4648]41! file re-formatted to follow the PALM coding standard
42!
43! 4548 2020-05-28 19:36:45Z suehring
[4548]44! Bugfix, move call for lsf_forcing_surf after lsf_init is called
[4648]45!
[4548]46! 4514 2020-04-30 16:29:59Z suehring
[4648]47! Add possibility to initialize surface sensible and latent heat fluxes via a static driver.
48!
[4514]49! 4493 2020-04-10 09:49:43Z pavelkrc
[4648]50! Overwrite u_init, v_init, pt_init, q_init and s_init with hom for all cyclic_fill-cases, not only
51! for turbulent_inflow = .TRUE.
52!
[4365]53! 4360 2020-01-07 11:25:50Z suehring
[4648]54! Introduction of wall_flags_total_0, which currently sets bits based on static topography
55! information used in wall_flags_static_0
56!
[4346]57! 4329 2019-12-10 15:46:36Z motisi
[4329]58! Renamed wall_flags_0 to wall_flags_static_0
[4648]59!
[4329]60! 4286 2019-10-30 16:01:14Z resler
[4227]61! implement new palm_date_time_mod
[4648]62!
[4227]63! 4223 2019-09-10 09:20:47Z gronemeier
[4648]64! Deallocate temporary string array since it may be re-used to read different input data in other
65! modules
66!
[4187]67! 4186 2019-08-23 16:06:14Z suehring
[4648]68! Design change, use variables defined in netcdf_data_input_mod to read netcd variables rather than
69! define local ones.
70!
[4186]71! 4185 2019-08-23 13:49:38Z oliver.maas
[4185]72! For initializing_actions = ' cyclic_fill':
[4648]73! Overwrite u_init, v_init, pt_init, q_init and s_init with the (temporally) and horizontally
74! averaged vertical profiles from the end of the prerun, because these profiles shall be used as the
75! basic state for the rayleigh damping and the pt_damping.
76!
[4185]77! 4182 2019-08-22 15:20:23Z scharf
[4182]78! Corrected "Former revisions" section
[4648]79!
[4182]80! 4168 2019-08-16 13:50:17Z suehring
[4168]81! Replace function get_topography_top_index by topo_top_ind
[4648]82!
[4168]83! 4151 2019-08-09 08:24:30Z suehring
[4151]84! Add netcdf directive around input calls (fix for last commit)
[4648]85!
[4151]86! 4150 2019-08-08 20:00:47Z suehring
[4648]87! Input of additional surface variables independent on land- or urban-surface model
88!
[4150]89! 4131 2019-08-02 11:06:18Z monakurppa
[4131]90! Allocate sums and sums_l to allow profile output for salsa variables.
[4648]91!
[4131]92! 4130 2019-08-01 13:04:13Z suehring
[4648]93! Effectively reduce 3D initialization to 1D initial profiles. This is because 3D initialization
94! produces structures in the w-component that are correlated with the processor grid for some
95! unknown reason
96!
[4130]97! 4090 2019-07-11 15:06:47Z Giersch
[4090]98! Unused variables removed
[4648]99!
[4090]100! 4088 2019-07-11 13:57:56Z Giersch
[4088]101! Pressure and density profile calculations revised using basic functions
[4648]102!
[4088]103! 4048 2019-06-21 21:00:21Z knoop
[4028]104! Further modularization of particle code components
[4648]105!
[4028]106! 4017 2019-06-06 12:16:46Z schwenkel
[4648]107! Convert most location messages to debug messages to reduce output in job logfile to a minimum
108!
[3939]109! unused variable removed
[4648]110!
[3939]111! 3937 2019-04-29 15:09:07Z suehring
[4648]112! Move initialization of synthetic turbulence generator behind initialization of offline nesting.
113! Remove call for stg_adjust, as this is now already done in stg_init.
114!
[3937]115! 3900 2019-04-16 15:17:43Z suehring
[3900]116! Fix problem with LOD = 2 initialization
[4648]117!
[3900]118! 3885 2019-04-11 11:29:34Z kanani
[4648]119! Changes related to global restructuring of location messages and introduction of additional debug
120! messages
121!
[3885]122! 3849 2019-04-01 16:35:16Z knoop
[3747]123! Move initialization of rmask before initializing user_init_arrays
[4648]124!
[3747]125! 3711 2019-01-31 13:44:26Z knoop
[3711]126! Introduced module_interface_init_checks for post-init checks in modules
[4648]127!
[3711]128! 3700 2019-01-26 17:03:42Z knoop
[3685]129! Some interface calls moved to module_interface + cleanup
[4648]130!
[3685]131! 3648 2019-01-02 16:35:46Z suehring
[3648]132! Rename subroutines for surface-data output
[3569]133!
[4182]134! Revision 1.1  1998/03/09 16:22:22  raasch
135! Initial revision
136!
137!
[1]138! Description:
139! ------------
[1682]140!> Allocation of arrays and initialization of the 3D model via
141!> a) pre-run the 1D model
142!> or
143!> b) pre-set constant linear profiles
144!> or
145!> c) read values of a previous run
[4648]146!--------------------------------------------------------------------------------------------------!
[1682]147 SUBROUTINE init_3d_model
[1]148
[3298]149
[667]150    USE advec_ws
[1320]151
[1]152    USE arrays_3d
[1849]153
[4648]154    USE basic_constants_and_equations_mod,                                                         &
155        ONLY:  barometric_formula, c_p, exner_function, exner_function_invers, g,                  &
156               ideal_gas_law_rho, ideal_gas_law_rho_pt, l_v, pi
[3274]157
[4648]158    USE bulk_cloud_model_mod,                                                                      &
[3685]159        ONLY:  bulk_cloud_model
[3274]160
[4648]161    USE chem_modules,                                                                              &
[3685]162        ONLY:  max_pr_cs ! ToDo: this dependency needs to be removed cause it is ugly #new_dom
[3298]163
[1]164    USE control_parameters
[3298]165
[4648]166    USE grid_variables,                                                                            &
[2037]167        ONLY:  dx, dy, ddx2_mg, ddy2_mg
[2817]168
[1]169    USE indices
[3469]170
[1320]171    USE kinds
[4648]172
173    USE lsf_nudging_mod,                                                                           &
[3685]174        ONLY:  ls_forcing_surf
[1849]175
[4648]176    USE model_1d_mod,                                                                              &
[3241]177        ONLY:  init_1d_model, l1d, u1d, v1d
[2338]178
[4648]179    USE module_interface,                                                                          &
180        ONLY:  module_interface_init_arrays,                                                       &
181               module_interface_init,                                                              &
[3711]182               module_interface_init_checks
[3685]183
[4648]184    USE multi_agent_system_mod,                                                                    &
[3159]185        ONLY:  agents_active, mas_init
186
[4648]187    USE netcdf_interface,                                                                          &
[3700]188        ONLY:  dots_max
[2696]189
[4648]190    USE netcdf_data_input_mod,                                                                     &
191        ONLY:  char_fill,                                                                          &
192               check_existence,                                                                    &
193               close_input_file,                                                                   &
194               get_attribute,                                                                      &
195               get_variable,                                                                       &
196               init_3d,                                                                            &
197               input_pids_static,                                                                  &
198               inquire_num_variables,                                                              &
199               inquire_variable_names,                                                             &
200               input_file_static,                                                                  &
201               netcdf_data_input_init_3d,                                                          &
202               num_var_pids,                                                                       &
203               open_read_file,                                                                     &
204               pids_id,                                                                            &
205               real_2d,                                                                            &
[4186]206               vars_pids
[4648]207
208    USE nesting_offl_mod,                                                                          &
[3347]209        ONLY:  nesting_offl_init
[3294]210
[4648]211    USE palm_date_time_mod,                                                                        &
[4680]212        ONLY:  init_date_time
[4227]213
[1]214    USE pegrid
[3298]215
[3524]216#if defined( __parallel )
[4648]217    USE pmc_interface,                                                                             &
[2934]218        ONLY:  nested_run
[3524]219#endif
[2934]220
[4648]221    USE random_function_mod
[3685]222
[4648]223    USE random_generator_parallel,                                                                 &
[2172]224        ONLY:  init_parallel_random_generator
[3685]225
[4648]226    USE read_restart_data_mod,                                                                     &
227        ONLY:  rrd_local, rrd_read_parts_of_global
[3685]228
[4648]229    USE statistics,                                                                                &
230        ONLY:  hom, hom_sum, mean_surface_level_height, pr_palm, rmask, statistic_regions, sums,   &
231               sums_divnew_l, sums_divold_l, sums_l, sums_l_l, sums_wsts_bc_l, ts_value,           &
[1833]232               weight_pres, weight_substep
[2259]233
[4648]234    USE synthetic_turbulence_generator_mod,                                                        &
[3939]235        ONLY:  stg_init, use_syn_turb_gen
[3685]236
[4648]237    USE surface_layer_fluxes_mod,                                                                  &
[1691]238        ONLY:  init_surface_layer_fluxes
[2232]239
[4648]240    USE surface_mod,                                                                               &
241        ONLY :  init_single_surface_properties,                                                    &
242                init_surface_arrays,                                                               &
243                init_surfaces,                                                                     &
244                surf_def_h,                                                                        &
245                surf_def_v,                                                                        &
246                surf_lsm_h,                                                                        &
[4168]247                surf_usm_h
[3685]248
[3849]249#if defined( _OPENACC )
[4648]250    USE surface_mod,                                                                               &
[3849]251        ONLY :  bc_h
252#endif
253
[4648]254    USE surface_data_output_mod,                                                                   &
[3648]255        ONLY:  surface_data_output_init
[3685]256
[2007]257    USE transpose_indices
[1]258
[4648]259
[1]260    IMPLICIT NONE
[4648]261
[4150]262    INTEGER(iwp) ::  i                    !< grid index in x direction
263    INTEGER(iwp) ::  ind_array(1)         !< dummy used to determine start index for external pressure forcing
264    INTEGER(iwp) ::  j                    !< grid index in y direction
265    INTEGER(iwp) ::  k                    !< grid index in z direction
266    INTEGER(iwp) ::  k_surf               !< surface level index
[4648]267    INTEGER(iwp) ::  l                    !< running index over surface orientation
268    INTEGER(iwp) ::  m                    !< index of surface element in surface data type
[4150]269    INTEGER(iwp) ::  nz_u_shift           !< topography-top index on u-grid, used to vertically shift initial profiles
270    INTEGER(iwp) ::  nz_v_shift           !< topography-top index on v-grid, used to vertically shift initial profiles
271    INTEGER(iwp) ::  nz_w_shift           !< topography-top index on w-grid, used to vertically shift initial profiles
272    INTEGER(iwp) ::  nz_s_shift           !< topography-top index on scalar-grid, used to vertically shift initial profiles
273    INTEGER(iwp) ::  nz_u_shift_l         !< topography-top index on u-grid, used to vertically shift initial profiles
274    INTEGER(iwp) ::  nz_v_shift_l         !< topography-top index on v-grid, used to vertically shift initial profiles
275    INTEGER(iwp) ::  nz_w_shift_l         !< topography-top index on w-grid, used to vertically shift initial profiles
276    INTEGER(iwp) ::  nz_s_shift_l         !< topography-top index on scalar-grid, used to vertically shift initial profiles
277    INTEGER(iwp) ::  nzt_l                !< index of top PE boundary for multigrid level
278    INTEGER(iwp) ::  sr                   !< index of statistic region
[1]279
[4648]280    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  ngp_2dh_l  !< toal number of horizontal grid points in statistical region on
281                                                             !< subdomain
[1]282
[4648]283    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l    !< number of horizontal non-wall bounded grid points on
284                                                                     !< subdomain
285    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner_l  !< number of horizontal non-topography grid points on
286                                                                     !< subdomain
[1]287
[2037]288    REAL(wp) ::  dx_l !< grid spacing along x on different multigrid level
289    REAL(wp) ::  dy_l !< grid spacing along y on different multigrid level
290
[4648]291    REAL(wp), DIMENSION(:), ALLOCATABLE ::  init_l                       !< dummy array used for averaging 3D data to obtain
292                                                                         !< inital profiles
293    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height_l  !< mean surface level height on subdomain
294    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l               !< total number of non-topography grid points on subdomain
295    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_tmp             !< total number of non-topography grid points
296    REAL(wp), DIMENSION(:), ALLOCATABLE ::  p_hydrostatic                !< hydrostatic pressure
297
[3547]298    REAL(wp), DIMENSION(1:3) ::  volume_flow_area_l     !< area of lateral and top model domain surface on local subdomain
299    REAL(wp), DIMENSION(1:3) ::  volume_flow_initial_l  !< initial volume flow into model domain
[1]300
[4648]301    TYPE(real_2d) ::  tmp_2d !< temporary variable to input additional surface-data from static file
[1]302
[3987]303    CALL location_message( 'model initialization', 'start' )
[4227]304!
305!-- Set reference date-time
[4680]306    CALL init_date_time( date_time_str=origin_date_time,                                           &
307                         use_fixed_date=use_fixed_date,                                            &
308                         use_fixed_time=use_fixed_time )
[3987]309
310    IF ( debug_output )  CALL debug_message( 'allocating arrays', 'start' )
[1]311!
312!-- Allocate arrays
[4648]313    ALLOCATE( mean_surface_level_height(0:statistic_regions),                                      &
314              mean_surface_level_height_l(0:statistic_regions),                                    &
315              ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions),                        &
316              ngp_3d(0:statistic_regions),                                                         &
317              ngp_3d_inner(0:statistic_regions),                                                   &
318              ngp_3d_inner_l(0:statistic_regions),                                                 &
319              ngp_3d_inner_tmp(0:statistic_regions),                                               &
320              sums_divnew_l(0:statistic_regions),                                                  &
[1]321              sums_divold_l(0:statistic_regions) )
[1195]322    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
[4648]323    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                                        &
324              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),                                      &
325              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),                                      &
326              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),                                    &
327              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),                                      &
328              sums(nzb:nzt+1,pr_palm+max_pr_user+max_pr_cs+max_pr_salsa),                          &
[4131]329              sums_l(nzb:nzt+1,pr_palm+max_pr_user+max_pr_cs+max_pr_salsa,0:threads_per_task-1),   &
[4648]330              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1),                        &
[3700]331              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions) )
332    ALLOCATE( ts_value(dots_max,0:statistic_regions) )
[978]333    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
[1]334
[4648]335    ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr),                                                        &
336              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                    &
[1010]337              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
338
[4648]339    ALLOCATE( pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                 &
340              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                 &
341              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
342              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
343              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
344              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
345              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
346              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
347              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
348              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                                  &
[667]349              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1788]350    IF (  .NOT.  neutral )  THEN
[1032]351       ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
352    ENDIF
[673]353!
[3747]354!-- Pre-set masks for regional statistics. Default is the total model domain.
[4648]355!-- Ghost points are excluded because counting values at the ghost boundaries would bias the
356!-- statistics.
[3747]357    rmask = 1.0_wp
358    rmask(:,nxlg:nxl-1,:) = 0.0_wp;  rmask(:,nxr+1:nxrg,:) = 0.0_wp
359    rmask(nysg:nys-1,:,:) = 0.0_wp;  rmask(nyn+1:nyng,:,:) = 0.0_wp
360!
[4648]361!-- Following array is required for perturbation pressure within the iterative pressure solvers. For
362!-- the multistep schemes (Runge-Kutta), array p holds the weighted average of the substeps and
363!-- cannot be used in the Poisson solver.
[707]364    IF ( psolver == 'sor' )  THEN
365       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1575]366    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
[707]367!
368!--    For performance reasons, multigrid is using one ghost layer only
369       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
[673]370    ENDIF
[1]371
[1111]372!
373!-- Array for storing constant coeffficients of the tridiagonal solver
374    IF ( psolver == 'poisfft' )  THEN
[1212]375       ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) )
[1111]376       ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) )
377    ENDIF
378
[1960]379    IF ( humidity )  THEN
[1]380!
[1960]381!--    3D-humidity
[4648]382       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                               &
383                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                               &
384                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                               &
385                 vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
386    ENDIF
387
[1960]388    IF ( passive_scalar )  THEN
[1]389
[1960]390!
391!--    3D scalar arrays
[4648]392       ALLOCATE( s_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                               &
393                 s_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                               &
[1960]394                 s_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[3636]395
[1960]396    ENDIF
397
[1]398!
[4648]399!-- Allocate and set 1d-profiles for Stokes drift velocity. It may be set to non-zero values later
400!-- in ocean_init.
401    ALLOCATE( u_stokes_zu(nzb:nzt+1), u_stokes_zw(nzb:nzt+1),                                      &
[3302]402              v_stokes_zu(nzb:nzt+1), v_stokes_zw(nzb:nzt+1) )
403    u_stokes_zu(:) = 0.0_wp
404    u_stokes_zw(:) = 0.0_wp
405    v_stokes_zu(:) = 0.0_wp
406    v_stokes_zw(:) = 0.0_wp
407
408!
[2037]409!-- Allocation of anelastic and Boussinesq approximation specific arrays
410    ALLOCATE( p_hydrostatic(nzb:nzt+1) )
411    ALLOCATE( rho_air(nzb:nzt+1) )
412    ALLOCATE( rho_air_zw(nzb:nzt+1) )
413    ALLOCATE( drho_air(nzb:nzt+1) )
414    ALLOCATE( drho_air_zw(nzb:nzt+1) )
415!
[4648]416!-- Density profile calculation for anelastic and Boussinesq approximation.
417!-- In case of a Boussinesq approximation, a constant density is calculated mainly for output
418!-- purposes. This density does not need to be considered in the model's system of equations.
[4088]419    IF ( TRIM( approximation ) == 'anelastic' )  THEN
[2037]420       DO  k = nzb, nzt+1
[4648]421          p_hydrostatic(k) = barometric_formula(zu(k), pt_surface *                                &
422                                                exner_function(surface_pressure * 100.0_wp),       &
423                                                surface_pressure * 100.0_wp)
424
[4088]425          rho_air(k) = ideal_gas_law_rho_pt(p_hydrostatic(k), pt_init(k))
[2037]426       ENDDO
[4648]427
[2037]428       DO  k = nzb, nzt
429          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
430       ENDDO
[4648]431
432       rho_air_zw(nzt+1)  = rho_air_zw(nzt) + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
433
[2037]434    ELSE
[2252]435       DO  k = nzb, nzt+1
[4648]436          p_hydrostatic(k) = barometric_formula(zu(nzb), pt_surface *                              &
437                                                exner_function(surface_pressure * 100.0_wp),       &
438                                                surface_pressure * 100.0_wp)
[4088]439
440          rho_air(k) = ideal_gas_law_rho_pt(p_hydrostatic(k), pt_init(nzb))
[2252]441       ENDDO
[4648]442
[2252]443       DO  k = nzb, nzt
444          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
445       ENDDO
[4648]446
447       rho_air_zw(nzt+1)  = rho_air_zw(nzt) + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
448
[2037]449    ENDIF
[2696]450!
[4648]451!-- Compute the inverse density array in order to avoid expencive divisions
[2037]452    drho_air    = 1.0_wp / rho_air
453    drho_air_zw = 1.0_wp / rho_air_zw
454
455!
456!-- Allocation of flux conversion arrays
457    ALLOCATE( heatflux_input_conversion(nzb:nzt+1) )
458    ALLOCATE( waterflux_input_conversion(nzb:nzt+1) )
459    ALLOCATE( momentumflux_input_conversion(nzb:nzt+1) )
460    ALLOCATE( heatflux_output_conversion(nzb:nzt+1) )
461    ALLOCATE( waterflux_output_conversion(nzb:nzt+1) )
462    ALLOCATE( momentumflux_output_conversion(nzb:nzt+1) )
463
464!
[4648]465!-- Calculate flux conversion factors according to approximation and in-/output mode
[2037]466    DO  k = nzb, nzt+1
467
468        IF ( TRIM( flux_input_mode ) == 'kinematic' )  THEN
469            heatflux_input_conversion(k)      = rho_air_zw(k)
470            waterflux_input_conversion(k)     = rho_air_zw(k)
471            momentumflux_input_conversion(k)  = rho_air_zw(k)
472        ELSEIF ( TRIM( flux_input_mode ) == 'dynamic' ) THEN
[3274]473            heatflux_input_conversion(k)      = 1.0_wp / c_p
[2037]474            waterflux_input_conversion(k)     = 1.0_wp / l_v
475            momentumflux_input_conversion(k)  = 1.0_wp
476        ENDIF
477
478        IF ( TRIM( flux_output_mode ) == 'kinematic' )  THEN
479            heatflux_output_conversion(k)     = drho_air_zw(k)
480            waterflux_output_conversion(k)    = drho_air_zw(k)
481            momentumflux_output_conversion(k) = drho_air_zw(k)
482        ELSEIF ( TRIM( flux_output_mode ) == 'dynamic' ) THEN
[3274]483            heatflux_output_conversion(k)     = c_p
[2037]484            waterflux_output_conversion(k)    = l_v
485            momentumflux_output_conversion(k) = 1.0_wp
486        ENDIF
487
488        IF ( .NOT. humidity ) THEN
489            waterflux_input_conversion(k)  = 1.0_wp
490            waterflux_output_conversion(k) = 1.0_wp
491        ENDIF
492
493    ENDDO
494
495!
[4648]496!-- In case of multigrid method, compute grid lengths and grid factors for the grid levels with
497!-- respective density on each grid.
[2037]498    IF ( psolver(1:9) == 'multigrid' )  THEN
499
500       ALLOCATE( ddx2_mg(maximum_grid_level) )
501       ALLOCATE( ddy2_mg(maximum_grid_level) )
502       ALLOCATE( dzu_mg(nzb+1:nzt+1,maximum_grid_level) )
503       ALLOCATE( dzw_mg(nzb+1:nzt+1,maximum_grid_level) )
504       ALLOCATE( f1_mg(nzb+1:nzt,maximum_grid_level) )
505       ALLOCATE( f2_mg(nzb+1:nzt,maximum_grid_level) )
506       ALLOCATE( f3_mg(nzb+1:nzt,maximum_grid_level) )
507       ALLOCATE( rho_air_mg(nzb:nzt+1,maximum_grid_level) )
508       ALLOCATE( rho_air_zw_mg(nzb:nzt+1,maximum_grid_level) )
509
510       dzu_mg(:,maximum_grid_level) = dzu
511       rho_air_mg(:,maximum_grid_level) = rho_air
[4648]512!
513!--    Next line to ensure an equally spaced grid.
[2037]514       dzu_mg(1,maximum_grid_level) = dzu(2)
[4648]515       rho_air_mg(nzb,maximum_grid_level) = rho_air(nzb) + (rho_air(nzb) - rho_air(nzb+1))
[2037]516
517       dzw_mg(:,maximum_grid_level) = dzw
518       rho_air_zw_mg(:,maximum_grid_level) = rho_air_zw
519       nzt_l = nzt
520       DO  l = maximum_grid_level-1, 1, -1
521           dzu_mg(nzb+1,l) = 2.0_wp * dzu_mg(nzb+1,l+1)
522           dzw_mg(nzb+1,l) = 2.0_wp * dzw_mg(nzb+1,l+1)
[4648]523           rho_air_mg(nzb,l)    = rho_air_mg(nzb,l+1)    + ( rho_air_mg(nzb,l+1)    -              &
524                                                             rho_air_mg(nzb+1,l+1)    )
525           rho_air_zw_mg(nzb,l) = rho_air_zw_mg(nzb,l+1) + ( rho_air_zw_mg(nzb,l+1) -              &
526                                                             rho_air_zw_mg(nzb+1,l+1) )
[2037]527           rho_air_mg(nzb+1,l)    = rho_air_mg(nzb+1,l+1)
528           rho_air_zw_mg(nzb+1,l) = rho_air_zw_mg(nzb+1,l+1)
529           nzt_l = nzt_l / 2
530           DO  k = 2, nzt_l+1
531              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
532              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
533              rho_air_mg(k,l)    = rho_air_mg(2*k-1,l+1)
534              rho_air_zw_mg(k,l) = rho_air_zw_mg(2*k-1,l+1)
535           ENDDO
536       ENDDO
537
538       nzt_l = nzt
539       dx_l  = dx
540       dy_l  = dy
541       DO  l = maximum_grid_level, 1, -1
542          ddx2_mg(l) = 1.0_wp / dx_l**2
543          ddy2_mg(l) = 1.0_wp / dy_l**2
544          DO  k = nzb+1, nzt_l
545             f2_mg(k,l) = rho_air_zw_mg(k,l) / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
546             f3_mg(k,l) = rho_air_zw_mg(k-1,l) / ( dzu_mg(k,l)   * dzw_mg(k,l) )
[4648]547             f1_mg(k,l) = 2.0_wp * ( ddx2_mg(l) + ddy2_mg(l) )                                     &
[2037]548                          * rho_air_mg(k,l) + f2_mg(k,l) + f3_mg(k,l)
549          ENDDO
550          nzt_l = nzt_l / 2
551          dx_l  = dx_l * 2.0_wp
552          dy_l  = dy_l * 2.0_wp
553       ENDDO
554
555    ENDIF
556
557!
[1299]558!-- 1D-array for large scale subsidence velocity
[1361]559    IF ( .NOT. ALLOCATED( w_subs ) )  THEN
560       ALLOCATE ( w_subs(nzb:nzt+1) )
561       w_subs = 0.0_wp
562    ENDIF
[1299]563
564!
[1]565!-- Initial assignment of the pointers
[1032]566    IF ( .NOT. neutral )  THEN
567       pt => pt_1;  pt_p => pt_2;  tpt_m => pt_3
568    ELSE
569       pt => pt_1;  pt_p => pt_1;  tpt_m => pt_3
570    ENDIF
[1001]571    u  => u_1;   u_p  => u_2;   tu_m  => u_3
572    v  => v_1;   v_p  => v_2;   tv_m  => v_3
573    w  => w_1;   w_p  => w_2;   tw_m  => w_3
[1]574
[1960]575    IF ( humidity )  THEN
[1001]576       q => q_1;  q_p => q_2;  tq_m => q_3
[3274]577       vpt  => vpt_1
[1001]578    ENDIF
[4648]579
[1960]580    IF ( passive_scalar )  THEN
581       s => s_1;  s_p => s_2;  ts_m => s_3
[4648]582    ENDIF
[1]583
584!
[2696]585!-- Initialize surface arrays
[2232]586    CALL init_surface_arrays
587!
[3294]588!-- Allocate arrays for other modules
[3685]589    CALL module_interface_init_arrays
[1551]590
[1914]591
[2320]592!
[4648]593!-- Allocate arrays containing the RK coefficient for calculation of perturbation pressure and
594!-- turbulent fluxes. At this point values are set for pressure calculation during initialization
595!-- (where no timestep is done). Further below the values needed within the timestep scheme will be
596!-- set.
597    ALLOCATE( weight_substep(1:intermediate_timestep_count_max),                                   &
[1878]598              weight_pres(1:intermediate_timestep_count_max) )
[1340]599    weight_substep = 1.0_wp
600    weight_pres    = 1.0_wp
[1918]601    intermediate_timestep_count = 0  ! needed when simulated_time = 0.0
[4648]602
[3987]603    IF ( debug_output )  CALL debug_message( 'allocating arrays', 'end' )
[1918]604
[673]605!
[3014]606!-- Initialize time series
607    ts_value = 0.0_wp
608
609!
[1918]610!-- Initialize local summation arrays for routine flow_statistics.
[4648]611!-- This is necessary because they may not yet have been initialized when they are called from
612!-- flow_statistics (or - depending on the chosen model run - are never initialized)
613    sums_divnew_l  = 0.0_wp
614    sums_divold_l  = 0.0_wp
615    sums_l_l       = 0.0_wp
616    sums_wsts_bc_l = 0.0_wp
617
[1918]618!
[1]619!-- Initialize model variables
[4648]620    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.                                &
[328]621         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
[1]622!
[2696]623!--    Initialization with provided input data derived from larger-scale model
624       IF ( INDEX( initializing_actions, 'inifor' ) /= 0 )  THEN
[3987]625          IF ( debug_output )  CALL debug_message( 'initializing with INIFOR', 'start' )
[2696]626!
[4648]627!--       Read initial 1D profiles or 3D data from NetCDF file, depending on the provided
628!--       level-of-detail.
629!--       At the moment, only u, v, w, pt and q are provided.
[2696]630          CALL netcdf_data_input_init_3d
631!
[4648]632!--       Please note, Inifor provides data from nzb+1 to nzt.
633!--       Bottom and top boundary conditions for Inifor profiles are already set (just after
634!--       reading), so that this is not necessary here.
635!--       Depending on the provided level-of-detail, initial Inifor data is either stored on data
636!--       type (lod=1), or directly on 3D arrays (lod=2).
637!--       In order to obtain also initial profiles in case of lod=2 (which is required for e.g.
638!--       damping), average over 3D data.
[3182]639          IF( init_3d%lod_u == 1 )  THEN
640             u_init = init_3d%u_init
[4648]641          ELSEIF( init_3d%lod_u == 2 )  THEN
642             ALLOCATE( init_l(nzb:nzt+1) )
[3182]643             DO  k = nzb, nzt+1
644                init_l(k) = SUM( u(k,nys:nyn,nxl:nxr) )
645             ENDDO
646             init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
[1384]647
[3182]648#if defined( __parallel )
[4648]649             CALL MPI_ALLREDUCE( init_l, u_init, nzt+1-nzb+1, MPI_REAL, MPI_SUM, comm2d, ierr )
[3182]650#else
651             u_init = init_l
652#endif
653             DEALLOCATE( init_l )
[3051]654
[2696]655          ENDIF
[4648]656
657          IF( init_3d%lod_v == 1 )  THEN
[3182]658             v_init = init_3d%v_init
[4648]659          ELSEIF( init_3d%lod_v == 2 )  THEN
660             ALLOCATE( init_l(nzb:nzt+1) )
[3182]661             DO  k = nzb, nzt+1
662                init_l(k) = SUM( v(k,nys:nyn,nxl:nxr) )
663             ENDDO
664             init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
[2696]665
[3182]666#if defined( __parallel )
[4648]667             CALL MPI_ALLREDUCE( init_l, v_init, nzt+1-nzb+1, MPI_REAL, MPI_SUM, comm2d, ierr )
[3182]668#else
669             v_init = init_l
670#endif
671             DEALLOCATE( init_l )
672          ENDIF
673          IF( .NOT. neutral )  THEN
674             IF( init_3d%lod_pt == 1 )  THEN
675                pt_init = init_3d%pt_init
[4648]676             ELSEIF( init_3d%lod_pt == 2 )  THEN
677                ALLOCATE( init_l(nzb:nzt+1) )
[3182]678                DO  k = nzb, nzt+1
679                   init_l(k) = SUM( pt(k,nys:nyn,nxl:nxr) )
680                ENDDO
681                init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
682
683#if defined( __parallel )
[4648]684                CALL MPI_ALLREDUCE( init_l, pt_init, nzt+1-nzb+1, MPI_REAL, MPI_SUM, comm2d, ierr )
[3182]685#else
686                pt_init = init_l
687#endif
688                DEALLOCATE( init_l )
689             ENDIF
690          ENDIF
691
692
693          IF( humidity )  THEN
694             IF( init_3d%lod_q == 1 )  THEN
695                q_init = init_3d%q_init
[4648]696             ELSEIF( init_3d%lod_q == 2 )  THEN
697                ALLOCATE( init_l(nzb:nzt+1) )
[3182]698                DO  k = nzb, nzt+1
699                   init_l(k) = SUM( q(k,nys:nyn,nxl:nxr) )
700                ENDDO
701                init_l = init_l / REAL( ( nx + 1 ) * ( ny + 1 ), KIND = wp )
702
703#if defined( __parallel )
[4648]704                CALL MPI_ALLREDUCE( init_l, q_init, nzt+1-nzb+1, MPI_REAL, MPI_SUM, comm2d, ierr )
[3182]705#else
706                q_init = init_l
707#endif
708                DEALLOCATE( init_l )
709             ENDIF
710          ENDIF
711
[2696]712!
[4648]713!--       Write initial profiles onto 3D arrays.
714!--       Work-around, 3D initialization of u,v,w creates artificial structures which correlate with
715!--       the processor grid. The reason for this is still unknown. To work-around this, 3D
716!--       initialization will be effectively reduce to a 1D initialization where no such artificial
717!--       structures appear.
[2696]718          DO  i = nxlg, nxrg
719             DO  j = nysg, nyng
[4648]720                IF( init_3d%lod_u == 1  .OR.  init_3d%lod_u == 2 )  u(:,j,i) = u_init(:)
721                IF( init_3d%lod_v == 1  .OR.  init_3d%lod_u == 2 )  v(:,j,i) = v_init(:)
722                IF( .NOT. neutral  .AND.  ( init_3d%lod_pt == 1  .OR.  init_3d%lod_pt == 2 ) )     &
[3051]723                   pt(:,j,i) = pt_init(:)
[4648]724                IF( humidity  .AND.  ( init_3d%lod_q == 1  .OR.  init_3d%lod_q == 2 ) )            &
[4130]725                   q(:,j,i) = q_init(:)
[2696]726             ENDDO
727          ENDDO
728!
[4648]729!--       Set geostrophic wind components.
[2938]730          IF ( init_3d%from_file_ug )  THEN
731             ug(:) = init_3d%ug_init(:)
732          ENDIF
733          IF ( init_3d%from_file_vg )  THEN
734             vg(:) = init_3d%vg_init(:)
735          ENDIF
[3404]736!
737!--       Set bottom and top boundary condition for geostrophic wind
[2938]738          ug(nzt+1) = ug(nzt)
739          vg(nzt+1) = vg(nzt)
[3404]740          ug(nzb)   = ug(nzb+1)
741          vg(nzb)   = vg(nzb+1)
[2696]742!
743!--       Set inital w to 0
744          w = 0.0_wp
745
746          IF ( passive_scalar )  THEN
747             DO  i = nxlg, nxrg
748                DO  j = nysg, nyng
749                   s(:,j,i) = s_init
750                ENDDO
751             ENDDO
752          ENDIF
753
754!
[4648]755!--       Set velocity components at non-atmospheric / oceanic grid points to zero.
[4346]756          u = MERGE( u, 0.0_wp, BTEST( wall_flags_total_0, 1 ) )
757          v = MERGE( v, 0.0_wp, BTEST( wall_flags_total_0, 2 ) )
758          w = MERGE( w, 0.0_wp, BTEST( wall_flags_total_0, 3 ) )
[2700]759!
[4648]760!--       Initialize surface variables, e.g. friction velocity, momentum fluxes, etc.
761          CALL  init_surfaces
[2696]762
[4648]763          IF ( debug_output )  CALL  debug_message( 'initializing with INIFOR', 'end' )
[2696]764!
765!--    Initialization via computed 1D-model profiles
766       ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
767
[4648]768          IF ( debug_output )  CALL  debug_message( 'initializing with 1D model profiles', 'start' )
[1]769!
770!--       Use solutions of the 1D model as initial profiles,
771!--       start 1D model
772          CALL init_1d_model
773!
774!--       Transfer initial profiles to the arrays of the 3D model
[667]775          DO  i = nxlg, nxrg
776             DO  j = nysg, nyng
[1]777                pt(:,j,i) = pt_init
778                u(:,j,i)  = u1d
779                v(:,j,i)  = v1d
780             ENDDO
781          ENDDO
782
[1960]783          IF ( humidity )  THEN
[667]784             DO  i = nxlg, nxrg
785                DO  j = nysg, nyng
[1]786                   q(:,j,i) = q_init
787                ENDDO
788             ENDDO
789          ENDIF
[2292]790
[1960]791          IF ( passive_scalar )  THEN
792             DO  i = nxlg, nxrg
793                DO  j = nysg, nyng
794                   s(:,j,i) = s_init
795                ENDDO
[4648]796             ENDDO
[1960]797          ENDIF
[1]798!
799!--          Store initial profiles for output purposes etc.
[2696]800          IF ( .NOT. constant_diffusion )  THEN
[1]801             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
802          ENDIF
803!
[2696]804!--       Set velocities back to zero
[4346]805          u = MERGE( u, 0.0_wp, BTEST( wall_flags_total_0, 1 ) )
[4648]806          v = MERGE( v, 0.0_wp, BTEST( wall_flags_total_0, 2 ) )
[1]807!
[4648]808!--       WARNING: The extra boundary conditions set after running the 1D model impose an error on
809!--       -------- the divergence one layer below the topography; need to correct later
810!--       ATTENTION: Provisional correction for Piacsek & Williams advection scheme: keep u and v
811!--       ---------- zero one layer below the topography.
[2696]812          IF ( ibc_uv_b == 1 )  THEN
[667]813!
[2696]814!--          Neumann condition
815             DO  i = nxl-1, nxr+1
816                DO  j = nys-1, nyn+1
817                   u(nzb,j,i) = u(nzb+1,j,i)
818                   v(nzb,j,i) = v(nzb+1,j,i)
[1]819                ENDDO
[2696]820             ENDDO
[1]821
822          ENDIF
[2618]823!
[4648]824!--       Initialize surface variables, e.g. friction velocity, momentum fluxes, etc.
[2618]825          CALL init_surfaces
[1]826
[4648]827          IF ( debug_output )  CALL  debug_message( 'initializing with 1D model profiles', 'end' )
[1384]828
[4648]829       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )  THEN
[1241]830
[4648]831          IF ( debug_output )  CALL  debug_message( 'initializing with constant profiles', 'start' )
[2259]832
833!
[4648]834!--       Use constructed initial profiles (velocity constant with height, temperature profile with
835!--       constant gradient)
[667]836          DO  i = nxlg, nxrg
837             DO  j = nysg, nyng
[1]838                pt(:,j,i) = pt_init
839                u(:,j,i)  = u_init
840                v(:,j,i)  = v_init
841             ENDDO
842          ENDDO
843!
[2758]844!--       Mask topography
[4346]845          u = MERGE( u, 0.0_wp, BTEST( wall_flags_total_0, 1 ) )
846          v = MERGE( v, 0.0_wp, BTEST( wall_flags_total_0, 2 ) )
[2758]847!
[4648]848!--       Set initial horizontal velocities at the lowest computational grid levels to zero in order
849!--       to avoid too small time steps caused by the diffusion limit in the initial phase of a run
850!--       (at k=1, dz/2 occurs in the limiting formula!).
851!--       Please note, in case land- or urban-surface model is used and a spinup is applied, masking
852!--       the lowest computational level is not possible as MOST as well as energy-balance
853!--       parametrizations will not work with zero wind velocity.
854          IF ( ibc_uv_b /= 1  .AND.  .NOT. spinup )  THEN
[1815]855             DO  i = nxlg, nxrg
856                DO  j = nysg, nyng
[2232]857                   DO  k = nzb, nzt
[4648]858                      u(k,j,i) = MERGE( u(k,j,i), 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 20 ) )
859                      v(k,j,i) = MERGE( v(k,j,i), 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 21 ) )
[2232]860                   ENDDO
[1815]861                ENDDO
862             ENDDO
863          ENDIF
[1]864
[1960]865          IF ( humidity )  THEN
[667]866             DO  i = nxlg, nxrg
867                DO  j = nysg, nyng
[1]868                   q(:,j,i) = q_init
869                ENDDO
870             ENDDO
871          ENDIF
[4648]872
[1960]873          IF ( passive_scalar )  THEN
874             DO  i = nxlg, nxrg
875                DO  j = nysg, nyng
876                   s(:,j,i) = s_init
877                ENDDO
878             ENDDO
879          ENDIF
[1]880
[1920]881!
[4648]882!--       Compute initial temperature field and other constants used in case of a sloping surface.
[1]883          IF ( sloping_surface )  CALL init_slope
[2618]884!
[4648]885!--       Initialize surface variables, e.g. friction velocity, momentum fluxes, etc.
[2618]886          CALL init_surfaces
[4648]887
[3987]888          IF ( debug_output )  CALL debug_message( 'initializing with constant profiles', 'end' )
[1384]889
[4648]890       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )  THEN
[1384]891
[3987]892          IF ( debug_output )  CALL debug_message( 'initializing by user', 'start' )
[46]893!
[4648]894!--       Pre-initialize surface variables, i.e. setting start- and end-indices at each
895!--       (j,i)-location. Please note, this does not supersede user-defined initialization of
896!--       surface quantities.
[2618]897          CALL init_surfaces
898!
[46]899!--       Initialization will completely be done by the user
900          CALL user_init_3d_model
901
[3987]902          IF ( debug_output )  CALL debug_message( 'initializing by user', 'end' )
[1384]903
[1]904       ENDIF
[1384]905
[4648]906       IF ( debug_output )  THEN
907          CALL debug_message( 'initializing statistics, boundary conditions, etc.', 'start' )
908       ENDIF
[1384]909
[667]910!
911!--    Bottom boundary
[4648]912       IF ( ibc_uv_b == 0  .OR.  ibc_uv_b == 2  )  THEN
[1340]913          u(nzb,:,:) = 0.0_wp
914          v(nzb,:,:) = 0.0_wp
[667]915       ENDIF
[1]916
917!
[151]918!--    Apply channel flow boundary condition
[132]919       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
[1340]920          u(nzt+1,:,:) = 0.0_wp
921          v(nzt+1,:,:) = 0.0_wp
[132]922       ENDIF
923
924!
[1]925!--    Calculate virtual potential temperature
[1960]926       IF ( humidity )  vpt = pt * ( 1.0_wp + 0.61_wp * q )
[1]927
928!
[4648]929!--    Store initial profiles for output purposes etc.. Please note, in case of initialization of u,
930!--    v, w, pt, and q via output data derived from larger scale models, data will not be
931!--    horizontally homogeneous. Actually, a mean profile should be calculated before.
[1]932       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
933       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
[667]934       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
[1340]935          hom(nzb,1,5,:) = 0.0_wp
936          hom(nzb,1,6,:) = 0.0_wp
[1]937       ENDIF
938       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
939
[75]940       IF ( humidity )  THEN
[1]941!
[4648]942!--       Store initial profile of total water content, virtual potential temperature
[1]943          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
944          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
[2696]945!
[4648]946!--       Store initial profile of mixing ratio and potential temperature
[3274]947          IF ( bulk_cloud_model  .OR.  cloud_droplets ) THEN
[1]948             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
949             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
950          ENDIF
951       ENDIF
952
[2696]953!
954!--    Store initial scalar profile
[1]955       IF ( passive_scalar )  THEN
[2513]956          hom(:,1,121,:) = SPREAD(  s(:,nys,nxl), 2, statistic_regions+1 )
[1]957       ENDIF
958
959!
[1400]960!--    Initialize the random number generators (from numerical recipes)
961       CALL random_function_ini
[4648]962
[1400]963       IF ( random_generator == 'random-parallel' )  THEN
[3241]964          CALL init_parallel_random_generator( nx, nys, nyn, nxl, nxr )
[1400]965       ENDIF
966!
[4648]967!--    Set the reference state to be used in the buoyancy terms (for ocean runs the reference state
968!--    will be set (overwritten) in init_ocean).
[1179]969       IF ( use_single_reference_value )  THEN
[4648]970          IF ( .NOT. humidity )  THEN
[1179]971             ref_state(:) = pt_reference
972          ELSE
973             ref_state(:) = vpt_reference
974          ENDIF
975       ELSE
[4648]976          IF ( .NOT. humidity )  THEN
[1179]977             ref_state(:) = pt_init(:)
978          ELSE
979             ref_state(:) = vpt(:,nys,nxl)
980          ENDIF
981       ENDIF
[152]982
983!
[707]984!--    For the moment, vertical velocity is zero
[1340]985       w = 0.0_wp
[1]986
987!
988!--    Initialize array sums (must be defined in first call of pres)
[1340]989       sums = 0.0_wp
[1]990
991!
[707]992!--    In case of iterative solvers, p must get an initial value
[1575]993       IF ( psolver(1:9) == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0_wp
[707]994!
[1]995!--    Impose vortex with vertical axis on the initial velocity profile
996       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
997          CALL init_rankine
998       ENDIF
999
1000!
[4648]1001!--    Impose temperature anomaly (advection test only) or warm air bubble close to surface.
1002       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0  .OR.                          &
[3035]1003            INDEX( initializing_actions, 'initialize_bubble' ) /= 0  )  THEN
[1]1004          CALL init_pt_anomaly
1005       ENDIF
[4648]1006
[1]1007!
1008!--    If required, change the surface temperature at the start of the 3D run
[1340]1009       IF ( pt_surface_initial_change /= 0.0_wp )  THEN
[1]1010          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
1011       ENDIF
1012
1013!
1014!--    If required, change the surface humidity/scalar at the start of the 3D
1015!--    run
[4648]1016       IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )                                  &
[1]1017          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
[4648]1018
1019       IF ( passive_scalar  .AND.  s_surface_initial_change /= 0.0_wp )                            &
[1960]1020          s(nzb,:,:) = s(nzb,:,:) + s_surface_initial_change
[1]1021
[4648]1022
[1]1023!
1024!--    Initialize old and new time levels.
[2696]1025       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1026       pt_p = pt; u_p = u; v_p = v; w_p = w
[1]1027
[1960]1028       IF ( humidity  )  THEN
[1340]1029          tq_m = 0.0_wp
[1]1030          q_p = q
1031       ENDIF
[4648]1032
[1960]1033       IF ( passive_scalar )  THEN
1034          ts_m = 0.0_wp
1035          s_p  = s
[4648]1036       ENDIF
[1]1037
[4648]1038       IF ( debug_output )  THEN
1039          CALL debug_message( 'initializing statistics, boundary conditions, etc.', 'end' )
1040       ENDIF
[94]1041
[1788]1042    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.         &
[2232]1043             TRIM( initializing_actions ) == 'cyclic_fill' )                   &
[1]1044    THEN
[1384]1045
[4648]1046       IF ( debug_output )  THEN
1047          CALL debug_message( 'initializing in case of restart / cyclic_fill', 'start' )
1048       ENDIF
[1]1049!
[4648]1050!--    Initialize surface elements and its attributes, e.g. heat- and momentumfluxes, roughness,
1051!--    scaling parameters. As number of surface elements might be different between runs, e.g. in
1052!--    case of cyclic fill, and not all surface elements are read, surface elements need to be
1053!--    initialized before.
1054!--    Please note, in case of cyclic fill, surfaces should be initialized after restart data is
1055!--    read, else, individual settings of surface parameters will be overwritten from data of
1056!--    precursor run, hence, init_surfaces is called a second time after reading the restart data.
1057       CALL init_surfaces
[3609]1058!
[4648]1059!--    When reading data for cyclic fill of 3D prerun data files, read some of the global variables
1060!--    from the restart file which are required for initializing the inflow
[328]1061       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[559]1062
[4783]1063!
1064!--       Blockwise I/O does not work together with MPI-I/O
1065          IF ( restart_data_format_input(1:3) == 'mpi' )  THEN
1066             CALL rrd_read_parts_of_global
1067          ELSE
1068             DO  i = 0, io_blocks-1
1069                IF ( i == io_group )  THEN
1070                   CALL rrd_read_parts_of_global
1071                ENDIF
1072#if defined( __parallel )
1073                CALL MPI_BARRIER( comm2d, ierr )
1074#endif
1075             ENDDO
1076          ENDIF
1077
1078       ENDIF
1079
1080!
1081!--    Read processor specific binary data from restart file.
1082!--    Blockwise I/O does not work together with MPI-I/O
1083       IF ( restart_data_format_input(1:3) == 'mpi' )  THEN
1084          CALL rrd_local
1085       ELSE
[759]1086          DO  i = 0, io_blocks-1
1087             IF ( i == io_group )  THEN
[4783]1088                CALL rrd_local
[759]1089             ENDIF
1090#if defined( __parallel )
1091             CALL MPI_BARRIER( comm2d, ierr )
1092#endif
1093          ENDDO
[767]1094       ENDIF
1095
[4648]1096
[4365]1097       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1098
[3608]1099!
[4648]1100!--       In case of cyclic fill, call init_surfaces a second time, so that surface properties such
1101!--       as heat fluxes are initialized as prescribed.
[3609]1102          CALL init_surfaces
[767]1103
[328]1104!
[4648]1105!--       Overwrite u_init, v_init, pt_init, q_init and s_init with the horizontally mean (hom)
1106!--       vertical profiles from the end of the prerun, because these profiles shall be used as the
1107!--       reference state for the rayleigh damping and the pt_damping. This is especially important
1108!--       for the use of large_scale_subsidence, because the reference temperature in the free
1109!--       atmosphere changes in time.
[4365]1110          u_init(:) = hom_sum(:,1,0)
1111          v_init(:) = hom_sum(:,2,0)
1112          pt_init(:) = hom_sum(:,4,0)
[4648]1113          IF ( humidity )  q_init(:) = hom_sum(:,41,0)
1114          IF ( passive_scalar )  s_init(:) = hom_sum(:,115,0)
[4365]1115       ENDIF
1116!
[4648]1117!--    In case of complex terrain and cyclic fill method as initialization, shift initial data in
1118!--    the vertical direction for each point in the x-y-plane depending on local surface height.
1119       IF ( complex_terrain  .AND.  TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
[2550]1120          DO  i = nxlg, nxrg
1121             DO  j = nysg, nyng
[4168]1122                nz_u_shift = topo_top_ind(j,i,1)
1123                nz_v_shift = topo_top_ind(j,i,2)
1124                nz_w_shift = topo_top_ind(j,i,3)
1125                nz_s_shift = topo_top_ind(j,i,0)
[2550]1126
[4648]1127                u(nz_u_shift:nzt+1,j,i)  = u(0:nzt+1-nz_u_shift,j,i)
[2550]1128
1129                v(nz_v_shift:nzt+1,j,i)  = v(0:nzt+1-nz_v_shift,j,i)
1130
1131                w(nz_w_shift:nzt+1,j,i)  = w(0:nzt+1-nz_w_shift,j,i)
1132
1133                p(nz_s_shift:nzt+1,j,i)  =  p(0:nzt+1-nz_s_shift,j,i)
1134                pt(nz_s_shift:nzt+1,j,i) = pt(0:nzt+1-nz_s_shift,j,i)
1135             ENDDO
1136          ENDDO
1137       ENDIF
1138!
[767]1139!--    Initialization of the turbulence recycling method
[4648]1140       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.  turbulent_inflow )  THEN
[767]1141!
1142!--       First store the profiles to be used at the inflow.
[4648]1143!--       These profiles are the (temporally) and horizontally averaged vertical profiles from the
1144!--       prerun. Alternatively, prescribed profiles for u,v-components can be used.
[3288]1145          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,1:num_mean_inflow_profiles) )
[151]1146
[767]1147          IF ( use_prescribed_profile_data )  THEN
1148             mean_inflow_profiles(:,1) = u_init            ! u
1149             mean_inflow_profiles(:,2) = v_init            ! v
1150          ELSE
[328]1151             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
1152             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
[767]1153          ENDIF
1154          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
[4648]1155          IF ( humidity )  mean_inflow_profiles(:,6) = hom_sum(:,41,0)          ! q
1156          IF ( passive_scalar )  mean_inflow_profiles(:,7) = hom_sum(:,115,0)   ! s
[4365]1157
[2550]1158!
[4648]1159!--       In case of complex terrain, determine vertical displacement at inflow boundary and adjust
1160!--       mean inflow profiles
[2550]1161          IF ( complex_terrain )  THEN
[4648]1162             IF ( nxlg <= 0  .AND.  nxrg >= 0  .AND.  nysg <= 0  .AND.  nyng >= 0 )  THEN
[4168]1163                nz_u_shift_l = topo_top_ind(j,i,1)
1164                nz_v_shift_l = topo_top_ind(j,i,2)
1165                nz_w_shift_l = topo_top_ind(j,i,3)
1166                nz_s_shift_l = topo_top_ind(j,i,0)
[2550]1167             ELSE
1168                nz_u_shift_l = 0
1169                nz_v_shift_l = 0
1170                nz_w_shift_l = 0
1171                nz_s_shift_l = 0
1172             ENDIF
[151]1173
[2550]1174#if defined( __parallel )
[4648]1175             CALL MPI_ALLREDUCE( nz_u_shift_l, nz_u_shift, 1, MPI_INTEGER, MPI_MAX, comm2d, ierr )
1176             CALL MPI_ALLREDUCE( nz_v_shift_l, nz_v_shift, 1, MPI_INTEGER, MPI_MAX, comm2d, ierr )
1177             CALL MPI_ALLREDUCE( nz_w_shift_l, nz_w_shift, 1, MPI_INTEGER, MPI_MAX, comm2d, ierr )
1178             CALL MPI_ALLREDUCE( nz_s_shift_l, nz_s_shift, 1, MPI_INTEGER, MPI_MAX, comm2d, ierr )
[2550]1179#else
1180             nz_u_shift = nz_u_shift_l
1181             nz_v_shift = nz_v_shift_l
1182             nz_w_shift = nz_w_shift_l
1183             nz_s_shift = nz_s_shift_l
1184#endif
1185
1186             mean_inflow_profiles(:,1) = 0.0_wp
1187             mean_inflow_profiles(nz_u_shift:nzt+1,1) = hom_sum(0:nzt+1-nz_u_shift,1,0)  ! u
1188
1189             mean_inflow_profiles(:,2) = 0.0_wp
1190             mean_inflow_profiles(nz_v_shift:nzt+1,2) = hom_sum(0:nzt+1-nz_v_shift,2,0)  ! v
1191
1192             mean_inflow_profiles(nz_s_shift:nzt+1,4) = hom_sum(0:nzt+1-nz_s_shift,4,0)  ! pt
1193
1194          ENDIF
1195
[151]1196!
[4648]1197!--       If necessary, adjust the horizontal flow field to the prescribed profiles
[767]1198          IF ( use_prescribed_profile_data )  THEN
1199             DO  i = nxlg, nxrg
[667]1200                DO  j = nysg, nyng
[328]1201                   DO  k = nzb, nzt+1
[767]1202                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1203                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
[328]1204                   ENDDO
[151]1205                ENDDO
[767]1206             ENDDO
1207          ENDIF
[151]1208
1209!
[4648]1210!--       Use these mean profiles at the inflow (provided that Dirichlet conditions are used)
[3182]1211          IF ( bc_dirichlet_l )  THEN
[767]1212             DO  j = nysg, nyng
1213                DO  k = nzb, nzt+1
1214                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1215                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
[1340]1216                   w(k,j,nxlg:-1)  = 0.0_wp
[767]1217                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
[4648]1218                   IF ( humidity )  q(k,j,nxlg:-1)  = mean_inflow_profiles(k,6)
1219                   IF ( passive_scalar )  s(k,j,nxlg:-1)  = mean_inflow_profiles(k,7)
[767]1220                ENDDO
1221             ENDDO
1222          ENDIF
1223
[151]1224!
[4648]1225!--       Calculate the damping factors to be used at the inflow. For a turbulent inflow the
1226!--       turbulent fluctuations have to be limited vertically because otherwise the turbulent
1227!--       inflow layer will grow in time.
[1340]1228          IF ( inflow_damping_height == 9999999.9_wp )  THEN
[767]1229!
[4648]1230!--          Default: use the inversion height calculated by the prerun; if this is zero,
1231!--          inflow_damping_height must be explicitly specified.
[1340]1232             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0_wp )  THEN
[767]1233                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1234             ELSE
[4648]1235                WRITE( message_string, * ) 'inflow_damping_height must be ',                       &
1236                                           'explicitly specified because&the inversion height ',   &
1237                                           'calculated by the prerun is zero.'
[767]1238                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
[292]1239             ENDIF
[151]1240
[767]1241          ENDIF
1242
[1340]1243          IF ( inflow_damping_width == 9999999.9_wp )  THEN
[151]1244!
[4648]1245!--          Default for the transition range: one tenth of the undamped layer
[1340]1246             inflow_damping_width = 0.1_wp * inflow_damping_height
[151]1247
[767]1248          ENDIF
[151]1249
[767]1250          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
[151]1251
[767]1252          DO  k = nzb, nzt+1
[151]1253
[767]1254             IF ( zu(k) <= inflow_damping_height )  THEN
[1340]1255                inflow_damping_factor(k) = 1.0_wp
[996]1256             ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) )  THEN
[4648]1257                inflow_damping_factor(k) = 1.0_wp -                                                &
1258                                           ( zu(k) - inflow_damping_height ) / inflow_damping_width
[767]1259             ELSE
[1340]1260                inflow_damping_factor(k) = 0.0_wp
[767]1261             ENDIF
[151]1262
[767]1263          ENDDO
[151]1264
[147]1265       ENDIF
1266
[152]1267!
[2696]1268!--    Inside buildings set velocities back to zero
[4648]1269       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.  topography /= 'flat' )  THEN
[359]1270!
[2696]1271!--       Inside buildings set velocities back to zero.
1272!--       Other scalars (pt, q, s, p, sa, ...) are ignored at present,
[359]1273!--       maybe revise later.
[1001]1274          DO  i = nxlg, nxrg
1275             DO  j = nysg, nyng
[2232]1276                DO  k = nzb, nzt
[4648]1277                   u(k,j,i) = MERGE( u(k,j,i), 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 1 ) )
1278                   v(k,j,i) = MERGE( v(k,j,i), 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 2 ) )
1279                   w(k,j,i) = MERGE( w(k,j,i), 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 3 ) )
[2232]1280                ENDDO
[359]1281             ENDDO
[1001]1282          ENDDO
[359]1283
1284       ENDIF
1285
1286!
[4648]1287!--    Calculate initial temperature field and other constants used in case of a sloping surface
[1]1288       IF ( sloping_surface )  CALL init_slope
1289
1290!
[4648]1291!--    Initialize new time levels (only done in order to set boundary values including ghost points)
[2696]1292       pt_p = pt; u_p = u; v_p = v; w_p = w
[1960]1293       IF ( humidity )  THEN
[1053]1294          q_p = q
1295       ENDIF
[1960]1296       IF ( passive_scalar )  s_p  = s
[181]1297!
[4648]1298!--    Allthough tendency arrays are set in prognostic_equations, they have have to be predefined
1299!--    here because they are used (but multiplied with 0) there before they are set.
[2696]1300       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
[1960]1301       IF ( humidity )  THEN
[1340]1302          tq_m = 0.0_wp
[1053]1303       ENDIF
[1960]1304       IF ( passive_scalar )  ts_m  = 0.0_wp
[181]1305
[4648]1306       IF ( debug_output )  THEN
1307          CALL debug_message( 'initializing in case of restart / cyclic_fill', 'end' )
1308       ENDIF
[1384]1309
[1]1310    ELSE
1311!
1312!--    Actually this part of the programm should not be reached
[254]1313       message_string = 'unknown initializing problem'
1314       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
[1]1315    ENDIF
1316
1317!
[667]1318!-- Calculate the initial volume flow at the right and north boundary
[709]1319    IF ( conserve_volume_flow )  THEN
[151]1320
[767]1321       IF ( use_prescribed_profile_data )  THEN
[667]1322
[1340]1323          volume_flow_initial_l = 0.0_wp
1324          volume_flow_area_l    = 0.0_wp
[732]1325
[667]1326          IF ( nxr == nx )  THEN
1327             DO  j = nys, nyn
[2232]1328                DO  k = nzb+1, nzt
[4648]1329                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +                           &
1330                                              u_init(k) * dzw(k)                                   &
1331                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1332                                                       BTEST( wall_flags_total_0(k,j,nxr), 1 )     &
1333                                                     )
[2232]1334
[4648]1335                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)                       &
1336                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1337                                                       BTEST( wall_flags_total_0(k,j,nxr), 1 )     &
1338                                                     )
[767]1339                ENDDO
1340             ENDDO
1341          ENDIF
[4648]1342
[767]1343          IF ( nyn == ny )  THEN
1344             DO  i = nxl, nxr
[2232]1345                DO  k = nzb+1, nzt
[4648]1346                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +                           &
1347                                              v_init(k) * dzw(k)                                   &
1348                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1349                                                       BTEST( wall_flags_total_0(k,nyn,i), 2 )     &
1350                                                     )
1351                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)                       &
1352                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1353                                                       BTEST( wall_flags_total_0(k,nyn,i), 2 )     &
1354                                                     )
[767]1355                ENDDO
1356             ENDDO
1357          ENDIF
1358
1359#if defined( __parallel )
[4648]1360          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1), 2, MPI_REAL,       &
1361                              MPI_SUM, comm2d, ierr )
1362          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), 2, MPI_REAL, MPI_SUM,    &
1363                              comm2d, ierr )
[767]1364
1365#else
1366          volume_flow_initial = volume_flow_initial_l
1367          volume_flow_area    = volume_flow_area_l
[4648]1368#endif
[767]1369
1370       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1371
[1340]1372          volume_flow_initial_l = 0.0_wp
1373          volume_flow_area_l    = 0.0_wp
[767]1374
1375          IF ( nxr == nx )  THEN
1376             DO  j = nys, nyn
[2232]1377                DO  k = nzb+1, nzt
[4648]1378                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +                           &
1379                                              hom_sum(k,1,0) * dzw(k)                              &
1380                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1381                                                       BTEST( wall_flags_total_0(k,j,nx), 1 )      &
1382                                                     )
1383                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)                       &
1384                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1385                                                       BTEST( wall_flags_total_0(k,j,nx), 1 )      &
1386                                                     )
[667]1387                ENDDO
1388             ENDDO
1389          ENDIF
[4648]1390
[667]1391          IF ( nyn == ny )  THEN
1392             DO  i = nxl, nxr
[2232]1393                DO  k = nzb+1, nzt
[4648]1394                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +                           &
1395                                              hom_sum(k,2,0) * dzw(k)                              &
1396                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1397                                                       BTEST( wall_flags_total_0(k,ny,i), 2 )      &
1398                                                     )
1399                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)                       &
1400                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1401                                                       BTEST( wall_flags_total_0(k,ny,i), 2 )      &
1402                                                     )
[667]1403                ENDDO
1404             ENDDO
1405          ENDIF
1406
[732]1407#if defined( __parallel )
[4648]1408          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1), 2, MPI_REAL,       &
1409                              MPI_SUM, comm2d, ierr )
1410          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), 2, MPI_REAL, MPI_SUM,    &
1411                              comm2d, ierr )
[732]1412
1413#else
1414          volume_flow_initial = volume_flow_initial_l
1415          volume_flow_area    = volume_flow_area_l
[4648]1416#endif
[732]1417
[667]1418       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1419
[1340]1420          volume_flow_initial_l = 0.0_wp
1421          volume_flow_area_l    = 0.0_wp
[732]1422
[667]1423          IF ( nxr == nx )  THEN
1424             DO  j = nys, nyn
[2232]1425                DO  k = nzb+1, nzt
[4648]1426                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +                           &
1427                                              u(k,j,nx) * dzw(k)                                   &
1428                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1429                                                       BTEST( wall_flags_total_0(k,j,nx), 1 )      &
1430                                                     )
1431                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)                       &
1432                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1433                                                   BTEST( wall_flags_total_0(k,j,nx), 1 )          &
1434                                                     )
[667]1435                ENDDO
1436             ENDDO
1437          ENDIF
[4648]1438
[667]1439          IF ( nyn == ny )  THEN
1440             DO  i = nxl, nxr
[2232]1441                DO  k = nzb+1, nzt
[4648]1442                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +                           &
1443                                              v(k,ny,i) * dzw(k)                                   &
1444                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1445                                                       BTEST( wall_flags_total_0(k,ny,i), 2 )      &
1446                                                     )
1447                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)                       &
1448                                              * MERGE( 1.0_wp, 0.0_wp,                             &
1449                                                       BTEST( wall_flags_total_0(k,ny,i), 2 )      &
1450                                                     )
[667]1451                ENDDO
1452             ENDDO
1453          ENDIF
1454
1455#if defined( __parallel )
[4648]1456          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1), 2, MPI_REAL,       &
1457                              MPI_SUM, comm2d, ierr )
1458          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1), 2, MPI_REAL, MPI_SUM,    &
1459                              comm2d, ierr )
[667]1460
1461#else
[732]1462          volume_flow_initial = volume_flow_initial_l
1463          volume_flow_area    = volume_flow_area_l
[4648]1464#endif
[667]1465
[732]1466       ENDIF
1467
[151]1468!
[4648]1469!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated from u|v_bulk instead
[680]1470       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
1471          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
1472          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
1473       ENDIF
[667]1474
[680]1475    ENDIF
[2232]1476!
[4648]1477!-- In the following, surface properties can be further initialized with input from static driver
1478!-- file.
1479!-- At the moment this affects only default surfaces. For example, roughness length or sensible /
1480!-- latent heat fluxes can be initialized heterogeneously for default surfaces. Therefore, a generic
1481!-- routine from netcdf_data_input_mod is called to read a 2D array.
[4150]1482    IF ( input_pids_static )  THEN
1483!
[4151]1484!--    Allocate memory for possible static input
1485       ALLOCATE( tmp_2d%var(nys:nyn,nxl:nxr) )
1486       tmp_2d%var = 0.0_wp
1487!
[4150]1488!--    Open the static input file
[4151]1489#if defined( __netcdf )
[4648]1490       CALL open_read_file( TRIM( input_file_static ) //                                           &
1491                            TRIM( coupling_char ), pids_id )
1492
[4186]1493       CALL inquire_num_variables( pids_id, num_var_pids )
[4150]1494!
1495!--    Allocate memory to store variable names and read them
[4186]1496       ALLOCATE( vars_pids(1:num_var_pids) )
1497       CALL inquire_variable_names( pids_id, vars_pids )
[4150]1498!
[4648]1499!--    Input roughness length.
[4186]1500       IF ( check_existence( vars_pids, 'z0' ) )  THEN
[4150]1501!
1502!--       Read _FillValue attribute
[4648]1503          CALL get_attribute( pids_id, char_fill, tmp_2d%fill, .FALSE., 'z0' )
1504!
1505!--       Read variable
1506          CALL get_variable( pids_id, 'z0', tmp_2d%var, nxl, nxr, nys, nyn )
1507!
1508!--       Initialize roughness length. Note, z0 will be only initialized at default-type surfaces.
1509!--       At natural or urban z0 is implicitly initialized by the respective parameter lists.
1510!--       Initialize horizontal surface elements.
1511          CALL init_single_surface_properties( surf_def_h(0)%z0, tmp_2d%var, surf_def_h(0)%ns,     &
1512                                               tmp_2d%fill, surf_def_h(0)%i, surf_def_h(0)%j )
1513!
1514!--       Initialize roughness also at vertical surface elements.
1515!--       Note, the actual 2D input arrays are only defined on the subdomain. Therefore, pass the
1516!--       index arrays with their respective offset values.
1517          DO  l = 0, 3
1518             CALL init_single_surface_properties( surf_def_v(l)%z0, tmp_2d%var, surf_def_v(l)%ns,  &
1519                                                  tmp_2d%fill, surf_def_v(l)%i+surf_def_v(l)%ioff, &
1520                                                  surf_def_v(l)%j+surf_def_v(l)%joff )
[4150]1521          ENDDO
[4648]1522
[4150]1523       ENDIF
1524!
[4648]1525!--    Input surface sensible heat flux.
[4514]1526       IF ( check_existence( vars_pids, 'shf' ) )  THEN
1527!
1528!--       Read _FillValue attribute
[4648]1529          CALL get_attribute( pids_id, char_fill, tmp_2d%fill, .FALSE., 'shf' )
[4514]1530!
1531!--       Read variable
[4648]1532          CALL get_variable( pids_id, 'shf', tmp_2d%var, nxl, nxr, nys, nyn )
[4514]1533!
[4648]1534!--       Initialize heat flux. Note, shf will be only initialized at default-type surfaces. At
1535!--       natural or urban shf is implicitly initialized by the respective parameter lists.
[4514]1536!--       Initialize horizontal surface elements.
[4648]1537          CALL init_single_surface_properties( surf_def_h(0)%shf, tmp_2d%var, surf_def_h(0)%ns,    &
1538                                               tmp_2d%fill, surf_def_h(0)%i, surf_def_h(0)%j )
[4514]1539!
1540!--       Initialize heat flux also at vertical surface elements.
[4648]1541!--       Note, the actual 2D input arrays are only defined on the subdomain. Therefore, pass the
1542!--       index arrays with their respective offset values.
[4514]1543          DO  l = 0, 3
[4648]1544             CALL init_single_surface_properties( surf_def_v(l)%shf, tmp_2d%var, surf_def_v(l)%ns, &
1545                                                  tmp_2d%fill, surf_def_v(l)%i+surf_def_v(l)%ioff, &
1546                                                  surf_def_v(l)%j+surf_def_v(l)%joff )
[4514]1547          ENDDO
1548
1549       ENDIF
1550!
[4791]1551!--    Input surface latent heat flux.
1552       IF ( humidity )  THEN
1553          IF ( check_existence( vars_pids, 'qsws' ) )  THEN
[4514]1554!
[4791]1555!--          Read _FillValue attribute
1556             CALL get_attribute( pids_id, char_fill, tmp_2d%fill, .FALSE., 'qsws' )
[4514]1557!
[4791]1558!--          Read variable
1559             CALL get_variable( pids_id, 'qsws', tmp_2d%var, nxl, nxr, nys, nyn )
[4514]1560!
[4791]1561!--          Initialize latent heat flux. Note, qsws will be only initialized at default-type surfaces.
1562!--          At natural or urban qsws is implicitly initialized by the respective parameter lists.
1563!--          Initialize horizontal surface elements.
1564             CALL init_single_surface_properties( surf_def_h(0)%qsws, tmp_2d%var, surf_def_h(0)%ns,&
1565                                                  tmp_2d%fill, surf_def_h(0)%i, surf_def_h(0)%j )
[4514]1566!
[4791]1567!--          Initialize latent heat flux also at vertical surface elements.
1568!--          Note, the actual 2D input arrays are only defined on the subdomain. Therefore, pass the
1569!--          index arrays with their respective offset values.
1570             DO  l = 0, 3
1571                CALL init_single_surface_properties( surf_def_v(l)%qsws, tmp_2d%var,               &
1572                                                     surf_def_v(l)%ns, tmp_2d%fill,                &
1573                                                     surf_def_v(l)%i+surf_def_v(l)%ioff,           &
1574                                                     surf_def_v(l)%j+surf_def_v(l)%joff )
1575             ENDDO
[4514]1576
[4791]1577          ENDIF
[4514]1578       ENDIF
1579!
[4791]1580!--    Input passive scalar flux.
1581       IF ( passive_scalar )  THEN
1582          IF ( check_existence( vars_pids, 'ssws' ) )  THEN
1583!
1584!--          Read _FillValue attribute
1585             CALL get_attribute( pids_id, char_fill, tmp_2d%fill, .FALSE., 'ssws' )
1586!
1587!--          Read variable
1588             CALL get_variable( pids_id, 'ssws', tmp_2d%var, nxl, nxr, nys, nyn )
1589!
1590!--          Initialize passive scalar flux. Initialize horizontal surface elements.
1591             CALL init_single_surface_properties( surf_def_h(0)%ssws, tmp_2d%var, surf_def_h(0)%ns,&
1592                                                  tmp_2d%fill, surf_def_h(0)%i, surf_def_h(0)%j )
1593
1594             DO  l = 0, 3
1595                CALL init_single_surface_properties( surf_def_v(l)%ssws, tmp_2d%var,               &
1596                                                     surf_def_v(l)%ns, tmp_2d%fill,                &
1597                                                     surf_def_v(l)%i+surf_def_v(l)%ioff,           &
1598                                                     surf_def_v(l)%j+surf_def_v(l)%joff )
1599             ENDDO
1600
1601          ENDIF
1602       ENDIF
1603!
[4648]1604!--    Additional variables, can be initialized the
[4150]1605!--    same way.
[4514]1606
[4150]1607!
[4187]1608!--    Finally, close the input file and deallocate temporary arrays
1609       DEALLOCATE( vars_pids )
[4648]1610
[4186]1611       CALL close_input_file( pids_id )
[4151]1612#endif
1613       DEALLOCATE( tmp_2d%var )
[4150]1614    ENDIF
1615!
[4648]1616!-- Finally, if random_heatflux is set, disturb shf at horizontal surfaces. Actually, this should be
1617!-- done in surface_mod, where all other initializations of surface quantities are done. However,
1618!-- this would create a ring dependency, hence, it is done here. Maybe delete disturb_heatflux and
1619!-- tranfer the respective code directly into the initialization in surface_mod.
1620    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.                                &
[2232]1621         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
[4648]1622
1623       IF ( use_surface_fluxes  .AND.  constant_heatflux  .AND.  random_heatflux )  THEN
[2232]1624          IF ( surf_def_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_def_h(0) )
[4671]1625          IF ( surf_lsm_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_lsm_h(0) )
1626          IF ( surf_usm_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_usm_h(0) )
[2232]1627       ENDIF
1628    ENDIF
[680]1629
[787]1630!
[4648]1631!-- Compute total sum of grid points and the mean surface level height for each statistic region.
1632!-- These are mainly used for horizontal averaging of turbulence statistics.
1633!-- ngp_2dh: number of grid points of a horizontal cross section through the respective statistic
1634!--          region
[3747]1635!-- ngp_3d:  number of grid points of the respective statistic region
[2696]1636    ngp_2dh_outer_l   = 0
1637    ngp_2dh_outer     = 0
1638    ngp_2dh_s_inner_l = 0
1639    ngp_2dh_s_inner   = 0
1640    ngp_2dh_l         = 0
1641    ngp_2dh           = 0
1642    ngp_3d_inner_l    = 0.0_wp
1643    ngp_3d_inner      = 0
1644    ngp_3d            = 0
1645    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
1646
1647    mean_surface_level_height   = 0.0_wp
1648    mean_surface_level_height_l = 0.0_wp
1649!
1650!-- To do: New concept for these non-topography grid points!
1651    DO  sr = 0, statistic_regions
1652       DO  i = nxl, nxr
1653          DO  j = nys, nyn
1654             IF ( rmask(j,i,sr) == 1.0_wp )  THEN
1655!
1656!--             All xy-grid points
1657                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
1658!
[4648]1659!--             Determine mean surface-level height. In case of downward-facing walls are present,
1660!--             more than one surface level exist.
1661!--             In this case, use the lowest surface-level height.
1662                IF ( surf_def_h(0)%start_index(j,i) <= surf_def_h(0)%end_index(j,i) )  THEN
[2696]1663                   m = surf_def_h(0)%start_index(j,i)
1664                   k = surf_def_h(0)%k(m)
[4648]1665                   mean_surface_level_height_l(sr) = mean_surface_level_height_l(sr) + zw(k-1)
[2696]1666                ENDIF
[4671]1667                IF ( surf_lsm_h(0)%start_index(j,i) <= surf_lsm_h(0)%end_index(j,i) )  THEN
1668                   m = surf_lsm_h(0)%start_index(j,i)
1669                   k = surf_lsm_h(0)%k(m)
[4648]1670                   mean_surface_level_height_l(sr) = mean_surface_level_height_l(sr) + zw(k-1)
[2696]1671                ENDIF
[4671]1672                IF ( surf_usm_h(0)%start_index(j,i) <= surf_usm_h(0)%end_index(j,i) )  THEN
1673                   m = surf_usm_h(0)%start_index(j,i)
1674                   k = surf_usm_h(0)%k(m)
[4648]1675                   mean_surface_level_height_l(sr) = mean_surface_level_height_l(sr) + zw(k-1)
[2696]1676                ENDIF
1677
1678                k_surf = k - 1
1679
1680                DO  k = nzb, nzt+1
1681!
1682!--                xy-grid points above topography
[4648]1683                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr)     +                             &
1684                                           MERGE( 1, 0, BTEST( wall_flags_total_0(k,j,i), 24 ) )
[2696]1685
[4648]1686                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) +                             &
1687                                             MERGE( 1, 0, BTEST( wall_flags_total_0(k,j,i), 22 ) )
[2696]1688
1689                ENDDO
1690!
1691!--             All grid points of the total domain above topography
1692                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + ( nz - k_surf + 2 )
1693
1694
1695
1696             ENDIF
1697          ENDDO
1698       ENDDO
1699    ENDDO
[3747]1700
[2696]1701    sr = statistic_regions + 1
1702#if defined( __parallel )
1703    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[4648]1704    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
1705    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1706    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr, MPI_INTEGER, MPI_SUM, &
[2696]1707                        comm2d, ierr )
1708    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[4648]1709    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0), (nz+2)*sr, MPI_INTEGER,      &
1710                        MPI_SUM, comm2d, ierr )
[2696]1711    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[4648]1712    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL, MPI_SUM, comm2d,     &
1713                        ierr )
[2696]1714    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
1715    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[4648]1716    CALL MPI_ALLREDUCE( mean_surface_level_height_l(0), mean_surface_level_height(0), sr, MPI_REAL,&
[2696]1717                        MPI_SUM, comm2d, ierr )
1718    mean_surface_level_height = mean_surface_level_height / REAL( ngp_2dh )
1719#else
1720    ngp_2dh         = ngp_2dh_l
1721    ngp_2dh_outer   = ngp_2dh_outer_l
1722    ngp_2dh_s_inner = ngp_2dh_s_inner_l
1723    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
1724    mean_surface_level_height = mean_surface_level_height_l / REAL( ngp_2dh_l )
1725#endif
1726
[4648]1727    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) *                                     &
[2696]1728             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
1729
1730!
[4648]1731!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics, buoyancy, etc. A
1732!-- zero value will occur for cases where all grid points of the respective subdomain lie below the
1733!-- surface topography
1734    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   )
1735    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )), ngp_3d_inner(:) )
1736    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) )
[2696]1737
[4648]1738    DEALLOCATE( mean_surface_level_height_l, ngp_2dh_l, ngp_2dh_outer_l, ngp_3d_inner_l,           &
1739                ngp_3d_inner_tmp )
[4548]1740
[2696]1741!
[4648]1742!-- Initializae 3D offline nesting in COSMO model and read data from
[3347]1743!-- external NetCDF file.
1744    IF ( nesting_offline )  CALL nesting_offl_init
1745!
[787]1746!-- Initialize quantities for special advections schemes
1747    CALL init_advec
[680]1748
[667]1749!
[680]1750!-- Impose random perturbation on the horizontal velocity field and then
1751!-- remove the divergences from the velocity field at the initial stage
[4648]1752    IF ( create_disturbances  .AND.  disturbance_energy_limit /= 0.0_wp  .AND.                     &
1753         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.                                &
[680]1754         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
1755
[4648]1756       IF ( debug_output )  THEN
1757          CALL debug_message( 'creating disturbances + applying pressure solver', 'start' )
1758       ENDIF
[3849]1759!
1760!--    Needed for both disturb_field and pres
1761!$ACC DATA &
1762!$ACC CREATE(tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
1763!$ACC COPY(u(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
1764!$ACC COPY(v(nzb:nzt+1,nysg:nyng,nxlg:nxrg))
1765
[2232]1766       CALL disturb_field( 'u', tend, u )
1767       CALL disturb_field( 'v', tend, v )
[1384]1768
[3849]1769!$ACC DATA &
1770!$ACC CREATE(d(nzb+1:nzt,nys:nyn,nxl:nxr)) &
1771!$ACC COPY(w(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
1772!$ACC COPY(p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
1773!$ACC COPYIN(rho_air(nzb:nzt+1), rho_air_zw(nzb:nzt+1)) &
1774!$ACC COPYIN(ddzu(1:nzt+1), ddzw(1:nzt+1)) &
[4346]1775!$ACC COPYIN(wall_flags_total_0(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
[3849]1776!$ACC COPYIN(bc_h(0:1)) &
1777!$ACC COPYIN(bc_h(0)%i(1:bc_h(0)%ns)) &
1778!$ACC COPYIN(bc_h(0)%j(1:bc_h(0)%ns)) &
1779!$ACC COPYIN(bc_h(0)%k(1:bc_h(0)%ns)) &
1780!$ACC COPYIN(bc_h(1)%i(1:bc_h(1)%ns)) &
1781!$ACC COPYIN(bc_h(1)%j(1:bc_h(1)%ns)) &
1782!$ACC COPYIN(bc_h(1)%k(1:bc_h(1)%ns))
1783
[680]1784       n_sor = nsor_ini
1785       CALL pres
1786       n_sor = nsor
[1384]1787
[3849]1788!$ACC END DATA
1789!$ACC END DATA
1790
[4648]1791       IF ( debug_output )  THEN
1792          CALL debug_message( 'creating disturbances + applying pressure solver', 'end' )
1793       ENDIF
[3987]1794
[680]1795    ENDIF
1796
[3294]1797    IF ( .NOT. ocean_mode )  THEN
[3274]1798
1799       ALLOCATE( hyp(nzb:nzt+1) )
1800       ALLOCATE( d_exner(nzb:nzt+1) )
1801       ALLOCATE( exner(nzb:nzt+1) )
1802       ALLOCATE( hyrho(nzb:nzt+1) )
[1849]1803!
[3274]1804!--    Check temperature in case of too large domain height
1805       DO  k = nzb, nzt+1
[4648]1806          IF ( ( pt_surface * exner_function( surface_pressure * 100.0_wp ) - g/c_p * zu(k) )      &
1807                 < 0.0_wp )  THEN
1808             WRITE( message_string, * )  'absolute temperature < 0.0 at zu(', k, ') = ', zu(k)
[3685]1809             CALL message( 'init_3d_model', 'PA0142', 1, 2, 0, 6, 0 )
[3274]1810          ENDIF
1811       ENDDO
1812
1813!
1814!--    Calculate vertical profile of the hydrostatic pressure (hyp)
[4648]1815       hyp    = barometric_formula( zu, pt_surface * exner_function( surface_pressure * 100.0_wp ),&
1816                                    surface_pressure * 100.0_wp )
1817       d_exner = exner_function_invers( hyp )
1818       exner = 1.0_wp / exner_function_invers( hyp )
1819       hyrho  = ideal_gas_law_rho_pt( hyp, pt_init )
[3274]1820!
1821!--    Compute reference density
[4648]1822       rho_surface = ideal_gas_law_rho( surface_pressure * 100.0_wp,                               &
1823                                        pt_surface * exner_function( surface_pressure * 100.0_wp ) )
[3274]1824
[96]1825    ENDIF
[1]1826
1827!
1828!-- If required, initialize particles
[3159]1829    IF ( agents_active )  CALL mas_init
1830!
[3937]1831!-- Initialization of synthetic turbulence generator
1832    IF ( use_syn_turb_gen )  CALL stg_init
[2696]1833!
[3685]1834!-- Initializing actions for all other modules
1835    CALL module_interface_init
[2696]1836!
[4648]1837!-- Initialize surface layer (done after LSM as roughness length are required for initialization
[3685]1838    IF ( constant_flux_layer )  CALL init_surface_layer_fluxes
[2977]1839!
[3421]1840!-- Initialize surface data output
[3685]1841    IF ( surface_output )  CALL surface_data_output_init
[3472]1842!
[4648]1843!-- Initialize the ws-scheme.
1844    IF ( ws_scheme_sca  .OR.  ws_scheme_mom )  CALL ws_init
[3711]1845!
1846!-- Perform post-initializing checks for all other modules
1847    CALL module_interface_init_checks
[1]1848
1849!
[4648]1850!-- Initialize surface forcing corresponding to large-scale forcing. Therein,
[4548]1851!-- initialize heat-fluxes, etc. via datatype. Revise it later!
[4648]1852    IF ( large_scale_forcing  .AND.  lsf_surf )  THEN
[4548]1853       IF ( use_surface_fluxes  .AND.  constant_heatflux )  THEN
[4648]1854          CALL ls_forcing_surf( simulated_time )
[4548]1855       ENDIF
1856    ENDIF
1857!
[4648]1858!-- Setting weighting factors for calculation of perturbation pressure and turbulent quantities from
1859!-- the RK substeps.
1860    IF ( TRIM( timestep_scheme ) == 'runge-kutta-3' )  THEN      ! for RK3-method
[709]1861
[1322]1862       weight_substep(1) = 1._wp/6._wp
1863       weight_substep(2) = 3._wp/10._wp
1864       weight_substep(3) = 8._wp/15._wp
[709]1865
[1322]1866       weight_pres(1)    = 1._wp/3._wp
1867       weight_pres(2)    = 5._wp/12._wp
1868       weight_pres(3)    = 1._wp/4._wp
[709]1869
[4648]1870    ELSEIF ( TRIM( timestep_scheme ) == 'runge-kutta-2' )  THEN  ! for RK2-method
[709]1871
[1322]1872       weight_substep(1) = 1._wp/2._wp
1873       weight_substep(2) = 1._wp/2._wp
[4648]1874
[1322]1875       weight_pres(1)    = 1._wp/2._wp
[4648]1876       weight_pres(2)    = 1._wp/2._wp
[709]1877
[1001]1878    ELSE                                     ! for Euler-method
[709]1879
[4648]1880       weight_substep(1) = 1.0_wp
1881       weight_pres(1)    = 1.0_wp
[709]1882
[673]1883    ENDIF
1884
1885!
[1]1886!-- Initialize Rayleigh damping factors
[1340]1887    rdf    = 0.0_wp
1888    rdf_sc = 0.0_wp
1889    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
[3294]1890
1891       IF (  .NOT.  ocean_mode )  THEN
[108]1892          DO  k = nzb+1, nzt
1893             IF ( zu(k) >= rayleigh_damping_height )  THEN
[4648]1894                rdf(k) = rayleigh_damping_factor *                                                 &
1895                         ( SIN( pi * 0.5_wp * ( zu(k) - rayleigh_damping_height )                  &
1896                                / ( zu(nzt) - rayleigh_damping_height ) )                          &
1897                         )**2
[108]1898             ENDIF
1899          ENDDO
1900       ELSE
[3294]1901!
[4648]1902!--       In ocean mode, rayleigh damping is applied in the lower part of the model domain
[108]1903          DO  k = nzt, nzb+1, -1
1904             IF ( zu(k) <= rayleigh_damping_height )  THEN
[4648]1905                rdf(k) = rayleigh_damping_factor *                                                 &
1906                         ( SIN( pi * 0.5_wp * ( rayleigh_damping_height - zu(k) )                  &
1907                                / ( rayleigh_damping_height - zu(nzb+1) ) )                        &
1908                         )**2
[108]1909             ENDIF
1910          ENDDO
1911       ENDIF
[3294]1912
[1]1913    ENDIF
[785]1914    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
[1]1915
1916!
[4648]1917!-- Initialize the starting level and the vertical smoothing factor used for the external pressure
1918!-- gradient
[1340]1919    dp_smooth_factor = 1.0_wp
[240]1920    IF ( dp_external )  THEN
1921!
[4648]1922!--    Set the starting level dp_level_ind_b only if it has not been set before (e.g. in init_grid).
[240]1923       IF ( dp_level_ind_b == 0 )  THEN
1924          ind_array = MINLOC( ABS( dp_level_b - zu ) )
[4648]1925          dp_level_ind_b = ind_array(1) - 1 + nzb
[240]1926                                        ! MINLOC uses lower array bound 1
1927       ENDIF
1928       IF ( dp_smooth )  THEN
[1340]1929          dp_smooth_factor(:dp_level_ind_b) = 0.0_wp
[240]1930          DO  k = dp_level_ind_b+1, nzt
[4648]1931             dp_smooth_factor(k) = 0.5_wp * ( 1.0_wp + SIN( pi *                                   &
1932                                             ( REAL( k - dp_level_ind_b, KIND=wp ) /               &
1933                                               REAL( nzt - dp_level_ind_b, KIND=wp ) - 0.5_wp ) ) )
[240]1934          ENDDO
1935       ENDIF
1936    ENDIF
1937
1938!
[4648]1939!-- Initialize damping zone for the potential temperature in case of non-cyclic lateral boundaries.
1940!-- The damping zone has the maximum value at the inflow boundary and decreases to zero at
1941!-- pt_damping_width.
[1340]1942    ptdf_x = 0.0_wp
1943    ptdf_y = 0.0_wp
[1159]1944    IF ( bc_lr_dirrad )  THEN
[996]1945       DO  i = nxl, nxr
[978]1946          IF ( ( i * dx ) < pt_damping_width )  THEN
[4648]1947             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5_wp *                                  &
1948                                                    REAL( pt_damping_width - i * dx, KIND=wp ) /   &
1949                                                    REAL( pt_damping_width, KIND=wp ) ) )**2
1950                                  ENDIF
[73]1951       ENDDO
[1159]1952    ELSEIF ( bc_lr_raddir )  THEN
[996]1953       DO  i = nxl, nxr
[978]1954          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
[4648]1955             ptdf_x(i) = pt_damping_factor * SIN( pi * 0.5_wp *                                    &
1956                                                  ( ( i - nx ) * dx + pt_damping_width ) /         &
1957                                                  REAL( pt_damping_width, KIND=wp ) )**2
[73]1958          ENDIF
[4648]1959       ENDDO
[1159]1960    ELSEIF ( bc_ns_dirrad )  THEN
[996]1961       DO  j = nys, nyn
[978]1962          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
[4648]1963             ptdf_y(j) = pt_damping_factor * SIN( pi * 0.5_wp *                                    &
1964                                                  ( ( j - ny ) * dy + pt_damping_width ) /         &
1965                                                  REAL( pt_damping_width, KIND=wp ) )**2
[1]1966          ENDIF
[4648]1967       ENDDO
[1159]1968    ELSEIF ( bc_ns_raddir )  THEN
[996]1969       DO  j = nys, nyn
[978]1970          IF ( ( j * dy ) < pt_damping_width )  THEN
[4648]1971             ptdf_y(j) = pt_damping_factor * SIN( pi * 0.5_wp *                                    &
1972                                                  ( pt_damping_width - j * dy ) /                  &
1973                                                  REAL( pt_damping_width, KIND=wp ) )**2
[1]1974          ENDIF
[73]1975       ENDDO
[1]1976    ENDIF
[51]1977
[1]1978!
[4648]1979!-- Input binary data file is not needed anymore. This line must be placed after call of user_init!
[1]1980    CALL close_file( 13 )
[2934]1981!
[4648]1982!-- In case of nesting, put an barrier to assure that all parent and child domains finished
1983!-- initialization.
[2934]1984#if defined( __parallel )
1985    IF ( nested_run )  CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
1986#endif
[1]1987
[2934]1988
[3987]1989    CALL location_message( 'model initialization', 'finished' )
[1]1990
1991 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.