source: palm/tags/release-5.0/SOURCE/init_3d_model.f90 @ 2977

Last change on this file since 2977 was 2701, checked in by suehring, 6 years ago

changes from last commit documented

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