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

Last change on this file since 3299 was 3298, checked in by kanani, 6 years ago

Merge chemistry branch at r3297 to trunk

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