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

Last change on this file since 3422 was 3421, checked in by gronemeier, 6 years ago

new surface-data output; renamed output variables (pt to theta, rho_air to rho, rho_ocean to rho_sea_water)

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