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

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

last commit documented

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