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

Last change on this file since 1992 was 1992, checked in by suehring, 8 years ago

Prescribing scalar flux at model top; several bugfixes concering data output of scalars and output of flight data

  • Property svn:keywords set to Id
  • Property svn:mergeinfo set to False
    /palm/branches/forwind/SOURCE/init_3d_model.f901564-1913
File size: 74.1 KB
Line 
1!> @file init_3d_model.f90
2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2016 Leibniz Universitaet Hannover
17!--------------------------------------------------------------------------------!
18!
19! Current revisions:
20! ------------------
21! Initializaton of scalarflux at model top
22! Bugfixes in initialization of surface and top salinity flux, top scalar and
23! humidity fluxes
24!
25! Former revisions:
26! -----------------
27! $Id: init_3d_model.f90 1992 2016-08-12 15:14:59Z suehring $
28!
29! 1960 2016-07-12 16:34:24Z suehring
30! Separate humidity and passive scalar
31! Increase dimension for mean_inflow_profiles
32! Remove inadvertent write-statement
33! Bugfix, large-scale forcing is still not implemented for passive scalars
34!
35! 1957 2016-07-07 10:43:48Z suehring
36! flight module added
37!
38! 1920 2016-05-30 10:50:15Z suehring
39! Initialize us with very small number to avoid segmentation fault during
40! calculation of Obukhov length
41!
42! 1918 2016-05-27 14:35:57Z raasch
43! intermediate_timestep_count is set 0 instead 1 for first call of pres,
44! bugfix: initialization of local sum arrays are moved to the beginning of the
45!         routine because otherwise results from pres are overwritten
46!
47! 1914 2016-05-26 14:44:07Z witha
48! Added initialization of the wind turbine model
49!
50! 1878 2016-04-19 12:30:36Z hellstea
51! The zeroth element of weight_pres removed as unnecessary
52!
53! 1849 2016-04-08 11:33:18Z hoffmann
54! Adapted for modularization of microphysics.
55! precipitation_amount, precipitation_rate, prr moved to arrays_3d.
56! Initialization of nc_1d, nr_1d, pt_1d, qc_1d, qr_1d, q_1d moved to
57! microphysics_init.
58!
59! 1845 2016-04-08 08:29:13Z raasch
60! nzb_2d replaced by nzb_u|v_inner
61!
62! 1833 2016-04-07 14:23:03Z raasch
63! initialization of spectra quantities moved to spectra_mod
64!
65! 1831 2016-04-07 13:15:51Z hoffmann
66! turbulence renamed collision_turbulence
67!
68! 1826 2016-04-07 12:01:39Z maronga
69! Renamed radiation calls.
70! Renamed canopy model calls.
71!
72! 1822 2016-04-07 07:49:42Z hoffmann
73! icloud_scheme replaced by microphysics_*
74!
75! 1817 2016-04-06 15:44:20Z maronga
76! Renamed lsm calls.
77!
78! 1815 2016-04-06 13:49:59Z raasch
79! zero-settings for velocities inside topography re-activated (was deactivated
80! in r1762)
81!
82! 1788 2016-03-10 11:01:04Z maronga
83! Added z0q.
84! Syntax layout improved.
85!
86! 1783 2016-03-06 18:36:17Z raasch
87! netcdf module name changed + related changes
88!
89! 1764 2016-02-28 12:45:19Z raasch
90! bugfix: increase size of volume_flow_area_l and volume_flow_initial_l by 1
91!
92! 1762 2016-02-25 12:31:13Z hellstea
93! Introduction of nested domain feature
94!
95! 1738 2015-12-18 13:56:05Z raasch
96! calculate mean surface level height for each statistic region
97!
98! 1734 2015-12-02 12:17:12Z raasch
99! no initial disturbances in case that the disturbance energy limit has been
100! set zero
101!
102! 1707 2015-11-02 15:24:52Z maronga
103! Bugfix: transfer of Richardson number from 1D model to Obukhov length caused
104! devision by zero in neutral stratification
105!
106! 1691 2015-10-26 16:17:44Z maronga
107! Call to init_surface_layer added. rif is replaced by ol and zeta.
108!
109! 1682 2015-10-07 23:56:08Z knoop
110! Code annotations made doxygen readable
111!
112! 1615 2015-07-08 18:49:19Z suehring
113! Enable turbulent inflow for passive_scalar and humidity
114!
115! 1585 2015-04-30 07:05:52Z maronga
116! Initialization of radiation code is now done after LSM initializtion
117!
118! 1575 2015-03-27 09:56:27Z raasch
119! adjustments for psolver-queries
120!
121! 1551 2015-03-03 14:18:16Z maronga
122! Allocation of land surface arrays is now done in the subroutine lsm_init_arrays,
123! which is part of land_surface_model.
124!
125! 1507 2014-12-10 12:14:18Z suehring
126! Bugfix: set horizontal velocity components to zero inside topography
127!
128! 1496 2014-12-02 17:25:50Z maronga
129! Added initialization of the land surface and radiation schemes
130!
131! 1484 2014-10-21 10:53:05Z kanani
132! Changes due to new module structure of the plant canopy model:
133! canopy-related initialization (e.g. lad and canopy_heat_flux) moved to new
134! subroutine init_plant_canopy within the module plant_canopy_model_mod,
135! call of subroutine init_plant_canopy added.
136!
137! 1431 2014-07-15 14:47:17Z suehring
138! var_d added, in order to normalize spectra.
139!
140! 1429 2014-07-15 12:53:45Z knoop
141! Ensemble run capability added to parallel random number generator
142!
143! 1411 2014-05-16 18:01:51Z suehring
144! Initial horizontal velocity profiles were not set to zero at the first vertical
145! grid level in case of non-cyclic lateral boundary conditions.
146!
147! 1406 2014-05-16 13:47:01Z raasch
148! bugfix: setting of initial velocities at k=1 to zero not in case of a
149! no-slip boundary condition for uv
150!
151! 1402 2014-05-09 14:25:13Z raasch
152! location messages modified
153!
154! 1400 2014-05-09 14:03:54Z knoop
155! Parallel random number generator added
156!
157! 1384 2014-05-02 14:31:06Z raasch
158! location messages added
159!
160! 1361 2014-04-16 15:17:48Z hoffmann
161! tend_* removed
162! Bugfix: w_subs is not allocated anymore if it is already allocated
163!
164! 1359 2014-04-11 17:15:14Z hoffmann
165! module lpm_init_mod added to use statements, because lpm_init has become a
166! module
167!
168! 1353 2014-04-08 15:21:23Z heinze
169! REAL constants provided with KIND-attribute
170!
171! 1340 2014-03-25 19:45:13Z kanani
172! REAL constants defined as wp-kind
173!
174! 1322 2014-03-20 16:38:49Z raasch
175! REAL constants defined as wp-kind
176! module interfaces removed
177!
178! 1320 2014-03-20 08:40:49Z raasch
179! ONLY-attribute added to USE-statements,
180! kind-parameters added to all INTEGER and REAL declaration statements,
181! kinds are defined in new module kinds,
182! revision history before 2012 removed,
183! comment fields (!:) to be used for variable explanations added to
184! all variable declaration statements
185!
186! 1316 2014-03-17 07:44:59Z heinze
187! Bugfix: allocation of w_subs
188!
189! 1299 2014-03-06 13:15:21Z heinze
190! Allocate w_subs due to extension of large scale subsidence in combination
191! with large scale forcing data (LSF_DATA)
192!
193! 1241 2013-10-30 11:36:58Z heinze
194! Overwrite initial profiles in case of nudging
195! Inititialize shf and qsws in case of large_scale_forcing
196!
197! 1221 2013-09-10 08:59:13Z raasch
198! +rflags_s_inner in copyin statement, use copyin for most arrays instead of
199! copy
200!
201! 1212 2013-08-15 08:46:27Z raasch
202! array tri is allocated and included in data copy statement
203!
204! 1195 2013-07-01 12:27:57Z heinze
205! Bugfix: move allocation of ref_state to parin.f90 and read_var_list.f90
206!
207! 1179 2013-06-14 05:57:58Z raasch
208! allocate and set ref_state to be used in buoyancy terms
209!
210! 1171 2013-05-30 11:27:45Z raasch
211! diss array is allocated with full size if accelerator boards are used
212!
213! 1159 2013-05-21 11:58:22Z fricke
214! -bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir
215!
216! 1153 2013-05-10 14:33:08Z raasch
217! diss array is allocated with dummy elements even if it is not needed
218! (required by PGI 13.4 / CUDA 5.0)
219!
220! 1115 2013-03-26 18:16:16Z hoffmann
221! unused variables removed
222!
223! 1113 2013-03-10 02:48:14Z raasch
224! openACC directive modified
225!
226! 1111 2013-03-08 23:54:10Z raasch
227! openACC directives added for pres
228! array diss allocated only if required
229!
230! 1092 2013-02-02 11:24:22Z raasch
231! unused variables removed
232!
233! 1065 2012-11-22 17:42:36Z hoffmann
234! allocation of diss (dissipation rate) in case of turbulence = .TRUE. added
235!
236! 1053 2012-11-13 17:11:03Z hoffmann
237! allocation and initialisation of necessary data arrays for the two-moment
238! cloud physics scheme the two new prognostic equations (nr, qr):
239! +dr, lambda_r, mu_r, sed_*, xr, *s, *sws, *swst, *, *_p, t*_m, *_1, *_2, *_3,
240! +tend_*, prr
241!
242! 1036 2012-10-22 13:43:42Z raasch
243! code put under GPL (PALM 3.9)
244!
245! 1032 2012-10-21 13:03:21Z letzel
246! save memory by not allocating pt_2 in case of neutral = .T.
247!
248! 1025 2012-10-07 16:04:41Z letzel
249! bugfix: swap indices of mask for ghost boundaries
250!
251! 1015 2012-09-27 09:23:24Z raasch
252! mask is set to zero for ghost boundaries
253!
254! 1010 2012-09-20 07:59:54Z raasch
255! cpp switch __nopointer added for pointer free version
256!
257! 1003 2012-09-14 14:35:53Z raasch
258! nxra,nyna, nzta replaced ny nxr, nyn, nzt
259!
260! 1001 2012-09-13 14:08:46Z raasch
261! all actions concerning leapfrog scheme removed
262!
263! 996 2012-09-07 10:41:47Z raasch
264! little reformatting
265!
266! 978 2012-08-09 08:28:32Z fricke
267! outflow damping layer removed
268! roughness length for scalar quantites z0h added
269! damping zone for the potential temperatur in case of non-cyclic lateral
270! boundaries added
271! initialization of ptdf_x, ptdf_y
272! initialization of c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l
273!
274! 849 2012-03-15 10:35:09Z raasch
275! init_particles renamed lpm_init
276!
277! 825 2012-02-19 03:03:44Z raasch
278! wang_collision_kernel renamed wang_kernel
279!
280! Revision 1.1  1998/03/09 16:22:22  raasch
281! Initial revision
282!
283!
284! Description:
285! ------------
286!> Allocation of arrays and initialization of the 3D model via
287!> a) pre-run the 1D model
288!> or
289!> b) pre-set constant linear profiles
290!> or
291!> c) read values of a previous run
292!------------------------------------------------------------------------------!
293 SUBROUTINE init_3d_model
294 
295
296    USE advec_ws
297
298    USE arrays_3d
299
300    USE constants,                                                             &
301        ONLY:  pi
302   
303    USE control_parameters
304   
305    USE flight_mod,                                                            &
306        ONLY:  flight_init
307   
308    USE grid_variables,                                                        &
309        ONLY:  dx, dy
310   
311    USE indices
312
313    USE lpm_init_mod,                                                          &
314        ONLY:  lpm_init
315   
316    USE kinds
317
318    USE land_surface_model_mod,                                                &
319        ONLY:  lsm_init, lsm_init_arrays, land_surface
320 
321    USE ls_forcing_mod
322
323    USE microphysics_mod,                                                      &
324        ONLY:  collision_turbulence, microphysics_init
325
326    USE model_1d,                                                              &
327        ONLY:  e1d, kh1d, km1d, l1d, rif1d, u1d, us1d, usws1d, v1d, vsws1d 
328   
329    USE netcdf_interface,                                                      &
330        ONLY:  dots_max, dots_num
331   
332    USE particle_attributes,                                                   &
333        ONLY:  particle_advection, use_sgs_for_particles, wang_kernel
334   
335    USE pegrid
336   
337    USE plant_canopy_model_mod,                                                &
338        ONLY:  pcm_init, plant_canopy
339
340    USE radiation_model_mod,                                                   &
341        ONLY:  radiation_init, radiation
342   
343    USE random_function_mod 
344   
345    USE random_generator_parallel,                                             &
346        ONLY:  random_number_parallel, random_seed_parallel, random_dummy,     &
347               id_random_array, seq_random_array
348   
349    USE statistics,                                                            &
350        ONLY:  hom, hom_sum, mean_surface_level_height, pr_palm, rmask,        &
351               statistic_regions, sums, sums_divnew_l, sums_divold_l, sums_l,  &
352               sums_l_l, sums_up_fraction_l, sums_wsts_bc_l, ts_value,         &
353               weight_pres, weight_substep
354 
355    USE surface_layer_fluxes_mod,                                              &
356        ONLY:  init_surface_layer_fluxes
357   
358    USE transpose_indices 
359
360    USE wind_turbine_model_mod,                                                &
361        ONLY:  wtm_init, wtm_init_arrays, wind_turbine
362
363    IMPLICIT NONE
364
365    INTEGER(iwp) ::  i             !<
366    INTEGER(iwp) ::  ind_array(1)  !<
367    INTEGER(iwp) ::  j             !<
368    INTEGER(iwp) ::  k             !<
369    INTEGER(iwp) ::  sr            !<
370
371    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  ngp_2dh_l  !<
372
373    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l    !<
374    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner_l  !<
375
376    REAL(wp), DIMENSION(1:3) ::  volume_flow_area_l     !<
377    REAL(wp), DIMENSION(1:3) ::  volume_flow_initial_l  !<
378
379    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height_l    !<
380    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l    !<
381    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_tmp  !<
382
383
384    CALL location_message( 'allocating arrays', .FALSE. )
385!
386!-- Allocate arrays
387    ALLOCATE( mean_surface_level_height(0:statistic_regions),                  &
388              mean_surface_level_height_l(0:statistic_regions),                &
389              ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions),    &
390              ngp_3d(0:statistic_regions),                                     &
391              ngp_3d_inner(0:statistic_regions),                               &
392              ngp_3d_inner_l(0:statistic_regions),                             &
393              ngp_3d_inner_tmp(0:statistic_regions),                           &
394              sums_divnew_l(0:statistic_regions),                              &
395              sums_divold_l(0:statistic_regions) )
396    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
397    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                    &
398              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),                  &
399              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),                  &
400              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),                &
401              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),                  &
402              sums(nzb:nzt+1,pr_palm+max_pr_user),                             &
403              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),      &
404              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1),    &
405              sums_up_fraction_l(10,3,0:statistic_regions),                    &
406              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                   &
407              ts_value(dots_max,0:statistic_regions) )
408    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
409
410    ALLOCATE( ol(nysg:nyng,nxlg:nxrg), shf(nysg:nyng,nxlg:nxrg),               &
411              ts(nysg:nyng,nxlg:nxrg), tswst(nysg:nyng,nxlg:nxrg),             &
412              us(nysg:nyng,nxlg:nxrg), usws(nysg:nyng,nxlg:nxrg),              &
413              uswst(nysg:nyng,nxlg:nxrg), vsws(nysg:nyng,nxlg:nxrg),           &
414              vswst(nysg:nyng,nxlg:nxrg), z0(nysg:nyng,nxlg:nxrg),             &
415              z0h(nysg:nyng,nxlg:nxrg), z0q(nysg:nyng,nxlg:nxrg) )
416
417    ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr),                                    &
418              kh(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
419              km(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
420              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
421              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
422
423#if defined( __nopointer )
424    ALLOCATE( e(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
425              e_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
426              pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
427              pt_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
428              u(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
429              u_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
430              v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
431              v_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
432              w(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
433              w_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
434              te_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
435              tpt_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
436              tu_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
437              tv_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
438              tw_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
439#else
440    ALLOCATE( e_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
441              e_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
442              e_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
443              pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
444              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
445              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
446              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
447              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
448              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
449              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
450              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
451              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
452              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
453              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
454    IF (  .NOT.  neutral )  THEN
455       ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
456    ENDIF
457#endif
458
459!
460!-- Following array is required for perturbation pressure within the iterative
461!-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds
462!-- the weighted average of the substeps and cannot be used in the Poisson
463!-- solver.
464    IF ( psolver == 'sor' )  THEN
465       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
466    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
467!
468!--    For performance reasons, multigrid is using one ghost layer only
469       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
470    ENDIF
471
472!
473!-- Array for storing constant coeffficients of the tridiagonal solver
474    IF ( psolver == 'poisfft' )  THEN
475       ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) )
476       ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) )
477    ENDIF
478
479    IF ( humidity )  THEN
480!
481!--    2D-humidity
482       ALLOCATE ( qs(nysg:nyng,nxlg:nxrg),                                     &
483                  qsws(nysg:nyng,nxlg:nxrg),                                   &
484                  qswst(nysg:nyng,nxlg:nxrg) )
485
486!
487!--    3D-humidity
488#if defined( __nopointer )
489       ALLOCATE( q(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
490                 q_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
491                 tq_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
492#else
493       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
494                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
495                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
496#endif
497
498!
499!--    3D-arrays needed for humidity
500       IF ( humidity )  THEN
501#if defined( __nopointer )
502          ALLOCATE( vpt(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
503#else
504          ALLOCATE( vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
505#endif
506
507          IF ( cloud_physics )  THEN
508
509!
510!--          Liquid water content
511#if defined( __nopointer )
512             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
513#else
514             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
515#endif
516!
517!--          Precipitation amount and rate (only needed if output is switched)
518             ALLOCATE( precipitation_amount(nysg:nyng,nxlg:nxrg),              &
519                       precipitation_rate(nysg:nyng,nxlg:nxrg) )
520
521!
522!--          3D-cloud water content
523#if defined( __nopointer )
524             ALLOCATE( qc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
525#else
526             ALLOCATE( qc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
527#endif
528!
529!--          3d-precipitation rate
530             ALLOCATE( prr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
531
532             IF ( microphysics_seifert )  THEN
533!
534!--             2D-rain water content and rain drop concentration arrays
535                ALLOCATE ( qrs(nysg:nyng,nxlg:nxrg),                        &
536                           qrsws(nysg:nyng,nxlg:nxrg),                      &
537                           qrswst(nysg:nyng,nxlg:nxrg),                     &
538                           nrs(nysg:nyng,nxlg:nxrg),                        &
539                           nrsws(nysg:nyng,nxlg:nxrg),                      &
540                           nrswst(nysg:nyng,nxlg:nxrg) )
541!
542!--             3D-rain water content, rain drop concentration arrays
543#if defined( __nopointer )
544                ALLOCATE( nr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
545                          nr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
546                          qr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
547                          qr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
548                          tnr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),             &
549                          tqr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
550#else
551                ALLOCATE( nr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
552                          nr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
553                          nr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
554                          qr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
555                          qr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
556                          qr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
557#endif
558             ENDIF
559
560          ENDIF
561
562          IF ( cloud_droplets )  THEN
563!
564!--          Liquid water content, change in liquid water content
565#if defined( __nopointer )
566             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
567                        ql_c(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
568#else
569             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
570                        ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
571#endif
572!
573!--          Real volume of particles (with weighting), volume of particles
574             ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
575                        ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
576          ENDIF
577
578       ENDIF
579
580    ENDIF
581   
582   
583    IF ( passive_scalar )  THEN
584!
585!--    2D-scalar arrays
586       ALLOCATE ( ss(nysg:nyng,nxlg:nxrg),                                     &
587                  ssws(nysg:nyng,nxlg:nxrg),                                   &
588                  sswst(nysg:nyng,nxlg:nxrg) )
589
590!
591!--    3D scalar arrays
592#if defined( __nopointer )
593       ALLOCATE( s(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
594                 s_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
595                 ts_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
596#else
597       ALLOCATE( s_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
598                 s_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
599                 s_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
600#endif
601    ENDIF
602
603    IF ( ocean )  THEN
604       ALLOCATE( saswsb(nysg:nyng,nxlg:nxrg),                                  &
605                 saswst(nysg:nyng,nxlg:nxrg) )
606#if defined( __nopointer )
607       ALLOCATE( prho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
608                 rho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
609                 sa(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
610                 sa_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
611                 tsa_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
612#else
613       ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
614                 rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                         &
615                 sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
616                 sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
617                 sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
618       prho => prho_1
619       rho  => rho_1  ! routines calc_mean_profile and diffusion_e require
620                      ! density to be apointer
621#endif
622       IF ( humidity_remote )  THEN
623          ALLOCATE( qswst_remote(nysg:nyng,nxlg:nxrg))
624          qswst_remote = 0.0_wp
625       ENDIF
626    ENDIF
627
628!
629!-- 3D-array for storing the dissipation, needed for calculating the sgs
630!-- particle velocities
631    IF ( use_sgs_for_particles  .OR.  wang_kernel  .OR.  collision_turbulence  &
632         .OR.  num_acc_per_node > 0 )  THEN
633       ALLOCATE( diss(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
634    ENDIF
635
636!
637!-- 1D-array for large scale subsidence velocity
638    IF ( .NOT. ALLOCATED( w_subs ) )  THEN
639       ALLOCATE ( w_subs(nzb:nzt+1) )
640       w_subs = 0.0_wp
641    ENDIF
642
643!
644!-- ID-array and state-space-array for the parallel random number generator
645    IF ( random_generator == 'random-parallel' )  THEN
646       ALLOCATE ( seq_random_array(5,nysg:nyng,nxlg:nxrg) )
647       ALLOCATE ( id_random_array(0:ny,0:nx) )
648       seq_random_array = 0
649       id_random_array  = 0
650    ENDIF
651   
652!
653!-- 4D-array for storing the Rif-values at vertical walls
654    IF ( topography /= 'flat' )  THEN
655       ALLOCATE( rif_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg,1:4) )
656       rif_wall = 0.0_wp
657    ENDIF
658
659!
660!-- Arrays to store velocity data from t-dt and the phase speeds which
661!-- are needed for radiation boundary conditions
662    IF ( outflow_l )  THEN
663       ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2),                               &
664                 v_m_l(nzb:nzt+1,nysg:nyng,0:1),                               &
665                 w_m_l(nzb:nzt+1,nysg:nyng,0:1) )
666    ENDIF
667    IF ( outflow_r )  THEN
668       ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
669                 v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
670                 w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) )
671    ENDIF
672    IF ( outflow_l  .OR.  outflow_r )  THEN
673       ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng),           &
674                 c_w(nzb:nzt+1,nysg:nyng) )
675    ENDIF
676    IF ( outflow_s )  THEN
677       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg),                               &
678                 v_m_s(nzb:nzt+1,1:2,nxlg:nxrg),                               &
679                 w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) )
680    ENDIF
681    IF ( outflow_n )  THEN
682       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
683                 v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
684                 w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) )
685    ENDIF
686    IF ( outflow_s  .OR.  outflow_n )  THEN
687       ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg),           &
688                 c_w(nzb:nzt+1,nxlg:nxrg) )
689    ENDIF
690    IF ( outflow_l  .OR.  outflow_r  .OR.  outflow_s  .OR.  outflow_n )  THEN
691       ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) )                   
692       ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) )
693    ENDIF
694
695
696#if ! defined( __nopointer )
697!
698!-- Initial assignment of the pointers
699    e  => e_1;   e_p  => e_2;   te_m  => e_3
700    IF ( .NOT. neutral )  THEN
701       pt => pt_1;  pt_p => pt_2;  tpt_m => pt_3
702    ELSE
703       pt => pt_1;  pt_p => pt_1;  tpt_m => pt_3
704    ENDIF
705    u  => u_1;   u_p  => u_2;   tu_m  => u_3
706    v  => v_1;   v_p  => v_2;   tv_m  => v_3
707    w  => w_1;   w_p  => w_2;   tw_m  => w_3
708
709    IF ( humidity )  THEN
710       q => q_1;  q_p => q_2;  tq_m => q_3
711       IF ( humidity )  THEN
712          vpt  => vpt_1   
713          IF ( cloud_physics )  THEN
714             ql => ql_1
715             qc => qc_1
716             IF ( microphysics_seifert )  THEN
717                qr => qr_1;  qr_p  => qr_2;  tqr_m  => qr_3
718                nr => nr_1;  nr_p  => nr_2;  tnr_m  => nr_3
719             ENDIF
720          ENDIF
721       ENDIF
722       IF ( cloud_droplets )  THEN
723          ql   => ql_1
724          ql_c => ql_2
725       ENDIF
726    ENDIF
727   
728    IF ( passive_scalar )  THEN
729       s => s_1;  s_p => s_2;  ts_m => s_3
730    ENDIF   
731
732    IF ( ocean )  THEN
733       sa => sa_1;  sa_p => sa_2;  tsa_m => sa_3
734    ENDIF
735#endif
736
737!
738!-- Allocate land surface model arrays
739    IF ( land_surface )  THEN
740       CALL lsm_init_arrays
741    ENDIF
742
743!
744!-- Allocate wind turbine model arrays
745    IF ( wind_turbine )  THEN
746       CALL wtm_init_arrays
747    ENDIF
748   
749!
750!-- Initialize virtual flight measurements
751    IF ( virtual_flight )  THEN
752       CALL flight_init
753    ENDIF
754
755!
756!-- Allocate arrays containing the RK coefficient for calculation of
757!-- perturbation pressure and turbulent fluxes. At this point values are
758!-- set for pressure calculation during initialization (where no timestep
759!-- is done). Further below the values needed within the timestep scheme
760!-- will be set.
761    ALLOCATE( weight_substep(1:intermediate_timestep_count_max),               &
762              weight_pres(1:intermediate_timestep_count_max) )
763    weight_substep = 1.0_wp
764    weight_pres    = 1.0_wp
765    intermediate_timestep_count = 0  ! needed when simulated_time = 0.0
766       
767    CALL location_message( 'finished', .TRUE. )
768
769!
770!-- Initialize local summation arrays for routine flow_statistics.
771!-- This is necessary because they may not yet have been initialized when they
772!-- are called from flow_statistics (or - depending on the chosen model run -
773!-- are never initialized)
774    sums_divnew_l      = 0.0_wp
775    sums_divold_l      = 0.0_wp
776    sums_l_l           = 0.0_wp
777    sums_up_fraction_l = 0.0_wp
778    sums_wsts_bc_l     = 0.0_wp
779
780
781!
782!-- Initialize model variables
783    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
784         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
785!
786!--    First model run of a possible job queue.
787!--    Initial profiles of the variables must be computes.
788       IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
789
790          CALL location_message( 'initializing with 1D model profiles', .FALSE. )
791!
792!--       Use solutions of the 1D model as initial profiles,
793!--       start 1D model
794          CALL init_1d_model
795!
796!--       Transfer initial profiles to the arrays of the 3D model
797          DO  i = nxlg, nxrg
798             DO  j = nysg, nyng
799                e(:,j,i)  = e1d
800                kh(:,j,i) = kh1d
801                km(:,j,i) = km1d
802                pt(:,j,i) = pt_init
803                u(:,j,i)  = u1d
804                v(:,j,i)  = v1d
805             ENDDO
806          ENDDO
807
808          IF ( humidity )  THEN
809             DO  i = nxlg, nxrg
810                DO  j = nysg, nyng
811                   q(:,j,i) = q_init
812                ENDDO
813             ENDDO
814             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
815                DO  i = nxlg, nxrg
816                   DO  j = nysg, nyng
817                      qr(:,j,i) = 0.0_wp
818                      nr(:,j,i) = 0.0_wp
819                   ENDDO
820                ENDDO
821
822             ENDIF
823          ENDIF
824          IF ( passive_scalar )  THEN
825             DO  i = nxlg, nxrg
826                DO  j = nysg, nyng
827                   s(:,j,i) = s_init
828                ENDDO
829             ENDDO   
830          ENDIF
831
832          IF ( .NOT. constant_diffusion )  THEN
833             DO  i = nxlg, nxrg
834                DO  j = nysg, nyng
835                   e(:,j,i)  = e1d
836                ENDDO
837             ENDDO
838!
839!--          Store initial profiles for output purposes etc.
840             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
841
842             IF ( constant_flux_layer )  THEN
843                ol   = ( zu(nzb+1) - zw(nzb) ) / ( rif1d(nzb+1) + 1.0E-20_wp )
844                ts   = 0.0_wp  ! could actually be computed more accurately in the
845                               ! 1D model. Update when opportunity arises.
846                us   = us1d
847                usws = usws1d
848                vsws = vsws1d
849             ELSE
850                ts   = 0.0_wp  ! must be set, because used in
851                ol   = ( zu(nzb+1) - zw(nzb) ) / zeta_min  ! flowste
852                us   = 0.0_wp
853                usws = 0.0_wp
854                vsws = 0.0_wp
855             ENDIF
856
857          ELSE
858             e    = 0.0_wp  ! must be set, because used in
859             ol   = ( zu(nzb+1) - zw(nzb) ) / zeta_min  ! flowste
860             ts   = 0.0_wp
861             us   = 0.0_wp
862             usws = 0.0_wp
863             vsws = 0.0_wp
864          ENDIF
865          uswst = top_momentumflux_u
866          vswst = top_momentumflux_v
867
868!
869!--       In every case qs = 0.0 (see also pt)
870!--       This could actually be computed more accurately in the 1D model.
871!--       Update when opportunity arises!
872          IF ( humidity )  THEN
873             qs = 0.0_wp
874             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
875                qrs = 0.0_wp
876                nrs = 0.0_wp
877             ENDIF
878          ENDIF
879!
880!--       Initialize scaling parameter for passive scalar
881          IF ( passive_scalar ) ss = 0.0_wp         
882
883!
884!--       Inside buildings set velocities back to zero
885          IF ( topography /= 'flat' )  THEN
886             DO  i = nxlg, nxrg
887                DO  j = nysg, nyng
888                   u(nzb:nzb_u_inner(j,i),j,i) = 0.0_wp
889                   v(nzb:nzb_v_inner(j,i),j,i) = 0.0_wp
890                ENDDO
891             ENDDO
892             
893!
894!--          WARNING: The extra boundary conditions set after running the
895!--          -------  1D model impose an error on the divergence one layer
896!--                   below the topography; need to correct later
897!--          ATTENTION: Provisional correction for Piacsek & Williams
898!--          ---------  advection scheme: keep u and v zero one layer below
899!--                     the topography.
900             IF ( ibc_uv_b == 1 )  THEN
901!
902!--             Neumann condition
903                DO  i = nxl-1, nxr+1
904                   DO  j = nys-1, nyn+1
905                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i)
906                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i)
907                   ENDDO
908                ENDDO
909
910             ENDIF
911
912          ENDIF
913
914          CALL location_message( 'finished', .TRUE. )
915
916       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )    &
917       THEN
918
919          CALL location_message( 'initializing with constant profiles', .FALSE. )
920!
921!--       Overwrite initial profiles in case of nudging
922          IF ( nudging )  THEN
923             pt_init = ptnudge(:,1)
924             u_init  = unudge(:,1)
925             v_init  = vnudge(:,1)
926             IF ( humidity  )  THEN ! is passive_scalar correct???
927                q_init = qnudge(:,1)
928             ENDIF
929
930             WRITE( message_string, * ) 'Initial profiles of u, v and ',       &
931                 'scalars from NUDGING_DATA are used.'
932             CALL message( 'init_3d_model', 'PA0370', 0, 0, 0, 6, 0 )
933          ENDIF
934
935!
936!--       Use constructed initial profiles (velocity constant with height,
937!--       temperature profile with constant gradient)
938          DO  i = nxlg, nxrg
939             DO  j = nysg, nyng
940                pt(:,j,i) = pt_init
941                u(:,j,i)  = u_init
942                v(:,j,i)  = v_init
943             ENDDO
944          ENDDO
945
946!
947!--       Set initial horizontal velocities at the lowest computational grid
948!--       levels to zero in order to avoid too small time steps caused by the
949!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
950!--       in the limiting formula!).
951          IF ( ibc_uv_b /= 1 )  THEN
952             DO  i = nxlg, nxrg
953                DO  j = nysg, nyng
954                   u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0_wp
955                   v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0_wp
956                ENDDO
957             ENDDO
958          ENDIF
959
960          IF ( humidity )  THEN
961             DO  i = nxlg, nxrg
962                DO  j = nysg, nyng
963                   q(:,j,i) = q_init
964                ENDDO
965             ENDDO
966             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
967
968                DO  i = nxlg, nxrg
969                   DO  j = nysg, nyng
970                      qr(:,j,i) = 0.0_wp
971                      nr(:,j,i) = 0.0_wp
972                   ENDDO
973                ENDDO
974
975             ENDIF
976          ENDIF
977         
978          IF ( passive_scalar )  THEN
979             DO  i = nxlg, nxrg
980                DO  j = nysg, nyng
981                   s(:,j,i) = s_init
982                ENDDO
983             ENDDO
984          ENDIF
985
986          IF ( ocean )  THEN
987             DO  i = nxlg, nxrg
988                DO  j = nysg, nyng
989                   sa(:,j,i) = sa_init
990                ENDDO
991             ENDDO
992          ENDIF
993         
994          IF ( constant_diffusion )  THEN
995             km   = km_constant
996             kh   = km / prandtl_number
997             e    = 0.0_wp
998          ELSEIF ( e_init > 0.0_wp )  THEN
999             DO  k = nzb+1, nzt
1000                km(k,:,:) = 0.1_wp * l_grid(k) * SQRT( e_init )
1001             ENDDO
1002             km(nzb,:,:)   = km(nzb+1,:,:)
1003             km(nzt+1,:,:) = km(nzt,:,:)
1004             kh   = km / prandtl_number
1005             e    = e_init
1006          ELSE
1007             IF ( .NOT. ocean )  THEN
1008                kh   = 0.01_wp   ! there must exist an initial diffusion, because
1009                km   = 0.01_wp   ! otherwise no TKE would be produced by the
1010                              ! production terms, as long as not yet
1011                              ! e = (u*/cm)**2 at k=nzb+1
1012             ELSE
1013                kh   = 0.00001_wp
1014                km   = 0.00001_wp
1015             ENDIF
1016             e    = 0.0_wp
1017          ENDIF
1018          ol    = ( zu(nzb+1) - zw(nzb) ) / zeta_min
1019          ts    = 0.0_wp
1020!
1021!--       Very small number is required for calculation of Obukhov length
1022!--       at first timestep     
1023          us    = 1E-30_wp 
1024          usws  = 0.0_wp
1025          uswst = top_momentumflux_u
1026          vsws  = 0.0_wp
1027          vswst = top_momentumflux_v
1028          IF ( humidity )       qs = 0.0_wp
1029          IF ( passive_scalar ) ss = 0.0_wp
1030
1031!
1032!--       Compute initial temperature field and other constants used in case
1033!--       of a sloping surface
1034          IF ( sloping_surface )  CALL init_slope
1035
1036          CALL location_message( 'finished', .TRUE. )
1037
1038       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )                  &
1039       THEN
1040
1041          CALL location_message( 'initializing by user', .FALSE. )
1042!
1043!--       Initialization will completely be done by the user
1044          CALL user_init_3d_model
1045
1046          CALL location_message( 'finished', .TRUE. )
1047
1048       ENDIF
1049
1050       CALL location_message( 'initializing statistics, boundary conditions, etc.', &
1051                              .FALSE. )
1052
1053!
1054!--    Bottom boundary
1055       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
1056          u(nzb,:,:) = 0.0_wp
1057          v(nzb,:,:) = 0.0_wp
1058       ENDIF
1059
1060!
1061!--    Apply channel flow boundary condition
1062       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
1063          u(nzt+1,:,:) = 0.0_wp
1064          v(nzt+1,:,:) = 0.0_wp
1065       ENDIF
1066
1067!
1068!--    Calculate virtual potential temperature
1069       IF ( humidity )  vpt = pt * ( 1.0_wp + 0.61_wp * q )
1070
1071!
1072!--    Store initial profiles for output purposes etc.
1073       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
1074       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
1075       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
1076          hom(nzb,1,5,:) = 0.0_wp
1077          hom(nzb,1,6,:) = 0.0_wp
1078       ENDIF
1079       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1080       hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 )
1081       hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 )
1082
1083       IF ( ocean )  THEN
1084!
1085!--       Store initial salinity profile
1086          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
1087       ENDIF
1088
1089       IF ( humidity )  THEN
1090!
1091!--       Store initial profile of total water content, virtual potential
1092!--       temperature
1093          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
1094          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
1095          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
1096!
1097!--          Store initial profile of specific humidity and potential
1098!--          temperature
1099             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1100             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1101          ENDIF
1102       ENDIF
1103
1104       IF ( passive_scalar )  THEN
1105!
1106!--       Store initial scalar profile
1107          hom(:,1,115,:) = SPREAD(  s(:,nys,nxl), 2, statistic_regions+1 )
1108       ENDIF
1109
1110!
1111!--    Initialize the random number generators (from numerical recipes)
1112       CALL random_function_ini
1113       
1114       IF ( random_generator == 'random-parallel' )  THEN
1115!--       Asigning an ID to every vertical gridpoint column
1116!--       dependig on the ensemble run number.
1117          random_dummy=1
1118          DO j=0,ny
1119             DO i=0,nx
1120                id_random_array(j,i) = random_dummy + 1E6                      &
1121                                       * ( ensemble_member_nr - 1000 )
1122                random_dummy = random_dummy + 1
1123             END DO
1124          ENDDO
1125!--       Initializing with random_seed_parallel for every vertical
1126!--       gridpoint column.
1127          random_dummy=0
1128          DO j = nysg, nyng
1129             DO i = nxlg, nxrg
1130                CALL random_seed_parallel (random_sequence=id_random_array(j, i))
1131                CALL random_number_parallel (random_dummy)
1132                CALL random_seed_parallel (get=seq_random_array(:, j, i))
1133             END DO
1134          ENDDO
1135       ENDIF
1136
1137!
1138!--    Initialize fluxes at bottom surface
1139       IF ( use_surface_fluxes )  THEN
1140
1141          IF ( constant_heatflux )  THEN
1142!
1143!--          Heat flux is prescribed
1144             IF ( random_heatflux )  THEN
1145                CALL disturb_heatflux
1146             ELSE
1147                shf = surface_heatflux
1148!
1149!--             Initialize shf with data from external file LSF_DATA
1150                IF ( large_scale_forcing .AND. lsf_surf )  THEN
1151                   CALL ls_forcing_surf ( simulated_time )
1152                ENDIF
1153
1154!
1155!--             Over topography surface_heatflux is replaced by wall_heatflux(0)
1156                IF ( TRIM( topography ) /= 'flat' )  THEN
1157                   DO  i = nxlg, nxrg
1158                      DO  j = nysg, nyng
1159                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
1160                            shf(j,i) = wall_heatflux(0)
1161                         ENDIF
1162                      ENDDO
1163                   ENDDO
1164                ENDIF
1165             ENDIF
1166          ENDIF
1167
1168!
1169!--       Determine the near-surface water flux
1170          IF ( humidity )  THEN
1171             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1172                qrsws = 0.0_wp
1173                nrsws = 0.0_wp
1174             ENDIF
1175             IF ( constant_waterflux )  THEN
1176                qsws   = surface_waterflux
1177!
1178!--             Over topography surface_waterflux is replaced by
1179!--             wall_humidityflux(0)
1180                IF ( TRIM( topography ) /= 'flat' )  THEN
1181                   wall_qflux = wall_humidityflux
1182                   DO  i = nxlg, nxrg
1183                      DO  j = nysg, nyng
1184                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
1185                            qsws(j,i) = wall_qflux(0)
1186                         ENDIF
1187                      ENDDO
1188                   ENDDO
1189                ENDIF
1190             ENDIF
1191          ENDIF
1192!
1193!--       Initialize the near-surface scalar flux
1194          IF ( passive_scalar )  THEN
1195             IF ( constant_scalarflux )  THEN
1196                ssws   = surface_scalarflux
1197!
1198!--             Over topography surface_scalarflux is replaced by
1199!--             wall_scalarflux(0)
1200                IF ( TRIM( topography ) /= 'flat' )  THEN
1201                   wall_sflux = wall_scalarflux
1202                   DO  i = nxlg, nxrg
1203                      DO  j = nysg, nyng
1204                         IF ( nzb_s_inner(j,i) /= 0 )  ssws(j,i) = wall_sflux(0)
1205                      ENDDO
1206                   ENDDO
1207                ENDIF
1208             ENDIF
1209          ENDIF   
1210!
1211!--       Initialize near-surface salinity flux
1212          IF ( ocean )  saswsb = bottom_salinityflux
1213
1214       ENDIF
1215
1216!
1217!--    Initialize fluxes at top surface
1218!--    Currently, only the heatflux and salinity flux can be prescribed.
1219!--    The latent flux is zero in this case!
1220       IF ( use_top_fluxes )  THEN
1221!
1222!--       Prescribe to heat flux
1223          IF ( constant_top_heatflux )  tswst = top_heatflux
1224!
1225!--       Prescribe zero latent flux at the top     
1226          IF ( humidity )  THEN
1227             qswst = 0.0_wp
1228             IF ( cloud_physics  .AND.  microphysics_seifert ) THEN
1229                nrswst = 0.0_wp
1230                qrswst = 0.0_wp
1231             ENDIF
1232          ENDIF
1233!
1234!--       Prescribe top scalar flux
1235          IF ( passive_scalar .AND. constant_top_scalarflux )                  &
1236             sswst = top_scalarflux
1237!
1238!--       Prescribe top salinity flux
1239          IF ( ocean .AND. constant_top_salinityflux)                          &
1240             saswst = top_salinityflux
1241!
1242!--       Initialization in case of a coupled model run
1243          IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1244             tswst = 0.0_wp
1245          ENDIF
1246
1247       ENDIF
1248
1249!
1250!--    Initialize Prandtl layer quantities
1251       IF ( constant_flux_layer )  THEN
1252
1253          z0 = roughness_length
1254          z0h = z0h_factor * z0
1255          z0q = z0h_factor * z0
1256
1257          IF ( .NOT. constant_heatflux )  THEN 
1258!
1259!--          Surface temperature is prescribed. Here the heat flux cannot be
1260!--          simply estimated, because therefore ol, u* and theta* would have
1261!--          to be computed by iteration. This is why the heat flux is assumed
1262!--          to be zero before the first time step. It approaches its correct
1263!--          value in the course of the first few time steps.
1264             shf   = 0.0_wp
1265          ENDIF
1266
1267          IF ( humidity  )  THEN
1268             IF (  .NOT.  constant_waterflux )  qsws   = 0.0_wp
1269             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1270                qrsws = 0.0_wp
1271                nrsws = 0.0_wp
1272             ENDIF
1273          ENDIF
1274          IF ( passive_scalar  .AND.  .NOT.  constant_scalarflux )  ssws = 0.0_wp
1275
1276       ENDIF
1277
1278!
1279!--    Set the reference state to be used in the buoyancy terms (for ocean runs
1280!--    the reference state will be set (overwritten) in init_ocean)
1281       IF ( use_single_reference_value )  THEN
1282          IF (  .NOT.  humidity )  THEN
1283             ref_state(:) = pt_reference
1284          ELSE
1285             ref_state(:) = vpt_reference
1286          ENDIF
1287       ELSE
1288          IF (  .NOT.  humidity )  THEN
1289             ref_state(:) = pt_init(:)
1290          ELSE
1291             ref_state(:) = vpt(:,nys,nxl)
1292          ENDIF
1293       ENDIF
1294
1295!
1296!--    For the moment, vertical velocity is zero
1297       w = 0.0_wp
1298
1299!
1300!--    Initialize array sums (must be defined in first call of pres)
1301       sums = 0.0_wp
1302
1303!
1304!--    In case of iterative solvers, p must get an initial value
1305       IF ( psolver(1:9) == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0_wp
1306
1307!
1308!--    Treating cloud physics, liquid water content and precipitation amount
1309!--    are zero at beginning of the simulation
1310       IF ( cloud_physics )  THEN
1311          ql = 0.0_wp
1312          qc = 0.0_wp
1313
1314          precipitation_amount = 0.0_wp
1315       ENDIF
1316!
1317!--    Impose vortex with vertical axis on the initial velocity profile
1318       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
1319          CALL init_rankine
1320       ENDIF
1321
1322!
1323!--    Impose temperature anomaly (advection test only)
1324       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
1325          CALL init_pt_anomaly
1326       ENDIF
1327
1328!
1329!--    If required, change the surface temperature at the start of the 3D run
1330       IF ( pt_surface_initial_change /= 0.0_wp )  THEN
1331          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
1332       ENDIF
1333
1334!
1335!--    If required, change the surface humidity/scalar at the start of the 3D
1336!--    run
1337       IF ( humidity  .AND.  q_surface_initial_change /= 0.0_wp )              &
1338          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
1339         
1340       IF ( passive_scalar .AND.  s_surface_initial_change /= 0.0_wp )         &
1341          s(nzb,:,:) = s(nzb,:,:) + s_surface_initial_change
1342       
1343
1344!
1345!--    Initialize old and new time levels.
1346       te_m = 0.0_wp; tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1347       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
1348
1349       IF ( humidity  )  THEN
1350          tq_m = 0.0_wp
1351          q_p = q
1352          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1353             tqr_m = 0.0_wp
1354             qr_p  = qr
1355             tnr_m = 0.0_wp
1356             nr_p  = nr
1357          ENDIF
1358       ENDIF
1359       
1360       IF ( passive_scalar )  THEN
1361          ts_m = 0.0_wp
1362          s_p  = s
1363       ENDIF       
1364
1365       IF ( ocean )  THEN
1366          tsa_m = 0.0_wp
1367          sa_p  = sa
1368       ENDIF
1369       
1370       CALL location_message( 'finished', .TRUE. )
1371
1372    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.         &
1373         TRIM( initializing_actions ) == 'cyclic_fill' )                       &
1374    THEN
1375
1376       CALL location_message( 'initializing in case of restart / cyclic_fill', &
1377                              .FALSE. )
1378!
1379!--    When reading data for cyclic fill of 3D prerun data files, read
1380!--    some of the global variables from the restart file which are required
1381!--    for initializing the inflow
1382       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1383
1384          DO  i = 0, io_blocks-1
1385             IF ( i == io_group )  THEN
1386                CALL read_parts_of_var_list
1387                CALL close_file( 13 )
1388             ENDIF
1389#if defined( __parallel )
1390             CALL MPI_BARRIER( comm2d, ierr )
1391#endif
1392          ENDDO
1393
1394       ENDIF
1395
1396!
1397!--    Read binary data from restart file
1398       DO  i = 0, io_blocks-1
1399          IF ( i == io_group )  THEN
1400             CALL read_3d_binary
1401          ENDIF
1402#if defined( __parallel )
1403          CALL MPI_BARRIER( comm2d, ierr )
1404#endif
1405       ENDDO
1406
1407!
1408!--    Initialization of the turbulence recycling method
1409       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.               &
1410            turbulent_inflow )  THEN
1411!
1412!--       First store the profiles to be used at the inflow.
1413!--       These profiles are the (temporally) and horizontally averaged vertical
1414!--       profiles from the prerun. Alternatively, prescribed profiles
1415!--       for u,v-components can be used.
1416          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,7) )
1417
1418          IF ( use_prescribed_profile_data )  THEN
1419             mean_inflow_profiles(:,1) = u_init            ! u
1420             mean_inflow_profiles(:,2) = v_init            ! v
1421          ELSE
1422             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
1423             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
1424          ENDIF
1425          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
1426          mean_inflow_profiles(:,5) = hom_sum(:,8,0)       ! e
1427          IF ( humidity )                                                      &
1428             mean_inflow_profiles(:,6) = hom_sum(:,41,0)   ! q
1429          IF ( passive_scalar )                                                &
1430             mean_inflow_profiles(:,7) = hom_sum(:,115,0)   ! s
1431
1432!
1433!--       If necessary, adjust the horizontal flow field to the prescribed
1434!--       profiles
1435          IF ( use_prescribed_profile_data )  THEN
1436             DO  i = nxlg, nxrg
1437                DO  j = nysg, nyng
1438                   DO  k = nzb, nzt+1
1439                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1440                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
1441                   ENDDO
1442                ENDDO
1443             ENDDO
1444          ENDIF
1445
1446!
1447!--       Use these mean profiles at the inflow (provided that Dirichlet
1448!--       conditions are used)
1449          IF ( inflow_l )  THEN
1450             DO  j = nysg, nyng
1451                DO  k = nzb, nzt+1
1452                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1453                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
1454                   w(k,j,nxlg:-1)  = 0.0_wp
1455                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
1456                   e(k,j,nxlg:-1)  = mean_inflow_profiles(k,5)
1457                   IF ( humidity )                                             &
1458                      q(k,j,nxlg:-1)  = mean_inflow_profiles(k,6)
1459                   IF ( passive_scalar )                                       &
1460                      s(k,j,nxlg:-1)  = mean_inflow_profiles(k,7)                     
1461                ENDDO
1462             ENDDO
1463          ENDIF
1464
1465!
1466!--       Calculate the damping factors to be used at the inflow. For a
1467!--       turbulent inflow the turbulent fluctuations have to be limited
1468!--       vertically because otherwise the turbulent inflow layer will grow
1469!--       in time.
1470          IF ( inflow_damping_height == 9999999.9_wp )  THEN
1471!
1472!--          Default: use the inversion height calculated by the prerun; if
1473!--          this is zero, inflow_damping_height must be explicitly
1474!--          specified.
1475             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0_wp )  THEN
1476                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1477             ELSE
1478                WRITE( message_string, * ) 'inflow_damping_height must be ',   &
1479                     'explicitly specified because&the inversion height ',     &
1480                     'calculated by the prerun is zero.'
1481                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
1482             ENDIF
1483
1484          ENDIF
1485
1486          IF ( inflow_damping_width == 9999999.9_wp )  THEN
1487!
1488!--          Default for the transition range: one tenth of the undamped
1489!--          layer
1490             inflow_damping_width = 0.1_wp * inflow_damping_height
1491
1492          ENDIF
1493
1494          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
1495
1496          DO  k = nzb, nzt+1
1497
1498             IF ( zu(k) <= inflow_damping_height )  THEN
1499                inflow_damping_factor(k) = 1.0_wp
1500             ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) )  THEN
1501                inflow_damping_factor(k) = 1.0_wp -                            &
1502                                           ( zu(k) - inflow_damping_height ) / &
1503                                           inflow_damping_width
1504             ELSE
1505                inflow_damping_factor(k) = 0.0_wp
1506             ENDIF
1507
1508          ENDDO
1509
1510       ENDIF
1511
1512!
1513!--    Inside buildings set velocities and TKE back to zero
1514       IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND.                &
1515            topography /= 'flat' )  THEN
1516!
1517!--       Inside buildings set velocities and TKE back to zero.
1518!--       Other scalars (pt, q, s, km, kh, p, sa, ...) are ignored at present,
1519!--       maybe revise later.
1520          DO  i = nxlg, nxrg
1521             DO  j = nysg, nyng
1522                u  (nzb:nzb_u_inner(j,i),j,i)   = 0.0_wp
1523                v  (nzb:nzb_v_inner(j,i),j,i)   = 0.0_wp
1524                w  (nzb:nzb_w_inner(j,i),j,i)   = 0.0_wp
1525                e  (nzb:nzb_w_inner(j,i),j,i)   = 0.0_wp
1526                tu_m(nzb:nzb_u_inner(j,i),j,i)  = 0.0_wp
1527                tv_m(nzb:nzb_v_inner(j,i),j,i)  = 0.0_wp
1528                tw_m(nzb:nzb_w_inner(j,i),j,i)  = 0.0_wp
1529                te_m(nzb:nzb_w_inner(j,i),j,i)  = 0.0_wp
1530                tpt_m(nzb:nzb_w_inner(j,i),j,i) = 0.0_wp
1531             ENDDO
1532          ENDDO
1533
1534       ENDIF
1535
1536!
1537!--    Calculate initial temperature field and other constants used in case
1538!--    of a sloping surface
1539       IF ( sloping_surface )  CALL init_slope
1540
1541!
1542!--    Initialize new time levels (only done in order to set boundary values
1543!--    including ghost points)
1544       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
1545       IF ( humidity )  THEN
1546          q_p = q
1547          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1548             qr_p = qr
1549             nr_p = nr
1550          ENDIF
1551       ENDIF
1552       IF ( passive_scalar )  s_p  = s
1553       IF ( ocean          )  sa_p = sa
1554
1555!
1556!--    Allthough tendency arrays are set in prognostic_equations, they have
1557!--    have to be predefined here because they are used (but multiplied with 0)
1558!--    there before they are set.
1559       te_m = 0.0_wp; tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1560       IF ( humidity )  THEN
1561          tq_m = 0.0_wp
1562          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1563             tqr_m = 0.0_wp
1564             tnr_m = 0.0_wp
1565          ENDIF
1566       ENDIF
1567       IF ( passive_scalar )  ts_m  = 0.0_wp
1568       IF ( ocean          )  tsa_m = 0.0_wp
1569
1570       CALL location_message( 'finished', .TRUE. )
1571
1572    ELSE
1573!
1574!--    Actually this part of the programm should not be reached
1575       message_string = 'unknown initializing problem'
1576       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
1577    ENDIF
1578
1579
1580    IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1581!
1582!--    Initialize old timelevels needed for radiation boundary conditions
1583       IF ( outflow_l )  THEN
1584          u_m_l(:,:,:) = u(:,:,1:2)
1585          v_m_l(:,:,:) = v(:,:,0:1)
1586          w_m_l(:,:,:) = w(:,:,0:1)
1587       ENDIF
1588       IF ( outflow_r )  THEN
1589          u_m_r(:,:,:) = u(:,:,nx-1:nx)
1590          v_m_r(:,:,:) = v(:,:,nx-1:nx)
1591          w_m_r(:,:,:) = w(:,:,nx-1:nx)
1592       ENDIF
1593       IF ( outflow_s )  THEN
1594          u_m_s(:,:,:) = u(:,0:1,:)
1595          v_m_s(:,:,:) = v(:,1:2,:)
1596          w_m_s(:,:,:) = w(:,0:1,:)
1597       ENDIF
1598       IF ( outflow_n )  THEN
1599          u_m_n(:,:,:) = u(:,ny-1:ny,:)
1600          v_m_n(:,:,:) = v(:,ny-1:ny,:)
1601          w_m_n(:,:,:) = w(:,ny-1:ny,:)
1602       ENDIF
1603       
1604    ENDIF
1605
1606!
1607!-- Calculate the initial volume flow at the right and north boundary
1608    IF ( conserve_volume_flow )  THEN
1609
1610       IF ( use_prescribed_profile_data )  THEN
1611
1612          volume_flow_initial_l = 0.0_wp
1613          volume_flow_area_l    = 0.0_wp
1614
1615          IF ( nxr == nx )  THEN
1616             DO  j = nys, nyn
1617                DO  k = nzb_u_inner(j,nx)+1, nzt
1618                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1619                                              u_init(k) * dzw(k)
1620                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1621                ENDDO
1622             ENDDO
1623          ENDIF
1624         
1625          IF ( nyn == ny )  THEN
1626             DO  i = nxl, nxr
1627                DO  k = nzb_v_inner(ny,i)+1, nzt
1628                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1629                                              v_init(k) * dzw(k)
1630                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1631                ENDDO
1632             ENDDO
1633          ENDIF
1634
1635#if defined( __parallel )
1636          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1637                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1638          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1639                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1640
1641#else
1642          volume_flow_initial = volume_flow_initial_l
1643          volume_flow_area    = volume_flow_area_l
1644#endif 
1645
1646       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1647
1648          volume_flow_initial_l = 0.0_wp
1649          volume_flow_area_l    = 0.0_wp
1650
1651          IF ( nxr == nx )  THEN
1652             DO  j = nys, nyn
1653                DO  k = nzb_u_inner(j,nx)+1, nzt
1654                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1655                                              hom_sum(k,1,0) * dzw(k)
1656                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1657                ENDDO
1658             ENDDO
1659          ENDIF
1660         
1661          IF ( nyn == ny )  THEN
1662             DO  i = nxl, nxr
1663                DO  k = nzb_v_inner(ny,i)+1, nzt
1664                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1665                                              hom_sum(k,2,0) * dzw(k)
1666                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1667                ENDDO
1668             ENDDO
1669          ENDIF
1670
1671#if defined( __parallel )
1672          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1673                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1674          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1675                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1676
1677#else
1678          volume_flow_initial = volume_flow_initial_l
1679          volume_flow_area    = volume_flow_area_l
1680#endif 
1681
1682       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1683
1684          volume_flow_initial_l = 0.0_wp
1685          volume_flow_area_l    = 0.0_wp
1686
1687          IF ( nxr == nx )  THEN
1688             DO  j = nys, nyn
1689                DO  k = nzb_u_inner(j,nx)+1, nzt
1690                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1691                                              u(k,j,nx) * dzw(k)
1692                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1693                ENDDO
1694             ENDDO
1695          ENDIF
1696         
1697          IF ( nyn == ny )  THEN
1698             DO  i = nxl, nxr
1699                DO  k = nzb_v_inner(ny,i)+1, nzt
1700                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1701                                              v(k,ny,i) * dzw(k)
1702                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1703                ENDDO
1704             ENDDO
1705          ENDIF
1706
1707#if defined( __parallel )
1708          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1709                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1710          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1711                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1712
1713#else
1714          volume_flow_initial = volume_flow_initial_l
1715          volume_flow_area    = volume_flow_area_l
1716#endif 
1717
1718       ENDIF
1719
1720!
1721!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated
1722!--    from u|v_bulk instead
1723       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
1724          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
1725          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
1726       ENDIF
1727
1728    ENDIF
1729
1730!
1731!-- Initialize quantities for special advections schemes
1732    CALL init_advec
1733
1734!
1735!-- Impose random perturbation on the horizontal velocity field and then
1736!-- remove the divergences from the velocity field at the initial stage
1737    IF ( create_disturbances  .AND.  disturbance_energy_limit /= 0.0_wp  .AND. &
1738         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
1739         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
1740
1741       CALL location_message( 'creating initial disturbances', .FALSE. )
1742       CALL disturb_field( nzb_u_inner, tend, u )
1743       CALL disturb_field( nzb_v_inner, tend, v )
1744       CALL location_message( 'finished', .TRUE. )
1745
1746       CALL location_message( 'calling pressure solver', .FALSE. )
1747       n_sor = nsor_ini
1748       !$acc data copyin( d, ddzu, ddzw, nzb_s_inner, nzb_u_inner )            &
1749       !$acc      copyin( nzb_v_inner, nzb_w_inner, p, rflags_s_inner, tend )  &
1750       !$acc      copyin( weight_pres, weight_substep )                        &
1751       !$acc      copy( tri, tric, u, v, w )
1752       CALL pres
1753       !$acc end data
1754       n_sor = nsor
1755       CALL location_message( 'finished', .TRUE. )
1756
1757    ENDIF
1758
1759!
1760!-- If required, initialize quantities needed for the plant canopy model
1761    CALL location_message( 'initializing plant canopy model', .FALSE. )
1762    IF ( plant_canopy )  CALL pcm_init
1763    CALL location_message( 'finished', .TRUE. )
1764
1765!
1766!-- If required, initialize dvrp-software
1767    IF ( dt_dvrp /= 9999999.9_wp )  CALL init_dvrp
1768
1769    IF ( ocean )  THEN
1770!
1771!--    Initialize quantities needed for the ocean model
1772       CALL init_ocean
1773
1774    ELSE
1775!
1776!--    Initialize quantities for handling cloud physics
1777!--    This routine must be called before lpm_init, because
1778!--    otherwise, array pt_d_t, needed in data_output_dvrp (called by
1779!--    lpm_init) is not defined.
1780       CALL init_cloud_physics
1781!
1782!--    Initialize bulk cloud microphysics
1783       CALL microphysics_init
1784    ENDIF
1785
1786!
1787!-- If required, initialize particles
1788    IF ( particle_advection )  CALL lpm_init
1789
1790!
1791!-- If required, initialize quantities needed for the LSM
1792    IF ( land_surface )  THEN
1793       CALL location_message( 'initializing land surface model', .FALSE. )
1794       CALL lsm_init
1795       CALL location_message( 'finished', .TRUE. )
1796    ENDIF
1797
1798!
1799!-- Initialize surface layer (done after LSM as roughness length are required
1800!-- for initialization
1801    IF ( constant_flux_layer )  THEN
1802       CALL location_message( 'initializing surface layer', .FALSE. )
1803       CALL init_surface_layer_fluxes
1804       CALL location_message( 'finished', .TRUE. )
1805    ENDIF
1806
1807!
1808!-- If required, initialize radiation model
1809    IF ( radiation )  THEN
1810       CALL location_message( 'initializing radiation model', .FALSE. )
1811       CALL radiation_init
1812       CALL location_message( 'finished', .TRUE. )
1813    ENDIF
1814   
1815!
1816!-- If required, initialize quantities needed for the wind turbine model
1817    IF ( wind_turbine )  THEN
1818       CALL location_message( 'initializing wind turbine model', .FALSE. )
1819       CALL wtm_init
1820       CALL location_message( 'finished', .TRUE. )
1821    ENDIF
1822
1823
1824!
1825!-- Initialize the ws-scheme.   
1826    IF ( ws_scheme_sca .OR. ws_scheme_mom )  CALL ws_init       
1827
1828!
1829!-- Setting weighting factors for calculation of perturbation pressure
1830!-- and turbulent quantities from the RK substeps
1831    IF ( TRIM(timestep_scheme) == 'runge-kutta-3' )  THEN      ! for RK3-method
1832
1833       weight_substep(1) = 1._wp/6._wp
1834       weight_substep(2) = 3._wp/10._wp
1835       weight_substep(3) = 8._wp/15._wp
1836
1837       weight_pres(1)    = 1._wp/3._wp
1838       weight_pres(2)    = 5._wp/12._wp
1839       weight_pres(3)    = 1._wp/4._wp
1840
1841    ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' )  THEN  ! for RK2-method
1842
1843       weight_substep(1) = 1._wp/2._wp
1844       weight_substep(2) = 1._wp/2._wp
1845         
1846       weight_pres(1)    = 1._wp/2._wp
1847       weight_pres(2)    = 1._wp/2._wp       
1848
1849    ELSE                                     ! for Euler-method
1850
1851       weight_substep(1) = 1.0_wp     
1852       weight_pres(1)    = 1.0_wp                   
1853
1854    ENDIF
1855
1856!
1857!-- Initialize Rayleigh damping factors
1858    rdf    = 0.0_wp
1859    rdf_sc = 0.0_wp
1860    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
1861       IF (  .NOT.  ocean )  THEN
1862          DO  k = nzb+1, nzt
1863             IF ( zu(k) >= rayleigh_damping_height )  THEN
1864                rdf(k) = rayleigh_damping_factor *                             &
1865                      ( SIN( pi * 0.5_wp * ( zu(k) - rayleigh_damping_height ) &
1866                             / ( zu(nzt) - rayleigh_damping_height ) )         &
1867                      )**2
1868             ENDIF
1869          ENDDO
1870       ELSE
1871          DO  k = nzt, nzb+1, -1
1872             IF ( zu(k) <= rayleigh_damping_height )  THEN
1873                rdf(k) = rayleigh_damping_factor *                             &
1874                      ( SIN( pi * 0.5_wp * ( rayleigh_damping_height - zu(k) ) &
1875                             / ( rayleigh_damping_height - zu(nzb+1) ) )       &
1876                      )**2
1877             ENDIF
1878          ENDDO
1879       ENDIF
1880    ENDIF
1881    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
1882
1883!
1884!-- Initialize the starting level and the vertical smoothing factor used for
1885!-- the external pressure gradient
1886    dp_smooth_factor = 1.0_wp
1887    IF ( dp_external )  THEN
1888!
1889!--    Set the starting level dp_level_ind_b only if it has not been set before
1890!--    (e.g. in init_grid).
1891       IF ( dp_level_ind_b == 0 )  THEN
1892          ind_array = MINLOC( ABS( dp_level_b - zu ) )
1893          dp_level_ind_b = ind_array(1) - 1 + nzb 
1894                                        ! MINLOC uses lower array bound 1
1895       ENDIF
1896       IF ( dp_smooth )  THEN
1897          dp_smooth_factor(:dp_level_ind_b) = 0.0_wp
1898          DO  k = dp_level_ind_b+1, nzt
1899             dp_smooth_factor(k) = 0.5_wp * ( 1.0_wp + SIN( pi *               &
1900                        ( REAL( k - dp_level_ind_b, KIND=wp ) /                &
1901                          REAL( nzt - dp_level_ind_b, KIND=wp ) - 0.5_wp ) ) )
1902          ENDDO
1903       ENDIF
1904    ENDIF
1905
1906!
1907!-- Initialize damping zone for the potential temperature in case of
1908!-- non-cyclic lateral boundaries. The damping zone has the maximum value
1909!-- at the inflow boundary and decreases to zero at pt_damping_width.
1910    ptdf_x = 0.0_wp
1911    ptdf_y = 0.0_wp
1912    IF ( bc_lr_dirrad )  THEN
1913       DO  i = nxl, nxr
1914          IF ( ( i * dx ) < pt_damping_width )  THEN
1915             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5_wp *              &
1916                            REAL( pt_damping_width - i * dx, KIND=wp ) / (     &
1917                            REAL( pt_damping_width, KIND=wp ) ) ) )**2 
1918          ENDIF
1919       ENDDO
1920    ELSEIF ( bc_lr_raddir )  THEN
1921       DO  i = nxl, nxr
1922          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
1923             ptdf_x(i) = pt_damping_factor *                                   &
1924                         SIN( pi * 0.5_wp *                                    &
1925                                 ( ( i - nx ) * dx + pt_damping_width ) /      &
1926                                 REAL( pt_damping_width, KIND=wp ) )**2
1927          ENDIF
1928       ENDDO 
1929    ELSEIF ( bc_ns_dirrad )  THEN
1930       DO  j = nys, nyn
1931          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
1932             ptdf_y(j) = pt_damping_factor *                                   &
1933                         SIN( pi * 0.5_wp *                                    &
1934                                 ( ( j - ny ) * dy + pt_damping_width ) /      &
1935                                 REAL( pt_damping_width, KIND=wp ) )**2
1936          ENDIF
1937       ENDDO 
1938    ELSEIF ( bc_ns_raddir )  THEN
1939       DO  j = nys, nyn
1940          IF ( ( j * dy ) < pt_damping_width )  THEN
1941             ptdf_y(j) = pt_damping_factor *                                   &
1942                         SIN( pi * 0.5_wp *                                    &
1943                                ( pt_damping_width - j * dy ) /                &
1944                                REAL( pt_damping_width, KIND=wp ) )**2
1945          ENDIF
1946       ENDDO
1947    ENDIF
1948
1949!
1950!-- Pre-set masks for regional statistics. Default is the total model domain.
1951!-- Ghost points are excluded because counting values at the ghost boundaries
1952!-- would bias the statistics
1953    rmask = 1.0_wp
1954    rmask(:,nxlg:nxl-1,:) = 0.0_wp;  rmask(:,nxr+1:nxrg,:) = 0.0_wp
1955    rmask(nysg:nys-1,:,:) = 0.0_wp;  rmask(nyn+1:nyng,:,:) = 0.0_wp
1956
1957!
1958!-- User-defined initializing actions. Check afterwards, if maximum number
1959!-- of allowed timeseries is exceeded
1960    CALL user_init
1961
1962    IF ( dots_num > dots_max )  THEN
1963       WRITE( message_string, * ) 'number of time series quantities exceeds',  &
1964                                  ' its maximum of dots_max = ', dots_max,     &
1965                                  ' &Please increase dots_max in modules.f90.'
1966       CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 )   
1967    ENDIF
1968
1969!
1970!-- Input binary data file is not needed anymore. This line must be placed
1971!-- after call of user_init!
1972    CALL close_file( 13 )
1973
1974!
1975!-- Compute total sum of active mask grid points
1976!-- and the mean surface level height for each statistic region
1977!-- ngp_2dh: number of grid points of a horizontal cross section through the
1978!--          total domain
1979!-- ngp_3d:  number of grid points of the total domain
1980    ngp_2dh_outer_l   = 0
1981    ngp_2dh_outer     = 0
1982    ngp_2dh_s_inner_l = 0
1983    ngp_2dh_s_inner   = 0
1984    ngp_2dh_l         = 0
1985    ngp_2dh           = 0
1986    ngp_3d_inner_l    = 0.0_wp
1987    ngp_3d_inner      = 0
1988    ngp_3d            = 0
1989    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
1990
1991    mean_surface_level_height   = 0.0_wp
1992    mean_surface_level_height_l = 0.0_wp
1993
1994    DO  sr = 0, statistic_regions
1995       DO  i = nxl, nxr
1996          DO  j = nys, nyn
1997             IF ( rmask(j,i,sr) == 1.0_wp )  THEN
1998!
1999!--             All xy-grid points
2000                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
2001                mean_surface_level_height_l(sr) = mean_surface_level_height_l(sr) &
2002                                                  + zw(nzb_s_inner(j,i))
2003!
2004!--             xy-grid points above topography
2005                DO  k = nzb_s_outer(j,i), nz + 1
2006                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1
2007                ENDDO
2008                DO  k = nzb_s_inner(j,i), nz + 1
2009                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) + 1
2010                ENDDO
2011!
2012!--             All grid points of the total domain above topography
2013                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr)                        &
2014                                     + ( nz - nzb_s_inner(j,i) + 2 )
2015             ENDIF
2016          ENDDO
2017       ENDDO
2018    ENDDO
2019
2020    sr = statistic_regions + 1
2021#if defined( __parallel )
2022    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2023    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,    &
2024                        comm2d, ierr )
2025    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2026    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr,   &
2027                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
2028    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2029    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),          &
2030                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
2031    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2032    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL,  &
2033                        MPI_SUM, comm2d, ierr )
2034    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
2035    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2036    CALL MPI_ALLREDUCE( mean_surface_level_height_l(0),                        &
2037                        mean_surface_level_height(0), sr, MPI_REAL,            &
2038                        MPI_SUM, comm2d, ierr )
2039    mean_surface_level_height = mean_surface_level_height / REAL( ngp_2dh )
2040#else
2041    ngp_2dh         = ngp_2dh_l
2042    ngp_2dh_outer   = ngp_2dh_outer_l
2043    ngp_2dh_s_inner = ngp_2dh_s_inner_l
2044    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
2045    mean_surface_level_height = mean_surface_level_height_l / REAL( ngp_2dh_l )
2046#endif
2047
2048    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * &
2049             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
2050
2051!
2052!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
2053!-- buoyancy, etc. A zero value will occur for cases where all grid points of
2054!-- the respective subdomain lie below the surface topography
2055    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   ) 
2056    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )),             &
2057                           ngp_3d_inner(:) )
2058    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) 
2059
2060    DEALLOCATE( mean_surface_level_height_l, ngp_2dh_l, ngp_2dh_outer_l,       &
2061                ngp_3d_inner_l, ngp_3d_inner_tmp )
2062
2063    CALL location_message( 'leaving init_3d_model', .TRUE. )
2064
2065 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.