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

Last change on this file since 2011 was 2011, checked in by kanani, 8 years ago

changes related to steering and formating of urban surface model

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