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

Last change on this file since 2773 was 2766, checked in by kanani, 6 years ago

Removal of chem directive, plus minor changes

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