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

Last change on this file since 4533 was 4514, checked in by suehring, 4 years ago

Bugfix in plant-canopy model for output of averaged transpiration rate after a restart; Revise check for output for plant heating rate and rename error message number; Surface-data output: enable output of mixing ratio and passive scalar concentration at the surface; Surface-data input: Add possibility to prescribe surface sensible and latent heat fluxes via static input file

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