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

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