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

Last change on this file since 3123 was 3083, checked in by gronemeier, 6 years ago

merge with branch rans

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