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

Last change on this file since 3552 was 3547, checked in by suehring, 6 years ago

variable description added some routines

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