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

Last change on this file since 3469 was 3469, checked in by kanani, 5 years ago

Implement indoor climate and energy demand model

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