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

Last change on this file since 2906 was 2906, checked in by Giersch, 6 years ago

new procedure for reading/writing svf data, initialization of local variable ids

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