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

Last change on this file since 4180 was 4180, checked in by scharf, 5 years ago

removed comments in 'Former revisions' section that are older than 01.01.2019

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