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

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

Bugfix, missing initialization of surface attributes in case of inifor-initialization branch

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