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

Last change on this file since 3254 was 3241, checked in by raasch, 6 years ago

various changes to avoid compiler warnings (mainly removal of unused variables)

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