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

Last change on this file since 3042 was 3042, checked in by schwenkel, 6 years ago

Changed the name specific humidity to mixing ratio

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