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

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

Bugfix, no initial masking of wind velocity at first prognostic grid level in case of land- or urban-surface spin-up.

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