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

Last change on this file since 2965 was 2938, checked in by suehring, 7 years ago

Nesting in RANS-LES and RANS-RANS mode enabled; synthetic turbulence generator at all lateral boundaries in nesting or non-cyclic forcing mode; revised Inifor initialization in nesting mode

  • Property svn:keywords set to Id
  • Property svn:mergeinfo set to False
    /palm/branches/forwind/SOURCE/init_3d_model.f901564-1913
    /palm/branches/palm4u/SOURCE/init_3d_model.f902540-2692
File size: 96.5 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-2018 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: init_3d_model.f90 2938 2018-03-27 15:52:42Z scharf $
27! - Revise Inifor initialization for geostrophic wind components
28! - Initialize synthetic turbulence generator in case of Inifor initialization 
29!
30! 2936 2018-03-27 14:49:27Z suehring
31! Synchronize parent and child models after initialization.
32! Remove obsolete masking of topography grid points for Runge-Kutta weighted
33! tendency arrays.
34!
35! 2920 2018-03-22 11:22:01Z kanani
36! Add call for precalculating apparent solar positions (moh.hefny)
37!
38! 2906 2018-03-19 08:56:40Z Giersch
39! The variables read/write_svf_on_init have been removed. Instead ENVIRONMENT
40! variables read/write_svf have been introduced. Location_message has been
41! added.
42!
43! 2894 2018-03-15 09:17:58Z Giersch
44! Renamed routines with respect to reading restart data, file 13 is closed in
45! rrd_read_parts_of_global now
46!
47! 2867 2018-03-09 09:40:23Z suehring
48! Further bugfix concerning call of user_init.
49!
50! 2864 2018-03-08 11:57:45Z suehring
51! Bugfix, move call of user_init in front of initialization of grid-point
52! arrays
53!
54! 2817 2018-02-19 16:32:21Z knoop
55! Preliminary gust module interface implemented
56!
57! 2776 2018-01-31 10:44:42Z Giersch
58! Variable use_synthetic_turbulence_generator has been abbreviated
59!
60! 2766 2018-01-22 17:17:47Z kanani
61! Removed preprocessor directive __chem
62!
63! 2758 2018-01-17 12:55:21Z suehring
64! In case of spinup of land- and urban-surface model, do not mask wind velocity
65! at first computational grid level
66!
67! 2746 2018-01-15 12:06:04Z suehring
68! Move flag plant canopy to modules
69!
70! 2718 2018-01-02 08:49:38Z maronga
71! Corrected "Former revisions" section
72!
73! 2705 2017-12-18 11:26:23Z maronga
74! Bugfix for reading initial profiles from ls/nuding file
75!
76! 2701 2017-12-15 15:40:50Z suehring
77! Changes from last commit documented
78!
79! 2700 2017-12-15 14:12:35Z suehring
80! Bugfix, missing initialization of surface attributes in case of
81! inifor-initialization branch
82!
83! 2698 2017-12-14 18:46:24Z suehring
84! Bugfix in get_topography_top_index
85!
86! 2696 2017-12-14 17:12:51Z kanani
87! Change in file header (GPL part)
88! Implementation of uv exposure model (FK)
89! Moved initialisation of diss, e, kh, km to turbulence_closure_mod (TG)
90! Added chemical emissions (FK)
91! Initialize masking arrays and number-of-grid-points arrays before initialize
92! LSM, USM and radiation module
93! Initialization with inifor (MS)
94!
95! 2618 2017-11-16 15:37:30Z suehring
96! Reorder calls of init_surfaces.
97!
98! 2564 2017-10-19 15:56:56Z Giersch
99! Variable wind_turbine was added to control_parameters.
100!
101! 2550 2017-10-16 17:12:01Z boeske
102! Modifications to cyclic fill method and turbulence recycling method in case of
103! complex terrain simulations
104!
105! 2513 2017-10-04 09:24:39Z kanani
106! Bugfix in storing initial scalar profile (wrong index)
107!
108! 2350 2017-08-15 11:48:26Z kanani
109! Bugfix in nopointer version
110!
111! 2339 2017-08-07 13:55:26Z gronemeier
112! corrected timestamp in header
113!
114! 2338 2017-08-07 12:15:38Z gronemeier
115! Modularize 1D model
116!
117! 2329 2017-08-03 14:24:56Z knoop
118! Removed temporary bugfix (r2327) as bug is properly resolved by this revision
119!
120! 2327 2017-08-02 07:40:57Z maronga
121! Temporary bugfix
122!
123! 2320 2017-07-21 12:47:43Z suehring
124! Modularize large-scale forcing and nudging
125!
126! 2292 2017-06-20 09:51:42Z schwenkel
127! Implementation of new microphysic scheme: cloud_scheme = 'morrison'
128! includes two more prognostic equations for cloud drop concentration (nc) 
129! and cloud water content (qc).
130!
131! 2277 2017-06-12 10:47:51Z kanani
132! Removed unused variable sums_up_fraction_l
133!
134! 2270 2017-06-09 12:18:47Z maronga
135! dots_num must be increased when LSM and/or radiation is used
136!
137! 2259 2017-06-08 09:09:11Z gronemeier
138! Implemented synthetic turbulence generator
139!
140! 2252 2017-06-07 09:35:37Z knoop
141! rho_air now depending on surface_pressure even in Boussinesq mode
142!
143! 2233 2017-05-30 18:08:54Z suehring
144!
145! 2232 2017-05-30 17:47:52Z suehring
146! Adjustments to new topography and surface concept:
147!   - Modify passed parameters for disturb_field
148!   - Topography representation via flags
149!   - Remove unused arrays.
150!   - Move initialization of surface-related quantities to surface_mod
151!
152! 2172 2017-03-08 15:55:25Z knoop
153! Bugfix: moved parallel random generator initialization into its module
154!
155! 2118 2017-01-17 16:38:49Z raasch
156! OpenACC directives removed
157!
158! 2037 2016-10-26 11:15:40Z knoop
159! Anelastic approximation implemented
160!
161! 2031 2016-10-21 15:11:58Z knoop
162! renamed variable rho to rho_ocean
163!
164! 2011 2016-09-19 17:29:57Z kanani
165! Flag urban_surface is now defined in module control_parameters.
166!
167! 2007 2016-08-24 15:47:17Z kanani
168! Added support for urban surface model,
169! adjusted location_message in case of plant_canopy
170!
171! 2000 2016-08-20 18:09:15Z knoop
172! Forced header and separation lines into 80 columns
173!
174! 1992 2016-08-12 15:14:59Z suehring
175! Initializaton of scalarflux at model top
176! Bugfixes in initialization of surface and top salinity flux, top scalar and
177! humidity fluxes
178!
179! 1960 2016-07-12 16:34:24Z suehring
180! Separate humidity and passive scalar
181! Increase dimension for mean_inflow_profiles
182! Remove inadvertent write-statement
183! Bugfix, large-scale forcing is still not implemented for passive scalars
184!
185! 1957 2016-07-07 10:43:48Z suehring
186! flight module added
187!
188! 1920 2016-05-30 10:50:15Z suehring
189! Initialize us with very small number to avoid segmentation fault during
190! calculation of Obukhov length
191!
192! 1918 2016-05-27 14:35:57Z raasch
193! intermediate_timestep_count is set 0 instead 1 for first call of pres,
194! bugfix: initialization of local sum arrays are moved to the beginning of the
195!         routine because otherwise results from pres are overwritten
196!
197! 1914 2016-05-26 14:44:07Z witha
198! Added initialization of the wind turbine model
199!
200! 1878 2016-04-19 12:30:36Z hellstea
201! The zeroth element of weight_pres removed as unnecessary
202!
203! 1849 2016-04-08 11:33:18Z hoffmann
204! Adapted for modularization of microphysics.
205! precipitation_amount, precipitation_rate, prr moved to arrays_3d.
206! Initialization of nc_1d, nr_1d, pt_1d, qc_1d, qr_1d, q_1d moved to
207! microphysics_init.
208!
209! 1845 2016-04-08 08:29:13Z raasch
210! nzb_2d replaced by nzb_u|v_inner
211!
212! 1833 2016-04-07 14:23:03Z raasch
213! initialization of spectra quantities moved to spectra_mod
214!
215! 1831 2016-04-07 13:15:51Z hoffmann
216! turbulence renamed collision_turbulence
217!
218! 1826 2016-04-07 12:01:39Z maronga
219! Renamed radiation calls.
220! Renamed canopy model calls.
221!
222! 1822 2016-04-07 07:49:42Z hoffmann
223! icloud_scheme replaced by microphysics_*
224!
225! 1817 2016-04-06 15:44:20Z maronga
226! Renamed lsm calls.
227!
228! 1815 2016-04-06 13:49:59Z raasch
229! zero-settings for velocities inside topography re-activated (was deactivated
230! in r1762)
231!
232! 1788 2016-03-10 11:01:04Z maronga
233! Added z0q.
234! Syntax layout improved.
235!
236! 1783 2016-03-06 18:36:17Z raasch
237! netcdf module name changed + related changes
238!
239! 1764 2016-02-28 12:45:19Z raasch
240! bugfix: increase size of volume_flow_area_l and volume_flow_initial_l by 1
241!
242! 1762 2016-02-25 12:31:13Z hellstea
243! Introduction of nested domain feature
244!
245! 1738 2015-12-18 13:56:05Z raasch
246! calculate mean surface level height for each statistic region
247!
248! 1734 2015-12-02 12:17:12Z raasch
249! no initial disturbances in case that the disturbance energy limit has been
250! set zero
251!
252! 1707 2015-11-02 15:24:52Z maronga
253! Bugfix: transfer of Richardson number from 1D model to Obukhov length caused
254! devision by zero in neutral stratification
255!
256! 1691 2015-10-26 16:17:44Z maronga
257! Call to init_surface_layer added. rif is replaced by ol and zeta.
258!
259! 1682 2015-10-07 23:56:08Z knoop
260! Code annotations made doxygen readable
261!
262! 1615 2015-07-08 18:49:19Z suehring
263! Enable turbulent inflow for passive_scalar and humidity
264!
265! 1585 2015-04-30 07:05:52Z maronga
266! Initialization of radiation code is now done after LSM initializtion
267!
268! 1575 2015-03-27 09:56:27Z raasch
269! adjustments for psolver-queries
270!
271! 1551 2015-03-03 14:18:16Z maronga
272! Allocation of land surface arrays is now done in the subroutine lsm_init_arrays,
273! which is part of land_surface_model.
274!
275! 1507 2014-12-10 12:14:18Z suehring
276! Bugfix: set horizontal velocity components to zero inside topography
277!
278! 1496 2014-12-02 17:25:50Z maronga
279! Added initialization of the land surface and radiation schemes
280!
281! 1484 2014-10-21 10:53:05Z kanani
282! Changes due to new module structure of the plant canopy model:
283! canopy-related initialization (e.g. lad and canopy_heat_flux) moved to new
284! subroutine init_plant_canopy within the module plant_canopy_model_mod,
285! call of subroutine init_plant_canopy added.
286!
287! 1431 2014-07-15 14:47:17Z suehring
288! var_d added, in order to normalize spectra.
289!
290! 1429 2014-07-15 12:53:45Z knoop
291! Ensemble run capability added to parallel random number generator
292!
293! 1411 2014-05-16 18:01:51Z suehring
294! Initial horizontal velocity profiles were not set to zero at the first vertical
295! grid level in case of non-cyclic lateral boundary conditions.
296!
297! 1406 2014-05-16 13:47:01Z raasch
298! bugfix: setting of initial velocities at k=1 to zero not in case of a
299! no-slip boundary condition for uv
300!
301! 1402 2014-05-09 14:25:13Z raasch
302! location messages modified
303!
304! 1400 2014-05-09 14:03:54Z knoop
305! Parallel random number generator added
306!
307! 1384 2014-05-02 14:31:06Z raasch
308! location messages added
309!
310! 1361 2014-04-16 15:17:48Z hoffmann
311! tend_* removed
312! Bugfix: w_subs is not allocated anymore if it is already allocated
313!
314! 1359 2014-04-11 17:15:14Z hoffmann
315! module lpm_init_mod added to use statements, because lpm_init has become a
316! module
317!
318! 1353 2014-04-08 15:21:23Z heinze
319! REAL constants provided with KIND-attribute
320!
321! 1340 2014-03-25 19:45:13Z kanani
322! REAL constants defined as wp-kind
323!
324! 1322 2014-03-20 16:38:49Z raasch
325! REAL constants defined as wp-kind
326! module interfaces removed
327!
328! 1320 2014-03-20 08:40:49Z raasch
329! ONLY-attribute added to USE-statements,
330! kind-parameters added to all INTEGER and REAL declaration statements,
331! kinds are defined in new module kinds,
332! revision history before 2012 removed,
333! comment fields (!:) to be used for variable explanations added to
334! all variable declaration statements
335!
336! 1316 2014-03-17 07:44:59Z heinze
337! Bugfix: allocation of w_subs
338!
339! 1299 2014-03-06 13:15:21Z heinze
340! Allocate w_subs due to extension of large scale subsidence in combination
341! with large scale forcing data (LSF_DATA)
342!
343! 1241 2013-10-30 11:36:58Z heinze
344! Overwrite initial profiles in case of nudging
345! Inititialize shf and qsws in case of large_scale_forcing
346!
347! 1221 2013-09-10 08:59:13Z raasch
348! +rflags_s_inner in copyin statement, use copyin for most arrays instead of
349! copy
350!
351! 1212 2013-08-15 08:46:27Z raasch
352! array tri is allocated and included in data copy statement
353!
354! 1195 2013-07-01 12:27:57Z heinze
355! Bugfix: move allocation of ref_state to parin.f90 and read_var_list.f90
356!
357! 1179 2013-06-14 05:57:58Z raasch
358! allocate and set ref_state to be used in buoyancy terms
359!
360! 1171 2013-05-30 11:27:45Z raasch
361! diss array is allocated with full size if accelerator boards are used
362!
363! 1159 2013-05-21 11:58:22Z fricke
364! -bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir
365!
366! 1153 2013-05-10 14:33:08Z raasch
367! diss array is allocated with dummy elements even if it is not needed
368! (required by PGI 13.4 / CUDA 5.0)
369!
370! 1115 2013-03-26 18:16:16Z hoffmann
371! unused variables removed
372!
373! 1113 2013-03-10 02:48:14Z raasch
374! openACC directive modified
375!
376! 1111 2013-03-08 23:54:10Z raasch
377! openACC directives added for pres
378! array diss allocated only if required
379!
380! 1092 2013-02-02 11:24:22Z raasch
381! unused variables removed
382!
383! 1065 2012-11-22 17:42:36Z hoffmann
384! allocation of diss (dissipation rate) in case of turbulence = .TRUE. added
385!
386! 1053 2012-11-13 17:11:03Z hoffmann
387! allocation and initialisation of necessary data arrays for the two-moment
388! cloud physics scheme the two new prognostic equations (nr, qr):
389! +dr, lambda_r, mu_r, sed_*, xr, *s, *sws, *swst, *, *_p, t*_m, *_1, *_2, *_3,
390! +tend_*, prr
391!
392! 1036 2012-10-22 13:43:42Z raasch
393! code put under GPL (PALM 3.9)
394!
395! 1032 2012-10-21 13:03:21Z letzel
396! save memory by not allocating pt_2 in case of neutral = .T.
397!
398! 1025 2012-10-07 16:04:41Z letzel
399! bugfix: swap indices of mask for ghost boundaries
400!
401! 1015 2012-09-27 09:23:24Z raasch
402! mask is set to zero for ghost boundaries
403!
404! 1010 2012-09-20 07:59:54Z raasch
405! cpp switch __nopointer added for pointer free version
406!
407! 1003 2012-09-14 14:35:53Z raasch
408! nxra,nyna, nzta replaced ny nxr, nyn, nzt
409!
410! 1001 2012-09-13 14:08:46Z raasch
411! all actions concerning leapfrog scheme removed
412!
413! 996 2012-09-07 10:41:47Z raasch
414! little reformatting
415!
416! 978 2012-08-09 08:28:32Z fricke
417! outflow damping layer removed
418! roughness length for scalar quantites z0h added
419! damping zone for the potential temperatur in case of non-cyclic lateral
420! boundaries added
421! initialization of ptdf_x, ptdf_y
422! initialization of c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l
423!
424! 849 2012-03-15 10:35:09Z raasch
425! init_particles renamed lpm_init
426!
427! 825 2012-02-19 03:03:44Z raasch
428! wang_collision_kernel renamed wang_kernel
429!
430! Revision 1.1  1998/03/09 16:22:22  raasch
431! Initial revision
432!
433!
434! Description:
435! ------------
436!> Allocation of arrays and initialization of the 3D model via
437!> a) pre-run the 1D model
438!> or
439!> b) pre-set constant linear profiles
440!> or
441!> c) read values of a previous run
442!------------------------------------------------------------------------------!
443 SUBROUTINE init_3d_model
444 
445
446    USE advec_ws
447
448    USE arrays_3d
449
450    USE chemistry_model_mod,                                                   &
451        ONLY:  chem_emissions
452
453    USE cloud_parameters,                                                      &
454        ONLY:  cp, l_v, r_d
455
456    USE constants,                                                             &
457        ONLY:  pi
458   
459    USE control_parameters
460   
461    USE flight_mod,                                                            &
462        ONLY:  flight_init
463   
464    USE grid_variables,                                                        &
465        ONLY:  dx, dy, ddx2_mg, ddy2_mg
466
467    USE gust_mod,                                                              &
468        ONLY:  gust_init, gust_init_arrays, gust_module_enabled
469   
470    USE indices
471
472    USE lpm_init_mod,                                                          &
473        ONLY:  lpm_init
474   
475    USE kinds
476
477    USE land_surface_model_mod,                                                &
478        ONLY:  lsm_init, lsm_init_arrays
479 
480    USE lsf_nudging_mod,                                                       &
481        ONLY:  lsf_init, ls_forcing_surf, nudge_init
482
483    USE microphysics_mod,                                                      &
484        ONLY:  collision_turbulence, microphysics_init
485
486    USE model_1d_mod,                                                          &
487        ONLY:  e1d, init_1d_model, kh1d, km1d, l1d, rif1d, u1d, us1d, usws1d,  &
488               v1d, vsws1d 
489
490    USE netcdf_interface,                                                      &
491        ONLY:  dots_max, dots_num, dots_unit, dots_label
492
493    USE netcdf_data_input_mod,                                                 &
494        ONLY:  init_3d, netcdf_data_input_interpolate, netcdf_data_input_init_3d
495   
496    USE particle_attributes,                                                   &
497        ONLY:  particle_advection, use_sgs_for_particles, wang_kernel
498   
499    USE pegrid
500   
501    USE plant_canopy_model_mod,                                                &
502        ONLY:  pcm_init
503
504    USE pmc_interface,                                                         &
505        ONLY:  nested_run
506
507    USE radiation_model_mod,                                                   &
508        ONLY:  radiation_init, radiation, radiation_control, radiation_scheme, &
509               radiation_calc_svf, radiation_write_svf,                        &
510               radiation_interaction, radiation_interactions,                  &
511               radiation_interaction_init, radiation_read_svf,                 &
512               radiation_presimulate_solar_pos
513   
514    USE random_function_mod 
515   
516    USE random_generator_parallel,                                             &
517        ONLY:  init_parallel_random_generator
518
519    USE read_restart_data_mod,                                                 &
520        ONLY:  rrd_read_parts_of_global, rrd_local                                     
521   
522    USE statistics,                                                            &
523        ONLY:  hom, hom_sum, mean_surface_level_height, pr_palm, rmask,        &
524               statistic_regions, sums, sums_divnew_l, sums_divold_l, sums_l,  &
525               sums_l_l, sums_wsts_bc_l, ts_value,                             &
526               weight_pres, weight_substep
527
528    USE synthetic_turbulence_generator_mod,                                    &
529        ONLY:  stg_init, use_syn_turb_gen
530
531    USE surface_layer_fluxes_mod,                                              &
532        ONLY:  init_surface_layer_fluxes
533
534    USE surface_mod,                                                           &
535        ONLY :  init_surface_arrays, init_surfaces, surf_def_h, surf_lsm_h,    &
536                surf_usm_h, get_topography_top_index_ji
537   
538    USE transpose_indices
539
540    USE turbulence_closure_mod,                                                &
541        ONLY:  tcm_init_arrays, tcm_init
542
543    USE urban_surface_mod,                                                     &
544        ONLY:  usm_init_urban_surface, usm_allocate_surface
545
546    USE uv_exposure_model_mod,                                                 &
547        ONLY:  uvem_init, uvem_init_arrays
548
549    USE wind_turbine_model_mod,                                                &
550        ONLY:  wtm_init, wtm_init_arrays
551
552    IMPLICIT NONE
553
554    INTEGER(iwp) ::  i             !<
555    INTEGER(iwp) ::  ind_array(1)  !<
556    INTEGER(iwp) ::  j             !<
557    INTEGER(iwp) ::  k             !<
558    INTEGER(iwp) ::  k_surf        !< surface level index
559    INTEGER(iwp) ::  m             !< index of surface element in surface data type
560    INTEGER(iwp) ::  sr            !< index of statistic region
561
562    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  ngp_2dh_l  !<
563
564    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l    !<
565    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner_l  !<
566
567    REAL(wp)     ::  t_surface !< air temperature at the surface
568
569    REAL(wp), DIMENSION(:), ALLOCATABLE ::  p_hydrostatic !< hydrostatic pressure
570
571    INTEGER(iwp) ::  l       !< loop variable
572    INTEGER(iwp) ::  nzt_l   !< index of top PE boundary for multigrid level
573    REAL(wp) ::  dx_l !< grid spacing along x on different multigrid level
574    REAL(wp) ::  dy_l !< grid spacing along y on different multigrid level
575
576    REAL(wp), DIMENSION(1:3) ::  volume_flow_area_l     !<
577    REAL(wp), DIMENSION(1:3) ::  volume_flow_initial_l  !<
578
579    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height_l    !<
580    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l    !<
581    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_tmp  !<
582
583    INTEGER(iwp) ::  nz_u_shift   !<
584    INTEGER(iwp) ::  nz_v_shift   !<
585    INTEGER(iwp) ::  nz_w_shift   !<
586    INTEGER(iwp) ::  nz_s_shift   !<
587    INTEGER(iwp) ::  nz_u_shift_l !<
588    INTEGER(iwp) ::  nz_v_shift_l !<
589    INTEGER(iwp) ::  nz_w_shift_l !<
590    INTEGER(iwp) ::  nz_s_shift_l !<
591
592    CALL location_message( 'allocating arrays', .FALSE. )
593!
594!-- Allocate arrays
595    ALLOCATE( mean_surface_level_height(0:statistic_regions),                  &
596              mean_surface_level_height_l(0:statistic_regions),                &
597              ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions),    &
598              ngp_3d(0:statistic_regions),                                     &
599              ngp_3d_inner(0:statistic_regions),                               &
600              ngp_3d_inner_l(0:statistic_regions),                             &
601              ngp_3d_inner_tmp(0:statistic_regions),                           &
602              sums_divnew_l(0:statistic_regions),                              &
603              sums_divold_l(0:statistic_regions) )
604    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
605    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                    &
606              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),                  &
607              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),                  &
608              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),                &
609              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),                  &
610              sums(nzb:nzt+1,pr_palm+max_pr_user),                             &
611              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),      &
612              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1),    &
613              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                   &
614              ts_value(dots_max,0:statistic_regions) )
615    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
616
617    ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr),                                    &
618              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
619              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
620
621#if defined( __nopointer )
622    ALLOCATE( pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
623              pt_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
624              u(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
625              u_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
626              v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
627              v_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
628              w(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
629              w_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
630              tpt_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
631              tu_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
632              tv_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
633              tw_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
634#else
635    ALLOCATE( pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
636              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
637              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
638              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
639              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
640              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
641              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
642              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
643              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
644              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
645              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
646    IF (  .NOT.  neutral )  THEN
647       ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
648    ENDIF
649#endif
650
651!
652!-- Following array is required for perturbation pressure within the iterative
653!-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds
654!-- the weighted average of the substeps and cannot be used in the Poisson
655!-- solver.
656    IF ( psolver == 'sor' )  THEN
657       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
658    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
659!
660!--    For performance reasons, multigrid is using one ghost layer only
661       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
662    ENDIF
663
664!
665!-- Array for storing constant coeffficients of the tridiagonal solver
666    IF ( psolver == 'poisfft' )  THEN
667       ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) )
668       ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) )
669    ENDIF
670
671    IF ( humidity )  THEN
672!
673!--    3D-humidity
674#if defined( __nopointer )
675       ALLOCATE( q(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
676                 q_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
677                 tq_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
678#else
679       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
680                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
681                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
682#endif
683
684!
685!--    3D-arrays needed for humidity
686       IF ( humidity )  THEN
687#if defined( __nopointer )
688          ALLOCATE( vpt(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
689#else
690          ALLOCATE( vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
691#endif
692
693          IF ( cloud_physics )  THEN
694!
695!--          Liquid water content
696#if defined( __nopointer )
697             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
698#else
699             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
700#endif
701
702!
703!--          3D-cloud water content
704             IF ( .NOT. microphysics_morrison )  THEN
705#if defined( __nopointer )
706                ALLOCATE( qc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
707#else
708                ALLOCATE( qc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
709#endif
710             ENDIF
711!
712!--          Precipitation amount and rate (only needed if output is switched)
713             ALLOCATE( precipitation_amount(nysg:nyng,nxlg:nxrg),              &
714                       precipitation_rate(nysg:nyng,nxlg:nxrg) )
715
716!
717!--          3d-precipitation rate
718             ALLOCATE( prr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
719
720             IF ( microphysics_morrison )  THEN
721!
722!--             3D-cloud drop water content, cloud drop concentration arrays
723#if defined( __nopointer )
724                ALLOCATE( nc(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
725                          nc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
726                          qc(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
727                          qc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
728                          tnc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),             & 
729                          tqc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
730#else
731                ALLOCATE( nc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
732                          nc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
733                          nc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
734                          qc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
735                          qc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
736                          qc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
737#endif
738             ENDIF
739
740             IF ( microphysics_seifert )  THEN
741!
742!--             3D-rain water content, rain drop concentration arrays
743#if defined( __nopointer )
744                ALLOCATE( nr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
745                          nr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
746                          qr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
747                          qr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
748                          tnr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),             &
749                          tqr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
750#else
751                ALLOCATE( nr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
752                          nr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
753                          nr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
754                          qr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
755                          qr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
756                          qr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
757#endif
758             ENDIF
759
760          ENDIF
761
762          IF ( cloud_droplets )  THEN
763!
764!--          Liquid water content, change in liquid water content
765#if defined( __nopointer )
766             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
767                        ql_c(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
768#else
769             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
770                        ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
771#endif
772!
773!--          Real volume of particles (with weighting), volume of particles
774             ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
775                        ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
776          ENDIF
777
778       ENDIF
779
780    ENDIF
781   
782   
783    IF ( passive_scalar )  THEN
784
785!
786!--    3D scalar arrays
787#if defined( __nopointer )
788       ALLOCATE( s(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
789                 s_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
790                 ts_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
791#else
792       ALLOCATE( s_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
793                 s_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
794                 s_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
795#endif
796    ENDIF
797
798    IF ( ocean )  THEN
799#if defined( __nopointer )
800       ALLOCATE( prho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
801                 rho_ocean(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
802                 sa(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
803                 sa_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
804                 tsa_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
805#else
806       ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
807                 rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                         &
808                 sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
809                 sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
810                 sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
811       prho => prho_1
812       rho_ocean  => rho_1  ! routines calc_mean_profile and diffusion_e require
813                      ! density to be apointer
814#endif
815    ENDIF
816
817!
818!-- Allocation of anelastic and Boussinesq approximation specific arrays
819    ALLOCATE( p_hydrostatic(nzb:nzt+1) )
820    ALLOCATE( rho_air(nzb:nzt+1) )
821    ALLOCATE( rho_air_zw(nzb:nzt+1) )
822    ALLOCATE( drho_air(nzb:nzt+1) )
823    ALLOCATE( drho_air_zw(nzb:nzt+1) )
824
825!
826!-- Density profile calculation for anelastic approximation
827    t_surface = pt_surface * ( surface_pressure / 1000.0_wp )**( r_d / cp )
828    IF ( TRIM( approximation ) == 'anelastic' ) THEN
829       DO  k = nzb, nzt+1
830          p_hydrostatic(k)    = surface_pressure * 100.0_wp *                  &
831                                ( 1 - ( g * zu(k) ) / ( cp * t_surface )       &
832                                )**( cp / r_d )
833          rho_air(k)          = ( p_hydrostatic(k) *                           &
834                                  ( 100000.0_wp / p_hydrostatic(k)             &
835                                  )**( r_d / cp )                              &
836                                ) / ( r_d * pt_init(k) )
837       ENDDO
838       DO  k = nzb, nzt
839          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
840       ENDDO
841       rho_air_zw(nzt+1)  = rho_air_zw(nzt)                                    &
842                            + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
843    ELSE
844       DO  k = nzb, nzt+1
845          p_hydrostatic(k)    = surface_pressure * 100.0_wp *                  &
846                                ( 1 - ( g * zu(nzb) ) / ( cp * t_surface )       &
847                                )**( cp / r_d )
848          rho_air(k)          = ( p_hydrostatic(k) *                           &
849                                  ( 100000.0_wp / p_hydrostatic(k)             &
850                                  )**( r_d / cp )                              &
851                                ) / ( r_d * pt_init(nzb) )
852       ENDDO
853       DO  k = nzb, nzt
854          rho_air_zw(k) = 0.5_wp * ( rho_air(k) + rho_air(k+1) )
855       ENDDO
856       rho_air_zw(nzt+1)  = rho_air_zw(nzt)                                    &
857                            + 2.0_wp * ( rho_air(nzt+1) - rho_air_zw(nzt)  )
858    ENDIF
859!
860!-- compute the inverse density array in order to avoid expencive divisions
861    drho_air    = 1.0_wp / rho_air
862    drho_air_zw = 1.0_wp / rho_air_zw
863
864!
865!-- Allocation of flux conversion arrays
866    ALLOCATE( heatflux_input_conversion(nzb:nzt+1) )
867    ALLOCATE( waterflux_input_conversion(nzb:nzt+1) )
868    ALLOCATE( momentumflux_input_conversion(nzb:nzt+1) )
869    ALLOCATE( heatflux_output_conversion(nzb:nzt+1) )
870    ALLOCATE( waterflux_output_conversion(nzb:nzt+1) )
871    ALLOCATE( momentumflux_output_conversion(nzb:nzt+1) )
872
873!
874!-- calculate flux conversion factors according to approximation and in-/output mode
875    DO  k = nzb, nzt+1
876
877        IF ( TRIM( flux_input_mode ) == 'kinematic' )  THEN
878            heatflux_input_conversion(k)      = rho_air_zw(k)
879            waterflux_input_conversion(k)     = rho_air_zw(k)
880            momentumflux_input_conversion(k)  = rho_air_zw(k)
881        ELSEIF ( TRIM( flux_input_mode ) == 'dynamic' ) THEN
882            heatflux_input_conversion(k)      = 1.0_wp / cp
883            waterflux_input_conversion(k)     = 1.0_wp / l_v
884            momentumflux_input_conversion(k)  = 1.0_wp
885        ENDIF
886
887        IF ( TRIM( flux_output_mode ) == 'kinematic' )  THEN
888            heatflux_output_conversion(k)     = drho_air_zw(k)
889            waterflux_output_conversion(k)    = drho_air_zw(k)
890            momentumflux_output_conversion(k) = drho_air_zw(k)
891        ELSEIF ( TRIM( flux_output_mode ) == 'dynamic' ) THEN
892            heatflux_output_conversion(k)     = cp
893            waterflux_output_conversion(k)    = l_v
894            momentumflux_output_conversion(k) = 1.0_wp
895        ENDIF
896
897        IF ( .NOT. humidity ) THEN
898            waterflux_input_conversion(k)  = 1.0_wp
899            waterflux_output_conversion(k) = 1.0_wp
900        ENDIF
901
902    ENDDO
903
904!
905!-- In case of multigrid method, compute grid lengths and grid factors for the
906!-- grid levels with respective density on each grid
907    IF ( psolver(1:9) == 'multigrid' )  THEN
908
909       ALLOCATE( ddx2_mg(maximum_grid_level) )
910       ALLOCATE( ddy2_mg(maximum_grid_level) )
911       ALLOCATE( dzu_mg(nzb+1:nzt+1,maximum_grid_level) )
912       ALLOCATE( dzw_mg(nzb+1:nzt+1,maximum_grid_level) )
913       ALLOCATE( f1_mg(nzb+1:nzt,maximum_grid_level) )
914       ALLOCATE( f2_mg(nzb+1:nzt,maximum_grid_level) )
915       ALLOCATE( f3_mg(nzb+1:nzt,maximum_grid_level) )
916       ALLOCATE( rho_air_mg(nzb:nzt+1,maximum_grid_level) )
917       ALLOCATE( rho_air_zw_mg(nzb:nzt+1,maximum_grid_level) )
918
919       dzu_mg(:,maximum_grid_level) = dzu
920       rho_air_mg(:,maximum_grid_level) = rho_air
921!       
922!--    Next line to ensure an equally spaced grid.
923       dzu_mg(1,maximum_grid_level) = dzu(2)
924       rho_air_mg(nzb,maximum_grid_level) = rho_air(nzb) +                     &
925                                             (rho_air(nzb) - rho_air(nzb+1))
926
927       dzw_mg(:,maximum_grid_level) = dzw
928       rho_air_zw_mg(:,maximum_grid_level) = rho_air_zw
929       nzt_l = nzt
930       DO  l = maximum_grid_level-1, 1, -1
931           dzu_mg(nzb+1,l) = 2.0_wp * dzu_mg(nzb+1,l+1)
932           dzw_mg(nzb+1,l) = 2.0_wp * dzw_mg(nzb+1,l+1)
933           rho_air_mg(nzb,l)    = rho_air_mg(nzb,l+1) + (rho_air_mg(nzb,l+1) - rho_air_mg(nzb+1,l+1))
934           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))
935           rho_air_mg(nzb+1,l)    = rho_air_mg(nzb+1,l+1)
936           rho_air_zw_mg(nzb+1,l) = rho_air_zw_mg(nzb+1,l+1)
937           nzt_l = nzt_l / 2
938           DO  k = 2, nzt_l+1
939              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
940              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
941              rho_air_mg(k,l)    = rho_air_mg(2*k-1,l+1)
942              rho_air_zw_mg(k,l) = rho_air_zw_mg(2*k-1,l+1)
943           ENDDO
944       ENDDO
945
946       nzt_l = nzt
947       dx_l  = dx
948       dy_l  = dy
949       DO  l = maximum_grid_level, 1, -1
950          ddx2_mg(l) = 1.0_wp / dx_l**2
951          ddy2_mg(l) = 1.0_wp / dy_l**2
952          DO  k = nzb+1, nzt_l
953             f2_mg(k,l) = rho_air_zw_mg(k,l) / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
954             f3_mg(k,l) = rho_air_zw_mg(k-1,l) / ( dzu_mg(k,l)   * dzw_mg(k,l) )
955             f1_mg(k,l) = 2.0_wp * ( ddx2_mg(l) + ddy2_mg(l) ) &
956                          * rho_air_mg(k,l) + f2_mg(k,l) + f3_mg(k,l)
957          ENDDO
958          nzt_l = nzt_l / 2
959          dx_l  = dx_l * 2.0_wp
960          dy_l  = dy_l * 2.0_wp
961       ENDDO
962
963    ENDIF
964
965!
966!-- 1D-array for large scale subsidence velocity
967    IF ( .NOT. ALLOCATED( w_subs ) )  THEN
968       ALLOCATE ( w_subs(nzb:nzt+1) )
969       w_subs = 0.0_wp
970    ENDIF
971
972!
973!-- Arrays to store velocity data from t-dt and the phase speeds which
974!-- are needed for radiation boundary conditions
975    IF ( outflow_l )  THEN
976       ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2),                               &
977                 v_m_l(nzb:nzt+1,nysg:nyng,0:1),                               &
978                 w_m_l(nzb:nzt+1,nysg:nyng,0:1) )
979    ENDIF
980    IF ( outflow_r )  THEN
981       ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
982                 v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
983                 w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) )
984    ENDIF
985    IF ( outflow_l  .OR.  outflow_r )  THEN
986       ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng),           &
987                 c_w(nzb:nzt+1,nysg:nyng) )
988    ENDIF
989    IF ( outflow_s )  THEN
990       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg),                               &
991                 v_m_s(nzb:nzt+1,1:2,nxlg:nxrg),                               &
992                 w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) )
993    ENDIF
994    IF ( outflow_n )  THEN
995       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
996                 v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
997                 w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) )
998    ENDIF
999    IF ( outflow_s  .OR.  outflow_n )  THEN
1000       ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg),           &
1001                 c_w(nzb:nzt+1,nxlg:nxrg) )
1002    ENDIF
1003    IF ( outflow_l  .OR.  outflow_r  .OR.  outflow_s  .OR.  outflow_n )  THEN
1004       ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) )                   
1005       ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) )
1006    ENDIF
1007
1008
1009#if ! defined( __nopointer )
1010!
1011!-- Initial assignment of the pointers
1012    IF ( .NOT. neutral )  THEN
1013       pt => pt_1;  pt_p => pt_2;  tpt_m => pt_3
1014    ELSE
1015       pt => pt_1;  pt_p => pt_1;  tpt_m => pt_3
1016    ENDIF
1017    u  => u_1;   u_p  => u_2;   tu_m  => u_3
1018    v  => v_1;   v_p  => v_2;   tv_m  => v_3
1019    w  => w_1;   w_p  => w_2;   tw_m  => w_3
1020
1021    IF ( humidity )  THEN
1022       q => q_1;  q_p => q_2;  tq_m => q_3
1023       IF ( humidity )  THEN
1024          vpt  => vpt_1   
1025          IF ( cloud_physics )  THEN
1026             ql => ql_1
1027             IF ( .NOT. microphysics_morrison )  THEN
1028                qc => qc_1
1029             ENDIF
1030             IF ( microphysics_morrison )  THEN
1031                qc => qc_1;  qc_p  => qc_2;  tqc_m  => qc_3
1032                nc => nc_1;  nc_p  => nc_2;  tnc_m  => nc_3
1033             ENDIF
1034             IF ( microphysics_seifert )  THEN
1035                qr => qr_1;  qr_p  => qr_2;  tqr_m  => qr_3
1036                nr => nr_1;  nr_p  => nr_2;  tnr_m  => nr_3
1037             ENDIF
1038          ENDIF
1039       ENDIF
1040       IF ( cloud_droplets )  THEN
1041          ql   => ql_1
1042          ql_c => ql_2
1043       ENDIF
1044    ENDIF
1045   
1046    IF ( passive_scalar )  THEN
1047       s => s_1;  s_p => s_2;  ts_m => s_3
1048    ENDIF   
1049
1050    IF ( ocean )  THEN
1051       sa => sa_1;  sa_p => sa_2;  tsa_m => sa_3
1052    ENDIF
1053#endif
1054!
1055!-- Initialize arrays for turbulence closure
1056    CALL tcm_init_arrays
1057!
1058!-- Initialize surface arrays
1059    CALL init_surface_arrays
1060!
1061!-- Allocate land surface model arrays
1062    IF ( land_surface )  THEN
1063       CALL lsm_init_arrays
1064    ENDIF
1065
1066!
1067!-- Allocate wind turbine model arrays
1068    IF ( wind_turbine )  THEN
1069       CALL wtm_init_arrays
1070    ENDIF
1071!
1072!-- Allocate gust module arrays
1073    IF ( gust_module_enabled )  THEN
1074       CALL gust_init_arrays
1075    ENDIF
1076
1077!
1078!-- Initialize virtual flight measurements
1079    IF ( virtual_flight )  THEN
1080       CALL flight_init
1081    ENDIF
1082
1083!
1084!-- Read uv exposure input data
1085    IF ( uv_exposure )  THEN
1086       CALL uvem_init
1087    ENDIF
1088!
1089!-- Allocate uv exposure arrays
1090    IF ( uv_exposure )  THEN
1091       CALL uvem_init_arrays
1092    ENDIF
1093
1094!
1095!-- Initialize nudging if required
1096    IF ( nudging )  THEN
1097       CALL nudge_init
1098    ENDIF
1099
1100!
1101!-- Initialize reading of large scale forcing from external file - if required
1102    IF ( large_scale_forcing  .OR.  forcing )  THEN
1103       CALL lsf_init
1104    ENDIF
1105
1106!
1107!-- Allocate arrays containing the RK coefficient for calculation of
1108!-- perturbation pressure and turbulent fluxes. At this point values are
1109!-- set for pressure calculation during initialization (where no timestep
1110!-- is done). Further below the values needed within the timestep scheme
1111!-- will be set.
1112    ALLOCATE( weight_substep(1:intermediate_timestep_count_max),               &
1113              weight_pres(1:intermediate_timestep_count_max) )
1114    weight_substep = 1.0_wp
1115    weight_pres    = 1.0_wp
1116    intermediate_timestep_count = 0  ! needed when simulated_time = 0.0
1117       
1118    CALL location_message( 'finished', .TRUE. )
1119
1120!
1121!-- Initialize local summation arrays for routine flow_statistics.
1122!-- This is necessary because they may not yet have been initialized when they
1123!-- are called from flow_statistics (or - depending on the chosen model run -
1124!-- are never initialized)
1125    sums_divnew_l      = 0.0_wp
1126    sums_divold_l      = 0.0_wp
1127    sums_l_l           = 0.0_wp
1128    sums_wsts_bc_l     = 0.0_wp
1129
1130
1131
1132!
1133!-- Initialize model variables
1134    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
1135         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
1136!
1137!--    Initialization with provided input data derived from larger-scale model
1138       IF ( INDEX( initializing_actions, 'inifor' ) /= 0 )  THEN
1139          CALL location_message( 'initializing with INIFOR', .FALSE. )
1140!
1141!--       Read initial 1D profiles from NetCDF file if available.
1142!--       At the moment, only u, v, w, pt and q are provided.
1143          CALL netcdf_data_input_init_3d
1144!
1145!--       Please note, at the moment INIFOR assumes only an equidistant vertical
1146!--       grid. In case of vertical grid stretching, input of inital data
1147!--       need to be inter- and/or extrapolated.
1148!--       Therefore, check if zu grid on file is identical to numeric zw grid.
1149!--       Please note 
1150          IF ( ANY( zu(1:nzt+1) /= init_3d%zu_atmos(1:init_3d%nzu) ) )  THEN
1151
1152             CALL netcdf_data_input_interpolate( init_3d%u_init(nzb+1:nzt+1),  &
1153                                                 zu(nzb+1:nzt+1),              &
1154                                                 init_3d%zu_atmos )
1155             CALL netcdf_data_input_interpolate( init_3d%v_init(nzb+1:nzt+1),  &
1156                                                 zu(nzb+1:nzt+1),              &
1157                                                 init_3d%zu_atmos )
1158!              CALL netcdf_data_input_interpolate( init_3d%w_init(nzb+1:nzt),    &
1159!                                                  zw(nzb+1:nzt),                &
1160!                                                  init_3d%zw_atmos )
1161             IF ( .NOT. neutral )                                              &
1162                CALL netcdf_data_input_interpolate(                            &
1163                                             init_3d%pt_init(nzb+1:nzt+1),     &
1164                                             zu(nzb+1:nzt+1),                  &
1165                                             init_3d%zu_atmos )
1166             IF ( humidity )                                                   &
1167                CALL netcdf_data_input_interpolate(                            &
1168                                             init_3d%q_init(nzb+1:nzt+1),      &
1169                                             zu(nzb+1:nzt+1),                  &
1170                                             init_3d%zu_atmos )
1171          ENDIF
1172
1173          u_init = init_3d%u_init
1174          v_init = init_3d%v_init   
1175          IF( .NOT. neutral )  pt_init = init_3d%pt_init
1176          IF( humidity      )  q_init  = init_3d%q_init
1177
1178!
1179!--       Please note, Inifor provides data from nzb+1 to nzt+1.
1180!--       Initialize pt and q with Neumann condition at nzb.
1181          IF( .NOT. neutral )  pt_init(nzb) = pt_init(nzb+1)
1182          IF( humidity      )  q_init(nzb)  = q_init(nzb+1)
1183          DO  i = nxlg, nxrg
1184             DO  j = nysg, nyng
1185                u(:,j,i) = u_init(:)
1186                v(:,j,i) = v_init(:)
1187                IF( .NOT. neutral )  pt(:,j,i) = pt_init(:)
1188                IF( humidity      )  q(:,j,i)  = q_init(:)
1189             ENDDO
1190          ENDDO
1191!
1192!--       MS: What about the geostrophic wind profiles? Actually these
1193!--           are not identical to the initial wind profiles in this case.
1194!--           This need to be further revised.
1195          IF ( init_3d%from_file_ug )  THEN
1196             ug(:) = init_3d%ug_init(:)
1197          ENDIF
1198          IF ( init_3d%from_file_vg )  THEN
1199             vg(:) = init_3d%vg_init(:)
1200          ENDIF
1201
1202          ug(nzt+1) = ug(nzt)
1203          vg(nzt+1) = vg(nzt)
1204
1205!
1206!--       Set inital w to 0
1207          w = 0.0_wp
1208!
1209!--       Initialize the remaining quantities
1210          IF ( humidity )  THEN
1211             IF ( cloud_physics  .AND.  microphysics_morrison )  THEN
1212                DO  i = nxlg, nxrg
1213                   DO  j = nysg, nyng
1214                      qc(:,j,i) = 0.0_wp
1215                      nc(:,j,i) = 0.0_wp
1216                   ENDDO
1217                ENDDO
1218             ENDIF
1219
1220             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1221                DO  i = nxlg, nxrg
1222                   DO  j = nysg, nyng
1223                      qr(:,j,i) = 0.0_wp
1224                      nr(:,j,i) = 0.0_wp
1225                   ENDDO
1226                ENDDO
1227             ENDIF
1228
1229          ENDIF
1230
1231          IF ( passive_scalar )  THEN
1232             DO  i = nxlg, nxrg
1233                DO  j = nysg, nyng
1234                   s(:,j,i) = s_init
1235                ENDDO
1236             ENDDO
1237          ENDIF
1238
1239          IF ( ocean )  THEN
1240             DO  i = nxlg, nxrg
1241                DO  j = nysg, nyng
1242                   sa(:,j,i) = sa_init
1243                ENDDO
1244             ENDDO
1245          ENDIF
1246
1247!
1248!--       Set velocity components at non-atmospheric / oceanic grid points to
1249!--       zero.
1250          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1251          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )
1252          w = MERGE( w, 0.0_wp, BTEST( wall_flags_0, 3 ) )
1253!
1254!--       Initialize surface variables, e.g. friction velocity, momentum
1255!--       fluxes, etc.
1256          CALL init_surfaces
1257!
1258!--       Initialize turbulence generator
1259          IF( use_syn_turb_gen )  CALL stg_init
1260
1261          CALL location_message( 'finished', .TRUE. )
1262!
1263!--    Initialization via computed 1D-model profiles
1264       ELSEIF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
1265
1266          CALL location_message( 'initializing with 1D model profiles', .FALSE. )
1267!
1268!--       Use solutions of the 1D model as initial profiles,
1269!--       start 1D model
1270          CALL init_1d_model
1271!
1272!--       Transfer initial profiles to the arrays of the 3D model
1273          DO  i = nxlg, nxrg
1274             DO  j = nysg, nyng
1275                pt(:,j,i) = pt_init
1276                u(:,j,i)  = u1d
1277                v(:,j,i)  = v1d
1278             ENDDO
1279          ENDDO
1280
1281          IF ( humidity )  THEN
1282             DO  i = nxlg, nxrg
1283                DO  j = nysg, nyng
1284                   q(:,j,i) = q_init
1285                ENDDO
1286             ENDDO
1287             IF ( cloud_physics  .AND.  microphysics_morrison )  THEN
1288                DO  i = nxlg, nxrg
1289                   DO  j = nysg, nyng
1290                      qc(:,j,i) = 0.0_wp
1291                      nc(:,j,i) = 0.0_wp
1292                   ENDDO
1293                ENDDO
1294             ENDIF
1295             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1296                DO  i = nxlg, nxrg
1297                   DO  j = nysg, nyng
1298                      qr(:,j,i) = 0.0_wp
1299                      nr(:,j,i) = 0.0_wp
1300                   ENDDO
1301                ENDDO
1302             ENDIF
1303          ENDIF
1304
1305          IF ( passive_scalar )  THEN
1306             DO  i = nxlg, nxrg
1307                DO  j = nysg, nyng
1308                   s(:,j,i) = s_init
1309                ENDDO
1310             ENDDO   
1311          ENDIF
1312!
1313!--          Store initial profiles for output purposes etc.
1314          IF ( .NOT. constant_diffusion )  THEN
1315             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
1316          ENDIF
1317!
1318!--       Set velocities back to zero
1319          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1320          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )         
1321!
1322!--       WARNING: The extra boundary conditions set after running the
1323!--       -------  1D model impose an error on the divergence one layer
1324!--                below the topography; need to correct later
1325!--       ATTENTION: Provisional correction for Piacsek & Williams
1326!--       ---------  advection scheme: keep u and v zero one layer below
1327!--                  the topography.
1328          IF ( ibc_uv_b == 1 )  THEN
1329!
1330!--          Neumann condition
1331             DO  i = nxl-1, nxr+1
1332                DO  j = nys-1, nyn+1
1333                   u(nzb,j,i) = u(nzb+1,j,i)
1334                   v(nzb,j,i) = v(nzb+1,j,i)
1335                ENDDO
1336             ENDDO
1337
1338          ENDIF
1339!
1340!--       Initialize surface variables, e.g. friction velocity, momentum
1341!--       fluxes, etc.
1342          CALL init_surfaces
1343
1344          CALL location_message( 'finished', .TRUE. )
1345
1346       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )    &
1347       THEN
1348
1349          CALL location_message( 'initializing with constant profiles', .FALSE. )
1350!
1351!--       Overwrite initial profiles in case of synthetic turbulence generator
1352          IF( use_syn_turb_gen )  CALL stg_init
1353
1354!
1355!--       Use constructed initial profiles (velocity constant with height,
1356!--       temperature profile with constant gradient)
1357          DO  i = nxlg, nxrg
1358             DO  j = nysg, nyng
1359                pt(:,j,i) = pt_init
1360                u(:,j,i)  = u_init
1361                v(:,j,i)  = v_init
1362             ENDDO
1363          ENDDO
1364!
1365!--       Mask topography
1366          u = MERGE( u, 0.0_wp, BTEST( wall_flags_0, 1 ) )
1367          v = MERGE( v, 0.0_wp, BTEST( wall_flags_0, 2 ) )
1368!
1369!--       Set initial horizontal velocities at the lowest computational grid
1370!--       levels to zero in order to avoid too small time steps caused by the
1371!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
1372!--       in the limiting formula!).
1373!--       Please note, in case land- or urban-surface model is used and a
1374!--       spinup is applied, masking the lowest computational level is not
1375!--       possible as MOST as well as energy-balance parametrizations will not
1376!--       work with zero wind velocity.
1377          IF ( ibc_uv_b /= 1  .AND.  .NOT.  spinup )  THEN
1378             DO  i = nxlg, nxrg
1379                DO  j = nysg, nyng
1380                   DO  k = nzb, nzt
1381                      u(k,j,i) = MERGE( u(k,j,i), 0.0_wp,                      &
1382                                        BTEST( wall_flags_0(k,j,i), 20 ) )
1383                      v(k,j,i) = MERGE( v(k,j,i), 0.0_wp,                      &
1384                                        BTEST( wall_flags_0(k,j,i), 21 ) )
1385                   ENDDO
1386                ENDDO
1387             ENDDO
1388          ENDIF
1389
1390          IF ( humidity )  THEN
1391             DO  i = nxlg, nxrg
1392                DO  j = nysg, nyng
1393                   q(:,j,i) = q_init
1394                ENDDO
1395             ENDDO
1396             IF ( cloud_physics  .AND.  microphysics_morrison )  THEN
1397                DO  i = nxlg, nxrg
1398                   DO  j = nysg, nyng
1399                      qc(:,j,i) = 0.0_wp
1400                      nc(:,j,i) = 0.0_wp
1401                   ENDDO
1402                ENDDO
1403             ENDIF
1404
1405             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1406                DO  i = nxlg, nxrg
1407                   DO  j = nysg, nyng
1408                      qr(:,j,i) = 0.0_wp
1409                      nr(:,j,i) = 0.0_wp
1410                   ENDDO
1411                ENDDO
1412             ENDIF
1413
1414          ENDIF
1415         
1416          IF ( passive_scalar )  THEN
1417             DO  i = nxlg, nxrg
1418                DO  j = nysg, nyng
1419                   s(:,j,i) = s_init
1420                ENDDO
1421             ENDDO
1422          ENDIF
1423
1424          IF ( ocean )  THEN
1425             DO  i = nxlg, nxrg
1426                DO  j = nysg, nyng
1427                   sa(:,j,i) = sa_init
1428                ENDDO
1429             ENDDO
1430          ENDIF
1431!
1432!--       Compute initial temperature field and other constants used in case
1433!--       of a sloping surface
1434          IF ( sloping_surface )  CALL init_slope
1435!
1436!--       Initialize surface variables, e.g. friction velocity, momentum
1437!--       fluxes, etc.
1438          CALL init_surfaces
1439
1440          CALL location_message( 'finished', .TRUE. )
1441
1442       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )                  &
1443       THEN
1444
1445          CALL location_message( 'initializing by user', .FALSE. )
1446!
1447!--       Pre-initialize surface variables, i.e. setting start- and end-indices
1448!--       at each (j,i)-location. Please note, this does not supersede
1449!--       user-defined initialization of surface quantities.
1450          CALL init_surfaces
1451!
1452!--       Initialization will completely be done by the user
1453          CALL user_init_3d_model
1454
1455          CALL location_message( 'finished', .TRUE. )
1456
1457       ENDIF
1458
1459       CALL location_message( 'initializing statistics, boundary conditions, etc.', &
1460                              .FALSE. )
1461
1462!
1463!--    Bottom boundary
1464       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
1465          u(nzb,:,:) = 0.0_wp
1466          v(nzb,:,:) = 0.0_wp
1467       ENDIF
1468
1469!
1470!--    Apply channel flow boundary condition
1471       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
1472          u(nzt+1,:,:) = 0.0_wp
1473          v(nzt+1,:,:) = 0.0_wp
1474       ENDIF
1475
1476!
1477!--    Calculate virtual potential temperature
1478       IF ( humidity )  vpt = pt * ( 1.0_wp + 0.61_wp * q )
1479
1480!
1481!--    Store initial profiles for output purposes etc.. Please note, in case of
1482!--    initialization of u, v, w, pt, and q via output data derived from larger
1483!--    scale models, data will not be horizontally homogeneous. Actually, a mean
1484!--    profile should be calculated before.   
1485       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
1486       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
1487       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
1488          hom(nzb,1,5,:) = 0.0_wp
1489          hom(nzb,1,6,:) = 0.0_wp
1490       ENDIF
1491       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1492
1493
1494!
1495!--    Store initial salinity profile
1496       IF ( ocean )  THEN
1497          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
1498       ENDIF
1499
1500       IF ( humidity )  THEN
1501!
1502!--       Store initial profile of total water content, virtual potential
1503!--       temperature
1504          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
1505          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
1506!
1507!--       Store initial profile of specific humidity and potential
1508!--       temperature
1509          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
1510             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1511             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1512          ENDIF
1513       ENDIF
1514
1515!
1516!--    Store initial scalar profile
1517       IF ( passive_scalar )  THEN
1518          hom(:,1,121,:) = SPREAD(  s(:,nys,nxl), 2, statistic_regions+1 )
1519       ENDIF
1520
1521!
1522!--    Initialize the random number generators (from numerical recipes)
1523       CALL random_function_ini
1524       
1525       IF ( random_generator == 'random-parallel' )  THEN
1526          CALL init_parallel_random_generator(nx, ny, nys, nyn, nxl, nxr)
1527       ENDIF
1528!
1529!--    Set the reference state to be used in the buoyancy terms (for ocean runs
1530!--    the reference state will be set (overwritten) in init_ocean)
1531       IF ( use_single_reference_value )  THEN
1532          IF (  .NOT.  humidity )  THEN
1533             ref_state(:) = pt_reference
1534          ELSE
1535             ref_state(:) = vpt_reference
1536          ENDIF
1537       ELSE
1538          IF (  .NOT.  humidity )  THEN
1539             ref_state(:) = pt_init(:)
1540          ELSE
1541             ref_state(:) = vpt(:,nys,nxl)
1542          ENDIF
1543       ENDIF
1544
1545!
1546!--    For the moment, vertical velocity is zero
1547       w = 0.0_wp
1548
1549!
1550!--    Initialize array sums (must be defined in first call of pres)
1551       sums = 0.0_wp
1552
1553!
1554!--    In case of iterative solvers, p must get an initial value
1555       IF ( psolver(1:9) == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0_wp
1556
1557!
1558!--    Treating cloud physics, liquid water content and precipitation amount
1559!--    are zero at beginning of the simulation
1560       IF ( cloud_physics )  THEN
1561          ql = 0.0_wp
1562          qc = 0.0_wp
1563
1564          precipitation_amount = 0.0_wp
1565       ENDIF
1566!
1567!--    Impose vortex with vertical axis on the initial velocity profile
1568       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
1569          CALL init_rankine
1570       ENDIF
1571
1572!
1573!--    Impose temperature anomaly (advection test only)
1574       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
1575          CALL init_pt_anomaly
1576       ENDIF
1577
1578!
1579!--    If required, change the surface temperature at the start of the 3D run
1580       IF ( pt_surface_initial_change /= 0.0_wp )  THEN
1581          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
1582       ENDIF
1583
1584!
1585!--    If required, change the surface humidity/scalar at the start of the 3D
1586!--    run
1587       IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )              &
1588          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
1589         
1590       IF ( passive_scalar .AND.  s_surface_initial_change /= 0.0_wp )         &
1591          s(nzb,:,:) = s(nzb,:,:) + s_surface_initial_change
1592       
1593
1594!
1595!--    Initialize old and new time levels.
1596       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1597       pt_p = pt; u_p = u; v_p = v; w_p = w
1598
1599       IF ( humidity  )  THEN
1600          tq_m = 0.0_wp
1601          q_p = q
1602          IF ( cloud_physics  .AND.  microphysics_morrison )  THEN
1603             tqc_m = 0.0_wp
1604             qc_p  = qc
1605             tnc_m = 0.0_wp
1606             nc_p  = nc
1607          ENDIF
1608          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1609             tqr_m = 0.0_wp
1610             qr_p  = qr
1611             tnr_m = 0.0_wp
1612             nr_p  = nr
1613          ENDIF
1614       ENDIF
1615       
1616       IF ( passive_scalar )  THEN
1617          ts_m = 0.0_wp
1618          s_p  = s
1619       ENDIF       
1620
1621       IF ( ocean )  THEN
1622          tsa_m = 0.0_wp
1623          sa_p  = sa
1624       ENDIF
1625       
1626       CALL location_message( 'finished', .TRUE. )
1627
1628    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.         &
1629             TRIM( initializing_actions ) == 'cyclic_fill' )                   &
1630    THEN
1631
1632       CALL location_message( 'initializing in case of restart / cyclic_fill', &
1633                              .FALSE. )
1634!
1635!--    Initialize surface elements and its attributes, e.g. heat- and
1636!--    momentumfluxes, roughness, scaling parameters. As number of surface
1637!--    elements might be different between runs, e.g. in case of cyclic fill,
1638!--    and not all surface elements are read, surface elements need to be
1639!--    initialized before.     
1640       CALL init_surfaces
1641!
1642!--    When reading data for cyclic fill of 3D prerun data files, read
1643!--    some of the global variables from the restart file which are required
1644!--    for initializing the inflow
1645       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1646
1647          DO  i = 0, io_blocks-1
1648             IF ( i == io_group )  THEN
1649                CALL rrd_read_parts_of_global
1650             ENDIF
1651#if defined( __parallel )
1652             CALL MPI_BARRIER( comm2d, ierr )
1653#endif
1654          ENDDO
1655
1656       ENDIF
1657
1658!
1659!--    Read processor specific binary data from restart file
1660       DO  i = 0, io_blocks-1
1661          IF ( i == io_group )  THEN
1662             CALL rrd_local
1663          ENDIF
1664#if defined( __parallel )
1665          CALL MPI_BARRIER( comm2d, ierr )
1666#endif
1667       ENDDO
1668
1669!
1670!--    In case of complex terrain and cyclic fill method as initialization,
1671!--    shift initial data in the vertical direction for each point in the
1672!--    x-y-plane depending on local surface height
1673       IF ( complex_terrain  .AND.                                             &
1674            TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1675          DO  i = nxlg, nxrg
1676             DO  j = nysg, nyng
1677                nz_u_shift = get_topography_top_index_ji( j, i, 'u' )
1678                nz_v_shift = get_topography_top_index_ji( j, i, 'v' )
1679                nz_w_shift = get_topography_top_index_ji( j, i, 'w' )
1680                nz_s_shift = get_topography_top_index_ji( j, i, 's' )
1681
1682                u(nz_u_shift:nzt+1,j,i)  = u(0:nzt+1-nz_u_shift,j,i)               
1683
1684                v(nz_v_shift:nzt+1,j,i)  = v(0:nzt+1-nz_v_shift,j,i)
1685
1686                w(nz_w_shift:nzt+1,j,i)  = w(0:nzt+1-nz_w_shift,j,i)
1687
1688                p(nz_s_shift:nzt+1,j,i)  =  p(0:nzt+1-nz_s_shift,j,i)
1689                pt(nz_s_shift:nzt+1,j,i) = pt(0:nzt+1-nz_s_shift,j,i)
1690             ENDDO
1691          ENDDO
1692       ENDIF
1693
1694!
1695!--    Initialization of the turbulence recycling method
1696       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.               &
1697            turbulent_inflow )  THEN
1698!
1699!--       First store the profiles to be used at the inflow.
1700!--       These profiles are the (temporally) and horizontally averaged vertical
1701!--       profiles from the prerun. Alternatively, prescribed profiles
1702!--       for u,v-components can be used.
1703          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,7) )
1704
1705          IF ( use_prescribed_profile_data )  THEN
1706             mean_inflow_profiles(:,1) = u_init            ! u
1707             mean_inflow_profiles(:,2) = v_init            ! v
1708          ELSE
1709             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
1710             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
1711          ENDIF
1712          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
1713          IF ( humidity )                                                      &
1714             mean_inflow_profiles(:,6) = hom_sum(:,41,0)   ! q
1715          IF ( passive_scalar )                                                &
1716             mean_inflow_profiles(:,7) = hom_sum(:,115,0)   ! s
1717!
1718!--       In case of complex terrain, determine vertical displacement at inflow
1719!--       boundary and adjust mean inflow profiles
1720          IF ( complex_terrain )  THEN
1721             IF ( nxlg <= 0 .AND. nxrg >= 0 .AND. nysg <= 0 .AND. nyng >= 0 )  THEN
1722                nz_u_shift_l = get_topography_top_index_ji( 0, 0, 'u' )
1723                nz_v_shift_l = get_topography_top_index_ji( 0, 0, 'v' )
1724                nz_w_shift_l = get_topography_top_index_ji( 0, 0, 'w' )
1725                nz_s_shift_l = get_topography_top_index_ji( 0, 0, 's' )
1726             ELSE
1727                nz_u_shift_l = 0
1728                nz_v_shift_l = 0
1729                nz_w_shift_l = 0
1730                nz_s_shift_l = 0
1731             ENDIF
1732
1733#if defined( __parallel )
1734             CALL MPI_ALLREDUCE(nz_u_shift_l, nz_u_shift, 1, MPI_INTEGER,      &
1735                                MPI_MAX, comm2d, ierr)
1736             CALL MPI_ALLREDUCE(nz_v_shift_l, nz_v_shift, 1, MPI_INTEGER,      &
1737                                MPI_MAX, comm2d, ierr)
1738             CALL MPI_ALLREDUCE(nz_w_shift_l, nz_w_shift, 1, MPI_INTEGER,      & 
1739                                MPI_MAX, comm2d, ierr)
1740             CALL MPI_ALLREDUCE(nz_s_shift_l, nz_s_shift, 1, MPI_INTEGER,      &
1741                                MPI_MAX, comm2d, ierr)
1742#else
1743             nz_u_shift = nz_u_shift_l
1744             nz_v_shift = nz_v_shift_l
1745             nz_w_shift = nz_w_shift_l
1746             nz_s_shift = nz_s_shift_l
1747#endif
1748
1749             mean_inflow_profiles(:,1) = 0.0_wp
1750             mean_inflow_profiles(nz_u_shift:nzt+1,1) = hom_sum(0:nzt+1-nz_u_shift,1,0)  ! u
1751
1752             mean_inflow_profiles(:,2) = 0.0_wp
1753             mean_inflow_profiles(nz_v_shift:nzt+1,2) = hom_sum(0:nzt+1-nz_v_shift,2,0)  ! v
1754
1755             mean_inflow_profiles(nz_s_shift:nzt+1,4) = hom_sum(0:nzt+1-nz_s_shift,4,0)  ! pt
1756
1757          ENDIF
1758
1759!
1760!--       If necessary, adjust the horizontal flow field to the prescribed
1761!--       profiles
1762          IF ( use_prescribed_profile_data )  THEN
1763             DO  i = nxlg, nxrg
1764                DO  j = nysg, nyng
1765                   DO  k = nzb, nzt+1
1766                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1767                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
1768                   ENDDO
1769                ENDDO
1770             ENDDO
1771          ENDIF
1772
1773!
1774!--       Use these mean profiles at the inflow (provided that Dirichlet
1775!--       conditions are used)
1776          IF ( inflow_l )  THEN
1777             DO  j = nysg, nyng
1778                DO  k = nzb, nzt+1
1779                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1780                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
1781                   w(k,j,nxlg:-1)  = 0.0_wp
1782                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
1783                   IF ( humidity )                                             &
1784                      q(k,j,nxlg:-1)  = mean_inflow_profiles(k,6)
1785                   IF ( passive_scalar )                                       &
1786                      s(k,j,nxlg:-1)  = mean_inflow_profiles(k,7)                     
1787                ENDDO
1788             ENDDO
1789          ENDIF
1790
1791!
1792!--       Calculate the damping factors to be used at the inflow. For a
1793!--       turbulent inflow the turbulent fluctuations have to be limited
1794!--       vertically because otherwise the turbulent inflow layer will grow
1795!--       in time.
1796          IF ( inflow_damping_height == 9999999.9_wp )  THEN
1797!
1798!--          Default: use the inversion height calculated by the prerun; if
1799!--          this is zero, inflow_damping_height must be explicitly
1800!--          specified.
1801             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0_wp )  THEN
1802                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1803             ELSE
1804                WRITE( message_string, * ) 'inflow_damping_height must be ',   &
1805                     'explicitly specified because&the inversion height ',     &
1806                     'calculated by the prerun is zero.'
1807                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
1808             ENDIF
1809
1810          ENDIF
1811
1812          IF ( inflow_damping_width == 9999999.9_wp )  THEN
1813!
1814!--          Default for the transition range: one tenth of the undamped
1815!--          layer
1816             inflow_damping_width = 0.1_wp * inflow_damping_height
1817
1818          ENDIF
1819
1820          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
1821
1822          DO  k = nzb, nzt+1
1823
1824             IF ( zu(k) <= inflow_damping_height )  THEN
1825                inflow_damping_factor(k) = 1.0_wp
1826             ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) )  THEN
1827                inflow_damping_factor(k) = 1.0_wp -                            &
1828                                           ( zu(k) - inflow_damping_height ) / &
1829                                           inflow_damping_width
1830             ELSE
1831                inflow_damping_factor(k) = 0.0_wp
1832             ENDIF
1833
1834          ENDDO
1835
1836       ENDIF
1837
1838!
1839!--    Inside buildings set velocities back to zero
1840       IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND.                &
1841            topography /= 'flat' )  THEN
1842!
1843!--       Inside buildings set velocities back to zero.
1844!--       Other scalars (pt, q, s, p, sa, ...) are ignored at present,
1845!--       maybe revise later.
1846          DO  i = nxlg, nxrg
1847             DO  j = nysg, nyng
1848                DO  k = nzb, nzt
1849                   u(k,j,i)     = MERGE( u(k,j,i), 0.0_wp,                     &
1850                                         BTEST( wall_flags_0(k,j,i), 1 ) )
1851                   v(k,j,i)     = MERGE( v(k,j,i), 0.0_wp,                     &
1852                                         BTEST( wall_flags_0(k,j,i), 2 ) )
1853                   w(k,j,i)     = MERGE( w(k,j,i), 0.0_wp,                     &
1854                                         BTEST( wall_flags_0(k,j,i), 3 ) )
1855                ENDDO
1856             ENDDO
1857          ENDDO
1858
1859       ENDIF
1860
1861!
1862!--    Calculate initial temperature field and other constants used in case
1863!--    of a sloping surface
1864       IF ( sloping_surface )  CALL init_slope
1865
1866!
1867!--    Initialize new time levels (only done in order to set boundary values
1868!--    including ghost points)
1869       pt_p = pt; u_p = u; v_p = v; w_p = w
1870       IF ( humidity )  THEN
1871          q_p = q
1872          IF ( cloud_physics  .AND.  microphysics_morrison )  THEN
1873             qc_p = qc
1874             nc_p = nc
1875          ENDIF
1876          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1877             qr_p = qr
1878             nr_p = nr
1879          ENDIF
1880       ENDIF
1881       IF ( passive_scalar )  s_p  = s
1882       IF ( ocean          )  sa_p = sa
1883
1884!
1885!--    Allthough tendency arrays are set in prognostic_equations, they have
1886!--    have to be predefined here because they are used (but multiplied with 0)
1887!--    there before they are set.
1888       tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1889       IF ( humidity )  THEN
1890          tq_m = 0.0_wp
1891          IF ( cloud_physics  .AND.  microphysics_morrison )  THEN
1892             tqc_m = 0.0_wp
1893             tnc_m = 0.0_wp
1894          ENDIF
1895          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1896             tqr_m = 0.0_wp
1897             tnr_m = 0.0_wp
1898          ENDIF
1899       ENDIF
1900       IF ( passive_scalar )  ts_m  = 0.0_wp
1901       IF ( ocean          )  tsa_m = 0.0_wp
1902!
1903!--    Initialize synthetic turbulence generator in case of restart.
1904       IF ( TRIM( initializing_actions ) == 'read_restart_data'  .AND.         &
1905            use_syn_turb_gen )  CALL stg_init
1906
1907       CALL location_message( 'finished', .TRUE. )
1908
1909    ELSE
1910!
1911!--    Actually this part of the programm should not be reached
1912       message_string = 'unknown initializing problem'
1913       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
1914    ENDIF
1915
1916!
1917!-- Initialize TKE, Kh and Km
1918    CALL tcm_init
1919
1920
1921    IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1922!
1923!--    Initialize old timelevels needed for radiation boundary conditions
1924       IF ( outflow_l )  THEN
1925          u_m_l(:,:,:) = u(:,:,1:2)
1926          v_m_l(:,:,:) = v(:,:,0:1)
1927          w_m_l(:,:,:) = w(:,:,0:1)
1928       ENDIF
1929       IF ( outflow_r )  THEN
1930          u_m_r(:,:,:) = u(:,:,nx-1:nx)
1931          v_m_r(:,:,:) = v(:,:,nx-1:nx)
1932          w_m_r(:,:,:) = w(:,:,nx-1:nx)
1933       ENDIF
1934       IF ( outflow_s )  THEN
1935          u_m_s(:,:,:) = u(:,0:1,:)
1936          v_m_s(:,:,:) = v(:,1:2,:)
1937          w_m_s(:,:,:) = w(:,0:1,:)
1938       ENDIF
1939       IF ( outflow_n )  THEN
1940          u_m_n(:,:,:) = u(:,ny-1:ny,:)
1941          v_m_n(:,:,:) = v(:,ny-1:ny,:)
1942          w_m_n(:,:,:) = w(:,ny-1:ny,:)
1943       ENDIF
1944       
1945    ENDIF
1946
1947!
1948!-- Calculate the initial volume flow at the right and north boundary
1949    IF ( conserve_volume_flow )  THEN
1950
1951       IF ( use_prescribed_profile_data )  THEN
1952
1953          volume_flow_initial_l = 0.0_wp
1954          volume_flow_area_l    = 0.0_wp
1955
1956          IF ( nxr == nx )  THEN
1957             DO  j = nys, nyn
1958                DO  k = nzb+1, nzt
1959                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1960                                              u_init(k) * dzw(k)               &
1961                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1962                                              BTEST( wall_flags_0(k,j,nxr), 1 )&
1963                                            )
1964
1965                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
1966                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1967                                              BTEST( wall_flags_0(k,j,nxr), 1 )&
1968                                            )
1969                ENDDO
1970             ENDDO
1971          ENDIF
1972         
1973          IF ( nyn == ny )  THEN
1974             DO  i = nxl, nxr
1975                DO  k = nzb+1, nzt
1976                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1977                                              v_init(k) * dzw(k)               &       
1978                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1979                                              BTEST( wall_flags_0(k,nyn,i), 2 )&
1980                                            )
1981                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
1982                                     * MERGE( 1.0_wp, 0.0_wp,                  &
1983                                              BTEST( wall_flags_0(k,nyn,i), 2 )&
1984                                            )
1985                ENDDO
1986             ENDDO
1987          ENDIF
1988
1989#if defined( __parallel )
1990          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1991                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1992          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1993                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1994
1995#else
1996          volume_flow_initial = volume_flow_initial_l
1997          volume_flow_area    = volume_flow_area_l
1998#endif 
1999
2000       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
2001
2002          volume_flow_initial_l = 0.0_wp
2003          volume_flow_area_l    = 0.0_wp
2004
2005          IF ( nxr == nx )  THEN
2006             DO  j = nys, nyn
2007                DO  k = nzb+1, nzt
2008                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
2009                                              hom_sum(k,1,0) * dzw(k)          &
2010                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2011                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
2012                                            )
2013                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
2014                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2015                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
2016                                            )
2017                ENDDO
2018             ENDDO
2019          ENDIF
2020         
2021          IF ( nyn == ny )  THEN
2022             DO  i = nxl, nxr
2023                DO  k = nzb+1, nzt
2024                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
2025                                              hom_sum(k,2,0) * dzw(k)          &       
2026                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2027                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
2028                                            )
2029                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
2030                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2031                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
2032                                            )
2033                ENDDO
2034             ENDDO
2035          ENDIF
2036
2037#if defined( __parallel )
2038          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
2039                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
2040          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
2041                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
2042
2043#else
2044          volume_flow_initial = volume_flow_initial_l
2045          volume_flow_area    = volume_flow_area_l
2046#endif 
2047
2048       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
2049
2050          volume_flow_initial_l = 0.0_wp
2051          volume_flow_area_l    = 0.0_wp
2052
2053          IF ( nxr == nx )  THEN
2054             DO  j = nys, nyn
2055                DO  k = nzb+1, nzt
2056                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
2057                                              u(k,j,nx) * dzw(k)               &
2058                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2059                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
2060                                            )
2061                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)   &
2062                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2063                                              BTEST( wall_flags_0(k,j,nx), 1 ) &
2064                                            )
2065                ENDDO
2066             ENDDO
2067          ENDIF
2068         
2069          IF ( nyn == ny )  THEN
2070             DO  i = nxl, nxr
2071                DO  k = nzb+1, nzt
2072                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
2073                                              v(k,ny,i) * dzw(k)               &       
2074                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2075                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
2076                                            )
2077                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)   &       
2078                                     * MERGE( 1.0_wp, 0.0_wp,                  &
2079                                              BTEST( wall_flags_0(k,ny,i), 2 ) &
2080                                            )
2081                ENDDO
2082             ENDDO
2083          ENDIF
2084
2085#if defined( __parallel )
2086          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
2087                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
2088          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
2089                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
2090
2091#else
2092          volume_flow_initial = volume_flow_initial_l
2093          volume_flow_area    = volume_flow_area_l
2094#endif 
2095
2096       ENDIF
2097
2098!
2099!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated
2100!--    from u|v_bulk instead
2101       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
2102          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
2103          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
2104       ENDIF
2105
2106    ENDIF
2107!
2108!-- Finally, if random_heatflux is set, disturb shf at horizontal
2109!-- surfaces. Actually, this should be done in surface_mod, where all other
2110!-- initializations of surface quantities are done. However, this
2111!-- would create a ring dependency, hence, it is done here. Maybe delete
2112!-- disturb_heatflux and tranfer the respective code directly into the
2113!-- initialization in surface_mod.         
2114    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
2115         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
2116 
2117       IF ( use_surface_fluxes  .AND.  constant_heatflux  .AND.                &
2118            random_heatflux )  THEN
2119          IF ( surf_def_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_def_h(0) )
2120          IF ( surf_lsm_h%ns    >= 1 )  CALL disturb_heatflux( surf_lsm_h    )
2121          IF ( surf_usm_h%ns    >= 1 )  CALL disturb_heatflux( surf_usm_h    )
2122       ENDIF
2123    ENDIF
2124
2125!
2126!-- Before initializing further modules, compute total sum of active mask
2127!-- grid points and the mean surface level height for each statistic region.
2128!-- ngp_2dh: number of grid points of a horizontal cross section through the
2129!--          total domain
2130!-- ngp_3d:  number of grid points of the total domain
2131    ngp_2dh_outer_l   = 0
2132    ngp_2dh_outer     = 0
2133    ngp_2dh_s_inner_l = 0
2134    ngp_2dh_s_inner   = 0
2135    ngp_2dh_l         = 0
2136    ngp_2dh           = 0
2137    ngp_3d_inner_l    = 0.0_wp
2138    ngp_3d_inner      = 0
2139    ngp_3d            = 0
2140    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
2141
2142    mean_surface_level_height   = 0.0_wp
2143    mean_surface_level_height_l = 0.0_wp
2144!
2145!-- Pre-set masks for regional statistics. Default is the total model domain.
2146!-- Ghost points are excluded because counting values at the ghost boundaries
2147!-- would bias the statistics
2148    rmask = 1.0_wp
2149    rmask(:,nxlg:nxl-1,:) = 0.0_wp;  rmask(:,nxr+1:nxrg,:) = 0.0_wp
2150    rmask(nysg:nys-1,:,:) = 0.0_wp;  rmask(nyn+1:nyng,:,:) = 0.0_wp
2151!
2152!-- User-defined initializing actions
2153    CALL user_init
2154!
2155!-- To do: New concept for these non-topography grid points!
2156    DO  sr = 0, statistic_regions
2157       DO  i = nxl, nxr
2158          DO  j = nys, nyn
2159             IF ( rmask(j,i,sr) == 1.0_wp )  THEN
2160!
2161!--             All xy-grid points
2162                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
2163!
2164!--             Determine mean surface-level height. In case of downward-
2165!--             facing walls are present, more than one surface level exist.
2166!--             In this case, use the lowest surface-level height.
2167                IF ( surf_def_h(0)%start_index(j,i) <=                         &
2168                     surf_def_h(0)%end_index(j,i) )  THEN
2169                   m = surf_def_h(0)%start_index(j,i)
2170                   k = surf_def_h(0)%k(m)
2171                   mean_surface_level_height_l(sr) =                           &
2172                                       mean_surface_level_height_l(sr) + zw(k-1)
2173                ENDIF
2174                IF ( surf_lsm_h%start_index(j,i) <=                            &
2175                     surf_lsm_h%end_index(j,i) )  THEN
2176                   m = surf_lsm_h%start_index(j,i)
2177                   k = surf_lsm_h%k(m)
2178                   mean_surface_level_height_l(sr) =                           &
2179                                       mean_surface_level_height_l(sr) + zw(k-1)
2180                ENDIF
2181                IF ( surf_usm_h%start_index(j,i) <=                            &
2182                     surf_usm_h%end_index(j,i) )  THEN
2183                   m = surf_usm_h%start_index(j,i)
2184                   k = surf_usm_h%k(m)
2185                   mean_surface_level_height_l(sr) =                           &
2186                                       mean_surface_level_height_l(sr) + zw(k-1)
2187                ENDIF
2188
2189                k_surf = k - 1
2190
2191                DO  k = nzb, nzt+1
2192!
2193!--                xy-grid points above topography
2194                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr)     +         &
2195                                  MERGE( 1, 0, BTEST( wall_flags_0(k,j,i), 24 ) )
2196
2197                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) +         &
2198                                  MERGE( 1, 0, BTEST( wall_flags_0(k,j,i), 22 ) )
2199
2200                ENDDO
2201!
2202!--             All grid points of the total domain above topography
2203                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + ( nz - k_surf + 2 )
2204
2205
2206
2207             ENDIF
2208          ENDDO
2209       ENDDO
2210    ENDDO
2211!
2212!-- Initialize arrays encompassing number of grid-points in inner and outer
2213!-- domains, statistic regions, etc. Mainly used for horizontal averaging
2214!-- of turbulence statistics. Please note, user_init must be called before
2215!-- doing this.   
2216    sr = statistic_regions + 1
2217#if defined( __parallel )
2218    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2219    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,    &
2220                        comm2d, ierr )
2221    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2222    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr,   &
2223                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
2224    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2225    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),          &
2226                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
2227    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2228    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL,  &
2229                        MPI_SUM, comm2d, ierr )
2230    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
2231    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2232    CALL MPI_ALLREDUCE( mean_surface_level_height_l(0),                        &
2233                        mean_surface_level_height(0), sr, MPI_REAL,            &
2234                        MPI_SUM, comm2d, ierr )
2235    mean_surface_level_height = mean_surface_level_height / REAL( ngp_2dh )
2236#else
2237    ngp_2dh         = ngp_2dh_l
2238    ngp_2dh_outer   = ngp_2dh_outer_l
2239    ngp_2dh_s_inner = ngp_2dh_s_inner_l
2240    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
2241    mean_surface_level_height = mean_surface_level_height_l / REAL( ngp_2dh_l )
2242#endif
2243
2244    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * &
2245             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
2246
2247!
2248!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
2249!-- buoyancy, etc. A zero value will occur for cases where all grid points of
2250!-- the respective subdomain lie below the surface topography
2251    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   ) 
2252    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )),             &
2253                           ngp_3d_inner(:) )
2254    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) 
2255
2256    DEALLOCATE( mean_surface_level_height_l, ngp_2dh_l, ngp_2dh_outer_l,       &
2257                ngp_3d_inner_l, ngp_3d_inner_tmp )
2258
2259!
2260!-- Initialize surface forcing corresponding to large-scale forcing. Therein,
2261!-- initialize heat-fluxes, etc. via datatype. Revise it later!
2262    IF ( large_scale_forcing .AND. lsf_surf )  THEN
2263       IF ( use_surface_fluxes  .AND.  constant_heatflux )  THEN
2264          CALL ls_forcing_surf ( simulated_time )
2265       ENDIF
2266    ENDIF
2267!
2268!-- Initialize quantities for special advections schemes
2269    CALL init_advec
2270
2271!
2272!-- Impose random perturbation on the horizontal velocity field and then
2273!-- remove the divergences from the velocity field at the initial stage
2274    IF ( create_disturbances  .AND.  disturbance_energy_limit /= 0.0_wp  .AND. &
2275         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
2276         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
2277
2278       CALL location_message( 'creating initial disturbances', .FALSE. )
2279       CALL disturb_field( 'u', tend, u )
2280       CALL disturb_field( 'v', tend, v )
2281       CALL location_message( 'finished', .TRUE. )
2282
2283       CALL location_message( 'calling pressure solver', .FALSE. )
2284       n_sor = nsor_ini
2285       CALL pres
2286       n_sor = nsor
2287       CALL location_message( 'finished', .TRUE. )
2288
2289    ENDIF
2290
2291!
2292!-- If required, initialize quantities needed for the plant canopy model
2293    IF ( plant_canopy )  THEN
2294       CALL location_message( 'initializing plant canopy model', .FALSE. )   
2295       CALL pcm_init
2296       CALL location_message( 'finished', .TRUE. )
2297    ENDIF
2298
2299!
2300!-- If required, initialize dvrp-software
2301    IF ( dt_dvrp /= 9999999.9_wp )  CALL init_dvrp
2302
2303    IF ( ocean )  THEN
2304!
2305!--    Initialize quantities needed for the ocean model
2306       CALL init_ocean
2307
2308    ELSE
2309!
2310!--    Initialize quantities for handling cloud physics
2311!--    This routine must be called before lpm_init, because
2312!--    otherwise, array pt_d_t, needed in data_output_dvrp (called by
2313!--    lpm_init) is not defined.
2314       CALL init_cloud_physics
2315!
2316!--    Initialize bulk cloud microphysics
2317       CALL microphysics_init
2318    ENDIF
2319
2320!
2321!-- If required, initialize particles
2322    IF ( particle_advection )  CALL lpm_init
2323
2324!
2325!-- If required, initialize quantities needed for the LSM
2326    IF ( land_surface )  THEN
2327       CALL location_message( 'initializing land surface model', .FALSE. )
2328       CALL lsm_init
2329       CALL location_message( 'finished', .TRUE. )
2330    ENDIF
2331
2332!
2333!-- If required, allocate USM and LSM surfaces
2334    IF ( urban_surface )  THEN
2335       CALL location_message( 'initializing and allocating urban surfaces', .FALSE. )
2336       CALL usm_allocate_surface
2337       CALL location_message( 'finished', .TRUE. )
2338    ENDIF
2339!
2340!-- If required, initialize urban surface model
2341    IF ( urban_surface )  THEN
2342       CALL location_message( 'initializing urban surface model', .FALSE. )
2343       CALL usm_init_urban_surface
2344       CALL location_message( 'finished', .TRUE. )
2345    ENDIF
2346
2347!
2348!-- Initialize surface layer (done after LSM as roughness length are required
2349!-- for initialization
2350    IF ( constant_flux_layer )  THEN
2351       CALL location_message( 'initializing surface layer', .FALSE. )
2352       CALL init_surface_layer_fluxes
2353       CALL location_message( 'finished', .TRUE. )
2354    ENDIF
2355
2356!
2357!-- If required, set chemical emissions
2358!-- (todo(FK): This should later on be CALLed time-dependently in init_3d_model)
2359    IF ( air_chemistry ) THEN
2360       CALL chem_emissions
2361    ENDIF
2362
2363!
2364!-- Initialize radiation processes
2365    IF ( radiation )  THEN
2366!
2367!--    If required, initialize radiation interactions between surfaces
2368!--    via sky-view factors. This must be done before radiation is initialized.
2369       IF ( radiation_interactions )  CALL radiation_interaction_init
2370
2371!
2372!--    Initialize radiation model
2373       CALL location_message( 'initializing radiation model', .FALSE. )
2374       CALL radiation_init
2375       CALL location_message( 'finished', .TRUE. )
2376
2377!
2378!--    Find all discretized apparent solar positions for radiation interaction.
2379!--    This must be done after radiation_init.
2380       IF ( radiation_interactions )  CALL radiation_presimulate_solar_pos
2381
2382!
2383!--    If required, read or calculate and write out the SVF
2384       IF ( radiation_interactions .AND. read_svf)  THEN
2385!
2386!--       Read sky-view factors and further required data from file
2387          CALL location_message( '    Start reading SVF from file', .FALSE. )
2388          CALL radiation_read_svf()
2389          CALL location_message( '    Reading SVF from file has finished', .TRUE. )
2390
2391       ELSEIF ( radiation_interactions .AND. .NOT. read_svf)  THEN
2392!
2393!--       calculate SFV and CSF
2394          CALL location_message( '    Start calculation of SVF', .FALSE. )
2395          CALL radiation_calc_svf()
2396          CALL location_message( '    Calculation of SVF has finished', .TRUE. )
2397       ENDIF
2398
2399       IF ( radiation_interactions .AND. write_svf)  THEN
2400!
2401!--       Write svf, csf svfsurf and csfsurf data to file
2402          CALL location_message( '    Start writing SVF in file', .FALSE. )
2403          CALL radiation_write_svf()
2404          CALL location_message( '    Writing SVF in file has finished', .TRUE. )
2405       ENDIF
2406
2407!
2408!--    Adjust radiative fluxes. In case of urban and land surfaces, also
2409!--    call an initial interaction.
2410       CALL radiation_control 
2411       IF ( radiation_interactions )  THEN
2412          CALL radiation_interaction
2413       ENDIF
2414    ENDIF
2415   
2416!
2417!-- Temporary solution to add LSM and radiation time series to the default
2418!-- output
2419    IF ( land_surface  .OR.  radiation )  THEN
2420       IF ( TRIM( radiation_scheme ) == 'rrtmg' )  THEN
2421          dots_num = dots_num + 15
2422       ELSE
2423          dots_num = dots_num + 11
2424       ENDIF
2425    ENDIF
2426   
2427
2428
2429!
2430!-- If required, initialize quantities needed for the wind turbine model
2431    IF ( wind_turbine )  THEN
2432       CALL location_message( 'initializing wind turbine model', .FALSE. )
2433       CALL wtm_init
2434       CALL location_message( 'finished', .TRUE. )
2435    ENDIF
2436
2437!
2438!-- If required, initialize quantities needed for the gust module
2439    IF ( gust_module_enabled )  THEN
2440       CALL gust_init( dots_label, dots_unit, dots_num, dots_max )
2441    ENDIF
2442
2443!
2444!-- Initialize the ws-scheme.   
2445    IF ( ws_scheme_sca .OR. ws_scheme_mom )  CALL ws_init       
2446
2447!
2448!-- Setting weighting factors for calculation of perturbation pressure
2449!-- and turbulent quantities from the RK substeps
2450    IF ( TRIM(timestep_scheme) == 'runge-kutta-3' )  THEN      ! for RK3-method
2451
2452       weight_substep(1) = 1._wp/6._wp
2453       weight_substep(2) = 3._wp/10._wp
2454       weight_substep(3) = 8._wp/15._wp
2455
2456       weight_pres(1)    = 1._wp/3._wp
2457       weight_pres(2)    = 5._wp/12._wp
2458       weight_pres(3)    = 1._wp/4._wp
2459
2460    ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' )  THEN  ! for RK2-method
2461
2462       weight_substep(1) = 1._wp/2._wp
2463       weight_substep(2) = 1._wp/2._wp
2464         
2465       weight_pres(1)    = 1._wp/2._wp
2466       weight_pres(2)    = 1._wp/2._wp       
2467
2468    ELSE                                     ! for Euler-method
2469
2470       weight_substep(1) = 1.0_wp     
2471       weight_pres(1)    = 1.0_wp                   
2472
2473    ENDIF
2474
2475!
2476!-- Initialize Rayleigh damping factors
2477    rdf    = 0.0_wp
2478    rdf_sc = 0.0_wp
2479    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
2480       IF (  .NOT.  ocean )  THEN
2481          DO  k = nzb+1, nzt
2482             IF ( zu(k) >= rayleigh_damping_height )  THEN
2483                rdf(k) = rayleigh_damping_factor *                             &
2484                      ( SIN( pi * 0.5_wp * ( zu(k) - rayleigh_damping_height ) &
2485                             / ( zu(nzt) - rayleigh_damping_height ) )         &
2486                      )**2
2487             ENDIF
2488          ENDDO
2489       ELSE
2490          DO  k = nzt, nzb+1, -1
2491             IF ( zu(k) <= rayleigh_damping_height )  THEN
2492                rdf(k) = rayleigh_damping_factor *                             &
2493                      ( SIN( pi * 0.5_wp * ( rayleigh_damping_height - zu(k) ) &
2494                             / ( rayleigh_damping_height - zu(nzb+1) ) )       &
2495                      )**2
2496             ENDIF
2497          ENDDO
2498       ENDIF
2499    ENDIF
2500    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
2501
2502!
2503!-- Initialize the starting level and the vertical smoothing factor used for
2504!-- the external pressure gradient
2505    dp_smooth_factor = 1.0_wp
2506    IF ( dp_external )  THEN
2507!
2508!--    Set the starting level dp_level_ind_b only if it has not been set before
2509!--    (e.g. in init_grid).
2510       IF ( dp_level_ind_b == 0 )  THEN
2511          ind_array = MINLOC( ABS( dp_level_b - zu ) )
2512          dp_level_ind_b = ind_array(1) - 1 + nzb 
2513                                        ! MINLOC uses lower array bound 1
2514       ENDIF
2515       IF ( dp_smooth )  THEN
2516          dp_smooth_factor(:dp_level_ind_b) = 0.0_wp
2517          DO  k = dp_level_ind_b+1, nzt
2518             dp_smooth_factor(k) = 0.5_wp * ( 1.0_wp + SIN( pi *               &
2519                        ( REAL( k - dp_level_ind_b, KIND=wp ) /                &
2520                          REAL( nzt - dp_level_ind_b, KIND=wp ) - 0.5_wp ) ) )
2521          ENDDO
2522       ENDIF
2523    ENDIF
2524
2525!
2526!-- Initialize damping zone for the potential temperature in case of
2527!-- non-cyclic lateral boundaries. The damping zone has the maximum value
2528!-- at the inflow boundary and decreases to zero at pt_damping_width.
2529    ptdf_x = 0.0_wp
2530    ptdf_y = 0.0_wp
2531    IF ( bc_lr_dirrad )  THEN
2532       DO  i = nxl, nxr
2533          IF ( ( i * dx ) < pt_damping_width )  THEN
2534             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5_wp *              &
2535                            REAL( pt_damping_width - i * dx, KIND=wp ) / (     &
2536                            REAL( pt_damping_width, KIND=wp ) ) ) )**2 
2537          ENDIF
2538       ENDDO
2539    ELSEIF ( bc_lr_raddir )  THEN
2540       DO  i = nxl, nxr
2541          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
2542             ptdf_x(i) = pt_damping_factor *                                   &
2543                         SIN( pi * 0.5_wp *                                    &
2544                                 ( ( i - nx ) * dx + pt_damping_width ) /      &
2545                                 REAL( pt_damping_width, KIND=wp ) )**2
2546          ENDIF
2547       ENDDO 
2548    ELSEIF ( bc_ns_dirrad )  THEN
2549       DO  j = nys, nyn
2550          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
2551             ptdf_y(j) = pt_damping_factor *                                   &
2552                         SIN( pi * 0.5_wp *                                    &
2553                                 ( ( j - ny ) * dy + pt_damping_width ) /      &
2554                                 REAL( pt_damping_width, KIND=wp ) )**2
2555          ENDIF
2556       ENDDO 
2557    ELSEIF ( bc_ns_raddir )  THEN
2558       DO  j = nys, nyn
2559          IF ( ( j * dy ) < pt_damping_width )  THEN
2560             ptdf_y(j) = pt_damping_factor *                                   &
2561                         SIN( pi * 0.5_wp *                                    &
2562                                ( pt_damping_width - j * dy ) /                &
2563                                REAL( pt_damping_width, KIND=wp ) )**2
2564          ENDIF
2565       ENDDO
2566    ENDIF
2567!
2568!-- Check if maximum number of allowed timeseries is exceeded
2569    IF ( dots_num > dots_max )  THEN
2570       WRITE( message_string, * ) 'number of time series quantities exceeds',  &
2571                                  ' its maximum of dots_max = ', dots_max,     &
2572                                  ' &Please increase dots_max in modules.f90.'
2573       CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 )   
2574    ENDIF
2575
2576!
2577!-- Input binary data file is not needed anymore. This line must be placed
2578!-- after call of user_init!
2579    CALL close_file( 13 )
2580!
2581!-- In case of nesting, put an barrier to assure that all parent and child
2582!-- domains finished initialization.
2583#if defined( __parallel )
2584    IF ( nested_run )  CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
2585#endif
2586
2587
2588    CALL location_message( 'leaving init_3d_model', .TRUE. )
2589
2590 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.