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

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