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

Last change on this file since 3542 was 3525, checked in by kanani, 6 years ago

Changes related to clean-up of biometeorology (by dom_dwd_user)

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