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

Last change on this file since 2750 was 2746, checked in by suehring, 7 years ago

Read information from statitic driver for resolved vegetation independently from land- or urban-surface model

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