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

Last change on this file since 2903 was 2894, checked in by Giersch, 7 years ago

Reading/Writing? data in case of restart runs revised

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