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

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

Code adjusted according to coding standards, renamed namelists, error messages revised until PA0347, output CASE 108 disabled

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