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

Last change on this file since 985 was 979, checked in by fricke, 12 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 59.6 KB
Line 
1#if defined( __ibmy_special )
2@PROCESS NOOPTimize
3#endif
4 SUBROUTINE init_3d_model
5
6!------------------------------------------------------------------------------!
7! Current revisions:
8! ------------------
9!
10!
11! Former revisions:
12! -----------------
13! $Id: init_3d_model.f90 979 2012-08-09 08:50:11Z maronga $
14!
15! 978 2012-08-09 08:28:32Z fricke
16! outflow damping layer removed
17! roughness length for scalar quantites z0h added
18! damping zone for the potential temperatur in case of non-cyclic lateral
19! boundaries added
20! initialization of ptdf_x, ptdf_y
21! initialization of c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l
22!
23! 849 2012-03-15 10:35:09Z raasch
24! init_particles renamed lpm_init
25!
26! 825 2012-02-19 03:03:44Z raasch
27! wang_collision_kernel renamed wang_kernel
28!
29! 790 2011-11-29 03:11:20Z raasch
30! diss is also allocated in case that the Wang kernel is used
31!
32! 787 2011-11-28 12:49:05Z heinze $
33! bugfix: call init_advec in every case - not only for inital runs
34!
35! 785 2011-11-28 09:47:19Z raasch
36! initialization of rdf_sc
37!
38! 767 2011-10-14 06:39:12Z raasch
39! adjustments concerning implementation of prescribed u,v-profiles
40! bugfix: dirichlet_0 conditions for ug/vg moved to check_parameters
41!
42! 759 2011-09-15 13:58:31Z raasch
43! Splitting of parallel I/O in blocks of PEs
44! Bugfix: No zero assignments to volume_flow_initial and volume_flow_area in
45! case of normal restart runs.
46!
47! 713 2011-03-30 14:21:21Z suehring
48! weight_substep and weight_pres are given as fractions.
49!
50! 709 2011-03-30 09:31:40Z raasch
51! formatting adjustments
52!
53! 707 2011-03-29 11:39:40Z raasch
54! p_sub renamed p_loc and allocated depending on the chosen pressure solver,
55! initial assignments of zero to array p for iterative solvers only,
56! bc_lr/ns replaced by bc_lr/ns_dirrad/raddir
57!
58! 680 2011-02-04 23:16:06Z gryschka
59! bugfix: volume_flow_control
60!
61! 673 2011-01-18 16:19:48Z suehring
62! weight_substep (moved from advec_ws) and weight_pres added.
63! Allocate p_sub when using Multigrid or SOR solver.
64! Call of ws_init moved behind the if requests.
65!
66! 667 2010-12-23 12:06:00Z suehring/gryschka
67! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in loops and
68! allocation of arrays. Calls of exchange_horiz are modified.
69! Call ws_init to initialize arrays needed for calculating statisticas and for
70! optimization when ws-scheme is used.
71! Initial volume flow is now calculated by using the variable hom_sum.
72! Therefore the correction of initial volume flow for non-flat topography
73! removed (removed u_nzb_p1_for_vfc and v_nzb_p1_for_vfc)
74! Changed surface boundary conditions for u and v in case of ibc_uv_b == 0 from
75! mirror to Dirichlet boundary conditions (u=v=0), so that k=nzb is
76! representative for the height z0.
77! Bugfix: type conversion of '1' to 64bit for the MAX function (ngp_3d_inner)
78!
79! 622 2010-12-10 08:08:13Z raasch
80! optional barriers included in order to speed up collective operations
81!
82! 560 2010-09-09 10:06:09Z weinreis
83! bugfix: correction of calculating ngp_3d for 64 bit
84!
85! 485 2010-02-05 10:57:51Z raasch
86! calculation of ngp_3d + ngp_3d_inner changed because they have now 64 bit
87!
88! 407 2009-12-01 15:01:15Z maronga
89! var_ts is replaced by dots_max
90! Enabled passive scalar/humidity wall fluxes for non-flat topography
91!
92! 388 2009-09-23 09:40:33Z raasch
93! Initialization of prho added.
94! bugfix: correction of initial volume flow for non-flat topography
95! bugfix: zero initialization of arrays within buildings for 'cyclic_fill'
96! bugfix: avoid that ngp_2dh_s_inner becomes zero
97! initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill', now
98! independent of turbulent_inflow
99! Output of messages replaced by message handling routine.
100! Set the starting level and the vertical smoothing factor used for
101! the external pressure gradient
102! +conserve_volume_flow_mode: 'default', 'initial_profiles', 'inflow_profile'
103! and 'bulk_velocity'
104! If the inversion height calculated by the prerun is zero,
105! inflow_damping_height must be explicitly specified.
106!
107! 181 2008-07-30 07:07:47Z raasch
108! bugfix: zero assignments to tendency arrays in case of restarts,
109! further extensions and modifications in the initialisation of the plant
110! canopy model,
111! allocation of hom_sum moved to parin, initialization of spectrum_x|y directly
112! after allocating theses arrays,
113! read data for recycling added as new initialization option,
114! dummy allocation for diss
115!
116! 138 2007-11-28 10:03:58Z letzel
117! New counter ngp_2dh_s_inner.
118! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
119! Corrected calculation of initial volume flow for 'set_1d-model_profiles' and
120! 'set_constant_profiles' in case of buildings in the reference cross-sections.
121!
122! 108 2007-08-24 15:10:38Z letzel
123! Flux initialization in case of coupled runs, +momentum fluxes at top boundary,
124! +arrays for phase speed c_u, c_v, c_w, indices for u|v|w_m_l|r changed
125! +qswst_remote in case of atmosphere model with humidity coupled to ocean
126! Rayleigh damping for ocean, optionally calculate km and kh from initial
127! TKE e_init
128!
129! 97 2007-06-21 08:23:15Z raasch
130! Initialization of salinity, call of init_ocean
131!
132! 87 2007-05-22 15:46:47Z raasch
133! var_hom and var_sum renamed pr_palm
134!
135! 75 2007-03-22 09:54:05Z raasch
136! Arrays for radiation boundary conditions are allocated (u_m_l, u_m_r, etc.),
137! bugfix for cases with the outflow damping layer extending over more than one
138! subdomain, moisture renamed humidity,
139! new initializing action "by_user" calls user_init_3d_model,
140! precipitation_amount/rate, ts_value are allocated, +module netcdf_control,
141! initial velocities at nzb+1 are regarded for volume
142! flow control in case they have been set zero before (to avoid small timesteps)
143! -uvmean_outflow, uxrp, vynp eliminated
144!
145! 19 2007-02-23 04:53:48Z raasch
146! +handling of top fluxes
147!
148! RCS Log replace by Id keyword, revision history cleaned up
149!
150! Revision 1.49  2006/08/22 15:59:07  raasch
151! No optimization of this file on the ibmy (Yonsei Univ.)
152!
153! Revision 1.1  1998/03/09 16:22:22  raasch
154! Initial revision
155!
156!
157! Description:
158! ------------
159! Allocation of arrays and initialization of the 3D model via
160! a) pre-run the 1D model
161! or
162! b) pre-set constant linear profiles
163! or
164! c) read values of a previous run
165!------------------------------------------------------------------------------!
166
167    USE advec_ws
168    USE arrays_3d
169    USE averaging
170    USE cloud_parameters
171    USE constants
172    USE control_parameters
173    USE cpulog
174    USE grid_variables
175    USE indices
176    USE interfaces
177    USE model_1d
178    USE netcdf_control
179    USE particle_attributes
180    USE pegrid
181    USE profil_parameter
182    USE random_function_mod
183    USE statistics
184
185    IMPLICIT NONE
186
187    INTEGER ::  i, ind_array(1), j, k, sr
188
189    INTEGER, DIMENSION(:), ALLOCATABLE ::  ngp_2dh_l
190
191    INTEGER, DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l,  &
192         ngp_2dh_s_inner_l
193
194    REAL ::  a, b
195
196    REAL, DIMENSION(1:2) ::  volume_flow_area_l, volume_flow_initial_l
197
198    REAL, DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l, ngp_3d_inner_tmp
199
200
201!
202!-- Allocate arrays
203    ALLOCATE( ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions), &
204              ngp_3d(0:statistic_regions),                                  &
205              ngp_3d_inner(0:statistic_regions),                            &
206              ngp_3d_inner_l(0:statistic_regions),                          &
207              ngp_3d_inner_tmp(0:statistic_regions),                        &
208              sums_divnew_l(0:statistic_regions),                           &
209              sums_divold_l(0:statistic_regions) )
210    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
211    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                 &
212              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),               &
213              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),               &
214              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),             &
215              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),           &
216              sums(nzb:nzt+1,pr_palm+max_pr_user),                          &
217              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),   &
218              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1), &
219              sums_up_fraction_l(10,3,0:statistic_regions),                 &
220              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                &
221              ts_value(dots_max,0:statistic_regions) )
222    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
223
224    ALLOCATE( rif_1(nysg:nyng,nxlg:nxrg), shf_1(nysg:nyng,nxlg:nxrg), &
225              ts(nysg:nyng,nxlg:nxrg), tswst_1(nysg:nyng,nxlg:nxrg),  &
226              us(nysg:nyng,nxlg:nxrg), usws_1(nysg:nyng,nxlg:nxrg),   &
227              uswst_1(nysg:nyng,nxlg:nxrg),                               &
228              vsws_1(nysg:nyng,nxlg:nxrg),                                &
229              vswst_1(nysg:nyng,nxlg:nxrg), z0(nysg:nyng,nxlg:nxrg),      &
230              z0h(nysg:nyng,nxlg:nxrg) )
231
232    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
233!
234!--    Leapfrog scheme needs two timelevels of diffusion quantities
235       ALLOCATE( rif_2(nysg:nyng,nxlg:nxrg),   &
236                 shf_2(nysg:nyng,nxlg:nxrg),   &
237                 tswst_2(nysg:nyng,nxlg:nxrg), &
238                 usws_2(nysg:nyng,nxlg:nxrg),  &
239                 uswst_2(nysg:nyng,nxlg:nxrg), &
240                 vswst_2(nysg:nyng,nxlg:nxrg), &
241                 vsws_2(nysg:nyng,nxlg:nxrg) )
242    ENDIF
243
244    ALLOCATE( d(nzb+1:nzta,nys:nyna,nxl:nxra),         &
245              e_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
246              e_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
247              e_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
248              kh_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
249              km_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
250              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),    &
251              pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
252              pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
253              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
254              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
255              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
256              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
257              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
258              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
259              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
260              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
261              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
262              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
263              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
264!
265!-- Following array is required for perturbation pressure within the iterative
266!-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds
267!-- the weighted average of the substeps and cannot be used in the Poisson
268!-- solver.
269    IF ( psolver == 'sor' )  THEN
270       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
271    ELSEIF ( psolver == 'multigrid' )  THEN
272!
273!--    For performance reasons, multigrid is using one ghost layer only
274       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
275    ENDIF
276
277    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
278       ALLOCATE( kh_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
279                 km_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
280    ENDIF
281
282    IF ( humidity  .OR.  passive_scalar ) THEN
283!
284!--    2D-humidity/scalar arrays
285       ALLOCATE ( qs(nysg:nyng,nxlg:nxrg),     &
286                  qsws_1(nysg:nyng,nxlg:nxrg), &
287                  qswst_1(nysg:nyng,nxlg:nxrg) )
288
289       IF ( timestep_scheme(1:5) /= 'runge' )  THEN
290          ALLOCATE( qsws_2(nysg:nyng,nxlg:nxrg), &
291                    qswst_2(nysg:nyng,nxlg:nxrg) )
292       ENDIF
293!
294!--    3D-humidity/scalar arrays
295       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
296                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
297                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
298
299!
300!--    3D-arrays needed for humidity only
301       IF ( humidity )  THEN
302          ALLOCATE( vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
303
304          IF ( timestep_scheme(1:5) /= 'runge' )  THEN
305             ALLOCATE( vpt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
306          ENDIF
307
308          IF ( cloud_physics ) THEN
309!
310!--          Liquid water content
311             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
312!
313!--          Precipitation amount and rate (only needed if output is switched)
314             ALLOCATE( precipitation_amount(nysg:nyng,nxlg:nxrg), &
315                       precipitation_rate(nysg:nyng,nxlg:nxrg) )
316          ENDIF
317
318          IF ( cloud_droplets )  THEN
319!
320!--          Liquid water content, change in liquid water content,
321!--          real volume of particles (with weighting), volume of particles
322             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
323                        ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
324                        ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
325                        ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
326          ENDIF
327
328       ENDIF
329
330    ENDIF
331
332    IF ( ocean )  THEN
333       ALLOCATE( saswsb_1(nysg:nyng,nxlg:nxrg), &
334                 saswst_1(nysg:nyng,nxlg:nxrg) )
335       ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
336                 rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
337                 sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),   &
338                 sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),   &
339                 sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
340       prho => prho_1
341       rho  => rho_1  ! routines calc_mean_profile and diffusion_e require
342                      ! density to be apointer
343       IF ( humidity_remote )  THEN
344          ALLOCATE( qswst_remote(nysg:nyng,nxlg:nxrg))
345          qswst_remote = 0.0
346       ENDIF
347    ENDIF
348
349!
350!-- 3D-array for storing the dissipation, needed for calculating the sgs
351!-- particle velocities
352    IF ( use_sgs_for_particles  .OR.  wang_kernel )  THEN
353       ALLOCATE ( diss(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
354    ELSE
355       ALLOCATE ( diss(2,2,2) )  ! required because diss is used as a
356                                 ! formal parameter
357    ENDIF
358
359    IF ( dt_dosp /= 9999999.9 )  THEN
360       ALLOCATE( spectrum_x( 1:nx/2, 1:10, 1:10 ), &
361                 spectrum_y( 1:ny/2, 1:10, 1:10 ) )
362       spectrum_x = 0.0
363       spectrum_y = 0.0
364    ENDIF
365
366!
367!-- 3D-arrays for the leaf area density and the canopy drag coefficient
368    IF ( plant_canopy ) THEN
369       ALLOCATE ( lad_s(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
370                  lad_u(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
371                  lad_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
372                  lad_w(nzb:nzt+1,nysg:nyng,nxlg:nxrg),  &
373                  cdc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
374
375       IF ( passive_scalar ) THEN
376          ALLOCATE ( sls(nzb:nzt+1,nysg:nyng,nxlg:nxrg),   &
377                     sec(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) 
378       ENDIF
379
380       IF ( cthf /= 0.0 ) THEN
381          ALLOCATE ( lai(nzb:nzt+1,nysg:nyng,nxlg:nxrg),   &
382                     canopy_heat_flux(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
383       ENDIF
384
385    ENDIF
386
387!
388!-- 4D-array for storing the Rif-values at vertical walls
389    IF ( topography /= 'flat' )  THEN
390       ALLOCATE( rif_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg,1:4) )
391       rif_wall = 0.0
392    ENDIF
393
394!
395!-- Arrays to store velocity data from t-dt and the phase speeds which
396!-- are needed for radiation boundary conditions
397    IF ( outflow_l )  THEN
398       ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2), &
399                 v_m_l(nzb:nzt+1,nysg:nyng,0:1), &
400                 w_m_l(nzb:nzt+1,nysg:nyng,0:1) )
401    ENDIF
402    IF ( outflow_r )  THEN
403       ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx), &
404                 v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx), &
405                 w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) )
406    ENDIF
407    IF ( outflow_l  .OR.  outflow_r )  THEN
408       ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng), &
409                 c_w(nzb:nzt+1,nysg:nyng) )
410    ENDIF
411    IF ( outflow_s )  THEN
412       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg), &
413                 v_m_s(nzb:nzt+1,1:2,nxlg:nxrg), &
414                 w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) )
415    ENDIF
416    IF ( outflow_n )  THEN
417       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg), &
418                 v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg), &
419                 w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) )
420    ENDIF
421    IF ( outflow_s  .OR.  outflow_n )  THEN
422       ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg), &
423                 c_w(nzb:nzt+1,nxlg:nxrg) )
424    ENDIF
425    IF ( outflow_l  .OR.  outflow_r .OR. outflow_s .OR. outflow_n )  THEN
426       ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) )                   
427       ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) )
428    ENDIF
429
430
431!
432!-- Initial assignment of the pointers
433    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
434
435       rif_m   => rif_1;    rif   => rif_2
436       shf_m   => shf_1;    shf   => shf_2
437       tswst_m => tswst_1;  tswst => tswst_2
438       usws_m  => usws_1;   usws  => usws_2
439       uswst_m => uswst_1;  uswst => uswst_2
440       vsws_m  => vsws_1;   vsws  => vsws_2
441       vswst_m => vswst_1;  vswst => vswst_2
442       e_m  => e_1;   e  => e_2;   e_p  => e_3;   te_m  => e_3
443       kh_m => kh_1;  kh => kh_2
444       km_m => km_1;  km => km_2
445       pt_m => pt_1;  pt => pt_2;  pt_p => pt_3;  tpt_m => pt_3
446       u_m  => u_1;   u  => u_2;   u_p  => u_3;   tu_m  => u_3
447       v_m  => v_1;   v  => v_2;   v_p  => v_3;   tv_m  => v_3
448       w_m  => w_1;   w  => w_2;   w_p  => w_3;   tw_m  => w_3
449
450       IF ( humidity  .OR.  passive_scalar )  THEN
451          qsws_m  => qsws_1;   qsws  => qsws_2
452          qswst_m => qswst_1;  qswst => qswst_2
453          q_m    => q_1;     q    => q_2;     q_p => q_3;     tq_m => q_3
454          IF ( humidity )        vpt_m  => vpt_1;   vpt  => vpt_2
455          IF ( cloud_physics )   ql   => ql_1
456          IF ( cloud_droplets )  THEN
457             ql   => ql_1
458             ql_c => ql_2
459          ENDIF
460       ENDIF
461
462    ELSE
463
464       rif   => rif_1
465       shf   => shf_1
466       tswst => tswst_1
467       usws  => usws_1
468       uswst => uswst_1
469       vsws  => vsws_1
470       vswst => vswst_1
471       e     => e_1;   e_p  => e_2;   te_m  => e_3;   e_m  => e_3
472       kh    => kh_1
473       km    => km_1
474       pt    => pt_1;  pt_p => pt_2;  tpt_m => pt_3;  pt_m => pt_3
475       u     => u_1;   u_p  => u_2;   tu_m  => u_3;   u_m  => u_3
476       v     => v_1;   v_p  => v_2;   tv_m  => v_3;   v_m  => v_3
477       w     => w_1;   w_p  => w_2;   tw_m  => w_3;   w_m  => w_3
478
479       IF ( humidity  .OR.  passive_scalar )  THEN
480          qsws   => qsws_1
481          qswst  => qswst_1
482          q      => q_1;     q_p  => q_2;     tq_m  => q_3;    q_m => q_3
483          IF ( humidity )        vpt  => vpt_1
484          IF ( cloud_physics )   ql   => ql_1
485          IF ( cloud_droplets )  THEN
486             ql   => ql_1
487             ql_c => ql_2
488          ENDIF
489       ENDIF
490
491       IF ( ocean )  THEN
492          saswsb => saswsb_1
493          saswst => saswst_1
494          sa     => sa_1;    sa_p => sa_2;    tsa_m => sa_3
495       ENDIF
496
497    ENDIF
498   
499!
500!-- Allocate arrays containing the RK coefficient for calculation of
501!-- perturbation pressure and turbulent fluxes. At this point values are
502!-- set for pressure calculation during initialization (where no timestep
503!-- is done). Further below the values needed within the timestep scheme
504!-- will be set.
505    ALLOCATE( weight_substep(1:intermediate_timestep_count_max), &
506              weight_pres(1:intermediate_timestep_count_max) )
507    weight_substep = 1.0
508    weight_pres    = 1.0
509    intermediate_timestep_count = 1  ! needed when simulated_time = 0.0
510       
511!
512!-- Initialize model variables
513    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.  &
514         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
515!
516!--    First model run of a possible job queue.
517!--    Initial profiles of the variables must be computes.
518       IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
519!
520!--       Use solutions of the 1D model as initial profiles,
521!--       start 1D model
522          CALL init_1d_model
523!
524!--       Transfer initial profiles to the arrays of the 3D model
525          DO  i = nxlg, nxrg
526             DO  j = nysg, nyng
527                e(:,j,i)  = e1d
528                kh(:,j,i) = kh1d
529                km(:,j,i) = km1d
530                pt(:,j,i) = pt_init
531                u(:,j,i)  = u1d
532                v(:,j,i)  = v1d
533             ENDDO
534          ENDDO
535
536          IF ( humidity  .OR.  passive_scalar )  THEN
537             DO  i = nxlg, nxrg
538                DO  j = nysg, nyng
539                   q(:,j,i) = q_init
540                ENDDO
541             ENDDO
542          ENDIF
543
544          IF ( .NOT. constant_diffusion )  THEN
545             DO  i = nxlg, nxrg
546                DO  j = nysg, nyng
547                   e(:,j,i)  = e1d
548                ENDDO
549             ENDDO
550!
551!--          Store initial profiles for output purposes etc.
552             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
553
554             IF ( prandtl_layer )  THEN
555                rif  = rif1d(nzb+1)
556                ts   = 0.0  ! could actually be computed more accurately in the
557                            ! 1D model. Update when opportunity arises.
558                us   = us1d
559                usws = usws1d
560                vsws = vsws1d
561             ELSE
562                ts   = 0.0  ! must be set, because used in
563                rif  = 0.0  ! flowste
564                us   = 0.0
565                usws = 0.0
566                vsws = 0.0
567             ENDIF
568
569          ELSE
570             e    = 0.0  ! must be set, because used in
571             rif  = 0.0  ! flowste
572             ts   = 0.0
573             us   = 0.0
574             usws = 0.0
575             vsws = 0.0
576          ENDIF
577          uswst = top_momentumflux_u
578          vswst = top_momentumflux_v
579
580!
581!--       In every case qs = 0.0 (see also pt)
582!--       This could actually be computed more accurately in the 1D model.
583!--       Update when opportunity arises!
584          IF ( humidity  .OR.  passive_scalar )  qs = 0.0
585
586!
587!--       inside buildings set velocities back to zero
588          IF ( topography /= 'flat' )  THEN
589             DO  i = nxl-1, nxr+1
590                DO  j = nys-1, nyn+1
591                   u(nzb:nzb_u_inner(j,i),j,i) = 0.0
592                   v(nzb:nzb_v_inner(j,i),j,i) = 0.0
593                ENDDO
594             ENDDO
595             
596!
597!--          WARNING: The extra boundary conditions set after running the
598!--          -------  1D model impose an error on the divergence one layer
599!--                   below the topography; need to correct later
600!--          ATTENTION: Provisional correction for Piacsek & Williams
601!--          ---------  advection scheme: keep u and v zero one layer below
602!--                     the topography.
603!
604!--           Following was removed, because mirror boundary condition are
605!--           replaced by dirichlet boundary conditions
606!
607!             IF ( ibc_uv_b == 0 )  THEN
608!!
609!!--             Satisfying the Dirichlet condition with an extra layer below
610!!--             the surface where the u and v component change their sign.
611!                DO  i = nxl-1, nxr+1
612!                   DO  j = nys-1, nyn+1
613!                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = -u(1,j,i)
614!                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = -v(1,j,i)
615!                   ENDDO
616!                ENDDO
617!
618!             ELSE
619             IF ( ibc_uv_b == 1 )  THEN
620!
621!--             Neumann condition
622                DO  i = nxl-1, nxr+1
623                   DO  j = nys-1, nyn+1
624                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i)
625                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i)
626                   ENDDO
627                ENDDO
628
629             ENDIF
630
631          ENDIF
632
633       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 ) &
634       THEN
635!
636!--       Use constructed initial profiles (velocity constant with height,
637!--       temperature profile with constant gradient)
638          DO  i = nxlg, nxrg
639             DO  j = nysg, nyng
640                pt(:,j,i) = pt_init
641                u(:,j,i)  = u_init
642                v(:,j,i)  = v_init
643             ENDDO
644          ENDDO
645
646!
647!--       Set initial horizontal velocities at the lowest computational grid
648!--       levels to zero in order to avoid too small time steps caused by the
649!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
650!--       in the limiting formula!). The original values are stored to be later
651!--       used for volume flow control.
652          DO  i = nxlg, nxrg
653             DO  j = nysg, nyng
654                u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0
655                v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0
656             ENDDO
657          ENDDO
658
659          IF ( humidity  .OR.  passive_scalar )  THEN
660             DO  i = nxlg, nxrg
661                DO  j = nysg, nyng
662                   q(:,j,i) = q_init
663                ENDDO
664             ENDDO
665          ENDIF
666
667          IF ( ocean )  THEN
668             DO  i = nxlg, nxrg
669                DO  j = nysg, nyng
670                   sa(:,j,i) = sa_init
671                ENDDO
672             ENDDO
673          ENDIF
674         
675          IF ( constant_diffusion )  THEN
676             km   = km_constant
677             kh   = km / prandtl_number
678             e    = 0.0
679          ELSEIF ( e_init > 0.0 )  THEN
680             DO  k = nzb+1, nzt
681                km(k,:,:) = 0.1 * l_grid(k) * SQRT( e_init )
682             ENDDO
683             km(nzb,:,:)   = km(nzb+1,:,:)
684             km(nzt+1,:,:) = km(nzt,:,:)
685             kh   = km / prandtl_number
686             e    = e_init
687          ELSE
688             IF ( .NOT. ocean )  THEN
689                kh   = 0.01   ! there must exist an initial diffusion, because
690                km   = 0.01   ! otherwise no TKE would be produced by the
691                              ! production terms, as long as not yet
692                              ! e = (u*/cm)**2 at k=nzb+1
693             ELSE
694                kh   = 0.00001
695                km   = 0.00001
696             ENDIF
697             e    = 0.0
698          ENDIF
699          rif   = 0.0
700          ts    = 0.0
701          us    = 0.0
702          usws  = 0.0
703          uswst = top_momentumflux_u
704          vsws  = 0.0
705          vswst = top_momentumflux_v
706          IF ( humidity  .OR.  passive_scalar ) qs = 0.0
707
708!
709!--       Compute initial temperature field and other constants used in case
710!--       of a sloping surface
711          IF ( sloping_surface )  CALL init_slope
712
713       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 ) &
714       THEN
715!
716!--       Initialization will completely be done by the user
717          CALL user_init_3d_model
718
719       ENDIF
720!
721!--    Bottom boundary
722       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
723          u(nzb,:,:) = 0.0
724          v(nzb,:,:) = 0.0
725       ENDIF
726
727!
728!--    Apply channel flow boundary condition
729       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
730          u(nzt+1,:,:) = 0.0
731          v(nzt+1,:,:) = 0.0
732       ENDIF
733
734!
735!--    Calculate virtual potential temperature
736       IF ( humidity ) vpt = pt * ( 1.0 + 0.61 * q )
737
738!
739!--    Store initial profiles for output purposes etc.
740       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
741       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
742       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
743          hom(nzb,1,5,:) = 0.0   
744          hom(nzb,1,6,:) = 0.0 
745       ENDIF
746       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
747       hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 )
748       hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 )
749
750       IF ( ocean )  THEN
751!
752!--       Store initial salinity profile
753          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
754       ENDIF
755
756       IF ( humidity )  THEN
757!
758!--       Store initial profile of total water content, virtual potential
759!--       temperature
760          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
761          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
762          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
763!
764!--          Store initial profile of specific humidity and potential
765!--          temperature
766             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
767             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
768          ENDIF
769       ENDIF
770
771       IF ( passive_scalar )  THEN
772!
773!--       Store initial scalar profile
774          hom(:,1,26,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
775       ENDIF
776
777!
778!--    Initialize fluxes at bottom surface
779       IF ( use_surface_fluxes )  THEN
780
781          IF ( constant_heatflux )  THEN
782!
783!--          Heat flux is prescribed
784             IF ( random_heatflux )  THEN
785                CALL disturb_heatflux
786             ELSE
787                shf = surface_heatflux
788!
789!--             Over topography surface_heatflux is replaced by wall_heatflux(0)
790                IF ( TRIM( topography ) /= 'flat' )  THEN
791                   DO  i = nxlg, nxrg
792                      DO  j = nysg, nyng
793                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
794                            shf(j,i) = wall_heatflux(0)
795                         ENDIF
796                      ENDDO
797                   ENDDO
798                ENDIF
799             ENDIF
800             IF ( ASSOCIATED( shf_m ) )  shf_m = shf
801          ENDIF
802
803!
804!--       Determine the near-surface water flux
805          IF ( humidity  .OR.  passive_scalar )  THEN
806             IF ( constant_waterflux )  THEN
807                qsws   = surface_waterflux
808!
809!--             Over topography surface_waterflux is replaced by
810!--             wall_humidityflux(0)
811                IF ( TRIM( topography ) /= 'flat' )  THEN
812                   wall_qflux = wall_humidityflux
813                   DO  i = nxlg, nxrg
814                      DO  j = nysg, nyng
815                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
816                            qsws(j,i) = wall_qflux(0)
817                         ENDIF
818                      ENDDO
819                   ENDDO
820                ENDIF
821                IF ( ASSOCIATED( qsws_m ) )  qsws_m = qsws
822             ENDIF
823          ENDIF
824
825       ENDIF
826
827!
828!--    Initialize fluxes at top surface
829!--    Currently, only the heatflux and salinity flux can be prescribed.
830!--    The latent flux is zero in this case!
831       IF ( use_top_fluxes )  THEN
832
833          IF ( constant_top_heatflux )  THEN
834!
835!--          Heat flux is prescribed
836             tswst = top_heatflux
837             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
838
839             IF ( humidity  .OR.  passive_scalar )  THEN
840                qswst = 0.0
841                IF ( ASSOCIATED( qswst_m ) )  qswst_m = qswst
842             ENDIF
843
844             IF ( ocean )  THEN
845                saswsb = bottom_salinityflux
846                saswst = top_salinityflux
847             ENDIF
848          ENDIF
849
850!
851!--       Initialization in case of a coupled model run
852          IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
853             tswst = 0.0
854             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
855          ENDIF
856
857       ENDIF
858
859!
860!--    Initialize Prandtl layer quantities
861       IF ( prandtl_layer )  THEN
862
863          z0 = roughness_length
864          z0h = z0h_factor * z0
865
866          IF ( .NOT. constant_heatflux )  THEN 
867!
868!--          Surface temperature is prescribed. Here the heat flux cannot be
869!--          simply estimated, because therefore rif, u* and theta* would have
870!--          to be computed by iteration. This is why the heat flux is assumed
871!--          to be zero before the first time step. It approaches its correct
872!--          value in the course of the first few time steps.
873             shf   = 0.0
874             IF ( ASSOCIATED( shf_m ) )  shf_m = 0.0
875          ENDIF
876
877          IF ( humidity  .OR.  passive_scalar )  THEN
878             IF ( .NOT. constant_waterflux )  THEN
879                qsws   = 0.0
880                IF ( ASSOCIATED( qsws_m ) )   qsws_m = 0.0
881             ENDIF
882          ENDIF
883
884       ENDIF
885
886
887!
888!--    For the moment, vertical velocity is zero
889       w = 0.0
890
891!
892!--    Initialize array sums (must be defined in first call of pres)
893       sums = 0.0
894
895!
896!--    In case of iterative solvers, p must get an initial value
897       IF ( psolver == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0
898
899!
900!--    Treating cloud physics, liquid water content and precipitation amount
901!--    are zero at beginning of the simulation
902       IF ( cloud_physics )  THEN
903          ql = 0.0
904          IF ( precipitation )  precipitation_amount = 0.0
905       ENDIF
906!
907!--    Impose vortex with vertical axis on the initial velocity profile
908       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
909          CALL init_rankine
910       ENDIF
911
912!
913!--    Impose temperature anomaly (advection test only)
914       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
915          CALL init_pt_anomaly
916       ENDIF
917
918!
919!--    If required, change the surface temperature at the start of the 3D run
920       IF ( pt_surface_initial_change /= 0.0 )  THEN
921          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
922       ENDIF
923
924!
925!--    If required, change the surface humidity/scalar at the start of the 3D
926!--    run
927       IF ( ( humidity .OR. passive_scalar ) .AND. &
928            q_surface_initial_change /= 0.0 )  THEN
929          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
930       ENDIF
931
932!
933!--    Initialize the random number generator (from numerical recipes)
934       CALL random_function_ini
935
936!
937!--    Initialize old and new time levels.
938       IF ( timestep_scheme(1:5) /= 'runge' )  THEN
939          e_m = e; pt_m = pt; u_m = u; v_m = v; w_m = w; kh_m = kh; km_m = km
940       ELSE
941          te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0
942       ENDIF
943       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
944
945       IF ( humidity  .OR.  passive_scalar )  THEN
946          IF ( ASSOCIATED( q_m ) )  q_m = q
947          IF ( timestep_scheme(1:5) == 'runge' )  tq_m = 0.0
948          q_p = q
949          IF ( humidity  .AND.  ASSOCIATED( vpt_m ) )  vpt_m = vpt
950       ENDIF
951
952       IF ( ocean )  THEN
953          tsa_m = 0.0
954          sa_p  = sa
955       ENDIF
956       
957
958    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.    &
959         TRIM( initializing_actions ) == 'cyclic_fill' )  &
960    THEN
961!
962!--    When reading data for cyclic fill of 3D prerun data files, read
963!--    some of the global variables from the restart file which are required
964!--    for initializing the inflow
965       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
966
967          DO  i = 0, io_blocks-1
968             IF ( i == io_group )  THEN
969                CALL read_parts_of_var_list
970                CALL close_file( 13 )
971             ENDIF
972#if defined( __parallel )
973             CALL MPI_BARRIER( comm2d, ierr )
974#endif
975          ENDDO
976
977       ENDIF
978
979!
980!--    Read binary data from restart file
981       DO  i = 0, io_blocks-1
982          IF ( i == io_group )  THEN
983             CALL read_3d_binary
984          ENDIF
985#if defined( __parallel )
986          CALL MPI_BARRIER( comm2d, ierr )
987#endif
988       ENDDO
989
990!
991!--    Initialization of the turbulence recycling method
992       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.  &
993            turbulent_inflow )  THEN
994!
995!--       First store the profiles to be used at the inflow.
996!--       These profiles are the (temporally) and horizontally averaged vertical
997!--       profiles from the prerun. Alternatively, prescribed profiles
998!--       for u,v-components can be used.
999          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,5) )
1000
1001          IF ( use_prescribed_profile_data )  THEN
1002             mean_inflow_profiles(:,1) = u_init            ! u
1003             mean_inflow_profiles(:,2) = v_init            ! v
1004          ELSE
1005             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
1006             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
1007          ENDIF
1008          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
1009          mean_inflow_profiles(:,5) = hom_sum(:,8,0)       ! e
1010
1011!
1012!--       If necessary, adjust the horizontal flow field to the prescribed
1013!--       profiles
1014          IF ( use_prescribed_profile_data )  THEN
1015             DO  i = nxlg, nxrg
1016                DO  j = nysg, nyng
1017                   DO  k = nzb, nzt+1
1018                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1019                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
1020                   ENDDO
1021                ENDDO
1022             ENDDO
1023          ENDIF
1024
1025!
1026!--       Use these mean profiles at the inflow (provided that Dirichlet
1027!--       conditions are used)
1028          IF ( inflow_l )  THEN
1029             DO  j = nysg, nyng
1030                DO  k = nzb, nzt+1
1031                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1032                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
1033                   w(k,j,nxlg:-1)  = 0.0
1034                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
1035                   e(k,j,nxlg:-1)  = mean_inflow_profiles(k,5)
1036                ENDDO
1037             ENDDO
1038          ENDIF
1039
1040!
1041!--       Calculate the damping factors to be used at the inflow. For a
1042!--       turbulent inflow the turbulent fluctuations have to be limited
1043!--       vertically because otherwise the turbulent inflow layer will grow
1044!--       in time.
1045          IF ( inflow_damping_height == 9999999.9 )  THEN
1046!
1047!--          Default: use the inversion height calculated by the prerun; if
1048!--          this is zero, inflow_damping_height must be explicitly
1049!--          specified.
1050             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0 )  THEN
1051                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1052             ELSE
1053                WRITE( message_string, * ) 'inflow_damping_height must be ',&
1054                     'explicitly specified because&the inversion height ', &
1055                     'calculated by the prerun is zero.'
1056                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
1057             ENDIF
1058
1059          ENDIF
1060
1061          IF ( inflow_damping_width == 9999999.9 )  THEN
1062!
1063!--          Default for the transition range: one tenth of the undamped
1064!--          layer
1065             inflow_damping_width = 0.1 * inflow_damping_height
1066
1067          ENDIF
1068
1069          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
1070
1071          DO  k = nzb, nzt+1
1072
1073             IF ( zu(k) <= inflow_damping_height )  THEN
1074                inflow_damping_factor(k) = 1.0
1075             ELSEIF ( zu(k) <= inflow_damping_height +  &
1076                               inflow_damping_width )  THEN
1077                inflow_damping_factor(k) = 1.0 -                            &
1078                                        ( zu(k) - inflow_damping_height ) / &
1079                                        inflow_damping_width
1080             ELSE
1081                inflow_damping_factor(k) = 0.0
1082             ENDIF
1083
1084          ENDDO
1085
1086       ENDIF
1087
1088!
1089!--    Inside buildings set velocities and TKE back to zero
1090       IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND.  &
1091            topography /= 'flat' )  THEN
1092!
1093!--       Inside buildings set velocities and TKE back to zero.
1094!--       Other scalars (pt, q, s, km, kh, p, sa, ...) are ignored at present,
1095!--       maybe revise later.
1096          IF ( timestep_scheme(1:5) == 'runge' )  THEN
1097             DO  i = nxlg, nxrg
1098                DO  j = nysg, nyng
1099                   u  (nzb:nzb_u_inner(j,i),j,i) = 0.0
1100                   v  (nzb:nzb_v_inner(j,i),j,i) = 0.0
1101                   w  (nzb:nzb_w_inner(j,i),j,i) = 0.0
1102                   e  (nzb:nzb_w_inner(j,i),j,i) = 0.0
1103                   u_m(nzb:nzb_u_inner(j,i),j,i) = 0.0
1104                   v_m(nzb:nzb_v_inner(j,i),j,i) = 0.0
1105                   w_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1106                   e_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1107                   tu_m(nzb:nzb_u_inner(j,i),j,i) = 0.0
1108                   tv_m(nzb:nzb_v_inner(j,i),j,i) = 0.0
1109                   tw_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1110                   te_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1111                   tpt_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1112                ENDDO
1113             ENDDO
1114          ELSE
1115             DO  i = nxlg, nxrg
1116                DO  j = nysg, nyng
1117                   u  (nzb:nzb_u_inner(j,i),j,i) = 0.0
1118                   v  (nzb:nzb_v_inner(j,i),j,i) = 0.0
1119                   w  (nzb:nzb_w_inner(j,i),j,i) = 0.0
1120                   e  (nzb:nzb_w_inner(j,i),j,i) = 0.0
1121                   u_m(nzb:nzb_u_inner(j,i),j,i) = 0.0
1122                   v_m(nzb:nzb_v_inner(j,i),j,i) = 0.0
1123                   w_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1124                   e_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1125                   u_p(nzb:nzb_u_inner(j,i),j,i) = 0.0
1126                   v_p(nzb:nzb_v_inner(j,i),j,i) = 0.0
1127                   w_p(nzb:nzb_w_inner(j,i),j,i) = 0.0
1128                   e_p(nzb:nzb_w_inner(j,i),j,i) = 0.0
1129                ENDDO
1130             ENDDO
1131          ENDIF
1132
1133       ENDIF
1134
1135!
1136!--    Calculate initial temperature field and other constants used in case
1137!--    of a sloping surface
1138       IF ( sloping_surface )  CALL init_slope
1139
1140!
1141!--    Initialize new time levels (only done in order to set boundary values
1142!--    including ghost points)
1143       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
1144       IF ( humidity  .OR.  passive_scalar )  q_p = q
1145       IF ( ocean )  sa_p = sa
1146
1147!
1148!--    Allthough tendency arrays are set in prognostic_equations, they have
1149!--    have to be predefined here because they are used (but multiplied with 0)
1150!--    there before they are set.
1151       IF ( timestep_scheme(1:5) == 'runge' )  THEN
1152          te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0
1153          IF ( humidity  .OR.  passive_scalar )  tq_m = 0.0
1154          IF ( ocean )  tsa_m = 0.0
1155       ENDIF
1156
1157    ELSE
1158!
1159!--    Actually this part of the programm should not be reached
1160       message_string = 'unknown initializing problem'
1161       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
1162    ENDIF
1163
1164
1165    IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1166!
1167!--    Initialize old timelevels needed for radiation boundary conditions
1168       IF ( outflow_l )  THEN
1169          u_m_l(:,:,:) = u(:,:,1:2)
1170          v_m_l(:,:,:) = v(:,:,0:1)
1171          w_m_l(:,:,:) = w(:,:,0:1)
1172       ENDIF
1173       IF ( outflow_r )  THEN
1174          u_m_r(:,:,:) = u(:,:,nx-1:nx)
1175          v_m_r(:,:,:) = v(:,:,nx-1:nx)
1176          w_m_r(:,:,:) = w(:,:,nx-1:nx)
1177       ENDIF
1178       IF ( outflow_s )  THEN
1179          u_m_s(:,:,:) = u(:,0:1,:)
1180          v_m_s(:,:,:) = v(:,1:2,:)
1181          w_m_s(:,:,:) = w(:,0:1,:)
1182       ENDIF
1183       IF ( outflow_n )  THEN
1184          u_m_n(:,:,:) = u(:,ny-1:ny,:)
1185          v_m_n(:,:,:) = v(:,ny-1:ny,:)
1186          w_m_n(:,:,:) = w(:,ny-1:ny,:)
1187       ENDIF
1188       
1189    ENDIF
1190
1191!
1192!-- Calculate the initial volume flow at the right and north boundary
1193    IF ( conserve_volume_flow )  THEN
1194
1195       IF ( use_prescribed_profile_data )  THEN
1196
1197          volume_flow_initial_l = 0.0
1198          volume_flow_area_l    = 0.0
1199
1200          IF ( nxr == nx )  THEN
1201             DO  j = nys, nyn
1202                DO  k = nzb_2d(j,nx)+1, nzt
1203                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1204                                              u_init(k) * dzw(k)
1205                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1206                ENDDO
1207             ENDDO
1208          ENDIF
1209         
1210          IF ( nyn == ny )  THEN
1211             DO  i = nxl, nxr
1212                DO  k = nzb_2d(ny,i)+1, nzt 
1213                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1214                                              v_init(k) * dzw(k)
1215                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1216                ENDDO
1217             ENDDO
1218          ENDIF
1219
1220#if defined( __parallel )
1221          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1222                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1223          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1224                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1225
1226#else
1227          volume_flow_initial = volume_flow_initial_l
1228          volume_flow_area    = volume_flow_area_l
1229#endif 
1230
1231       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1232
1233          volume_flow_initial_l = 0.0
1234          volume_flow_area_l    = 0.0
1235
1236          IF ( nxr == nx )  THEN
1237             DO  j = nys, nyn
1238                DO  k = nzb_2d(j,nx)+1, nzt
1239                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1240                                              hom_sum(k,1,0) * dzw(k)
1241                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1242                ENDDO
1243             ENDDO
1244          ENDIF
1245         
1246          IF ( nyn == ny )  THEN
1247             DO  i = nxl, nxr
1248                DO  k = nzb_2d(ny,i)+1, nzt 
1249                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1250                                              hom_sum(k,2,0) * dzw(k)
1251                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1252                ENDDO
1253             ENDDO
1254          ENDIF
1255
1256#if defined( __parallel )
1257          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1258                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1259          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1260                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1261
1262#else
1263          volume_flow_initial = volume_flow_initial_l
1264          volume_flow_area    = volume_flow_area_l
1265#endif 
1266
1267       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1268
1269          volume_flow_initial_l = 0.0
1270          volume_flow_area_l    = 0.0
1271
1272          IF ( nxr == nx )  THEN
1273             DO  j = nys, nyn
1274                DO  k = nzb_2d(j,nx)+1, nzt
1275                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1276                                              u(k,j,nx) * dzw(k)
1277                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1278                ENDDO
1279             ENDDO
1280          ENDIF
1281         
1282          IF ( nyn == ny )  THEN
1283             DO  i = nxl, nxr
1284                DO  k = nzb_2d(ny,i)+1, nzt 
1285                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1286                                              v(k,ny,i) * dzw(k)
1287                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1288                ENDDO
1289             ENDDO
1290          ENDIF
1291
1292#if defined( __parallel )
1293          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1294                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1295          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1296                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1297
1298#else
1299          volume_flow_initial = volume_flow_initial_l
1300          volume_flow_area    = volume_flow_area_l
1301#endif 
1302
1303       ENDIF
1304
1305!
1306!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated
1307!--    from u|v_bulk instead
1308       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
1309          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
1310          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
1311       ENDIF
1312
1313    ENDIF
1314
1315!
1316!-- Initialize quantities for special advections schemes
1317    CALL init_advec
1318
1319!
1320!-- Impose random perturbation on the horizontal velocity field and then
1321!-- remove the divergences from the velocity field at the initial stage
1322    IF ( create_disturbances .AND. &
1323         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.  &
1324         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
1325
1326       CALL disturb_field( nzb_u_inner, tend, u )
1327       CALL disturb_field( nzb_v_inner, tend, v )
1328       n_sor = nsor_ini
1329       CALL pres
1330       n_sor = nsor
1331    ENDIF
1332
1333!
1334!-- Initialization of the leaf area density
1335    IF ( plant_canopy )  THEN
1336 
1337       SELECT CASE ( TRIM( canopy_mode ) )
1338
1339          CASE( 'block' )
1340
1341             DO  i = nxlg, nxrg
1342                DO  j = nysg, nyng
1343                   lad_s(:,j,i) = lad(:)
1344                   cdc(:,j,i)   = drag_coefficient
1345                   IF ( passive_scalar )  THEN
1346                      sls(:,j,i) = leaf_surface_concentration
1347                      sec(:,j,i) = scalar_exchange_coefficient
1348                   ENDIF
1349                ENDDO
1350             ENDDO
1351
1352          CASE DEFAULT
1353
1354!
1355!--          The DEFAULT case is reached either if the parameter
1356!--          canopy mode contains a wrong character string or if the
1357!--          user has coded a special case in the user interface.
1358!--          There, the subroutine user_init_plant_canopy checks
1359!--          which of these two conditions applies.
1360             CALL user_init_plant_canopy
1361 
1362          END SELECT
1363
1364       CALL exchange_horiz( lad_s, nbgp )
1365       CALL exchange_horiz( cdc, nbgp )
1366
1367       IF ( passive_scalar )  THEN
1368          CALL exchange_horiz( sls, nbgp )
1369          CALL exchange_horiz( sec, nbgp )
1370       ENDIF
1371
1372!
1373!--    Sharp boundaries of the plant canopy in horizontal directions
1374!--    In vertical direction the interpolation is retained, as the leaf
1375!--    area density is initialised by prescribing a vertical profile
1376!--    consisting of piecewise linear segments. The upper boundary
1377!--    of the plant canopy is now defined by lad_w(pch_index,:,:) = 0.0.
1378
1379       DO  i = nxl, nxr
1380          DO  j = nys, nyn
1381             DO  k = nzb, nzt+1 
1382                IF ( lad_s(k,j,i) > 0.0 )  THEN
1383                   lad_u(k,j,i)   = lad_s(k,j,i) 
1384                   lad_u(k,j,i+1) = lad_s(k,j,i)
1385                   lad_v(k,j,i)   = lad_s(k,j,i)
1386                   lad_v(k,j+1,i) = lad_s(k,j,i)
1387                ENDIF
1388             ENDDO
1389             DO  k = nzb, nzt
1390                lad_w(k,j,i) = 0.5 * ( lad_s(k+1,j,i) + lad_s(k,j,i) )
1391             ENDDO
1392          ENDDO
1393       ENDDO
1394
1395       lad_w(pch_index,:,:) = 0.0
1396       lad_w(nzt+1,:,:)     = lad_w(nzt,:,:)
1397
1398       CALL exchange_horiz( lad_u, nbgp )
1399       CALL exchange_horiz( lad_v, nbgp )
1400       CALL exchange_horiz( lad_w, nbgp )
1401
1402!
1403!--    Initialisation of the canopy heat source distribution
1404       IF ( cthf /= 0.0 )  THEN
1405!
1406!--       Piecewise evaluation of the leaf area index by
1407!--       integration of the leaf area density
1408          lai(:,:,:) = 0.0
1409          DO  i = nxlg, nxrg
1410             DO  j = nysg, nyng
1411                DO  k = pch_index-1, 0, -1
1412                   lai(k,j,i) = lai(k+1,j,i) +                   &
1413                                ( 0.5 * ( lad_w(k+1,j,i) +       &
1414                                          lad_s(k+1,j,i) ) *     &
1415                                  ( zw(k+1) - zu(k+1) ) )  +     &
1416                                ( 0.5 * ( lad_w(k,j,i)   +       &
1417                                          lad_s(k+1,j,i) ) *     &
1418                                  ( zu(k+1) - zw(k) ) )
1419                ENDDO
1420             ENDDO
1421          ENDDO
1422
1423!
1424!--       Evaluation of the upward kinematic vertical heat flux within the
1425!--       canopy
1426          DO  i = nxlg, nxrg
1427             DO  j = nysg, nyng
1428                DO  k = 0, pch_index
1429                   canopy_heat_flux(k,j,i) = cthf *                    &
1430                                             exp( -0.6 * lai(k,j,i) )
1431                ENDDO
1432             ENDDO
1433          ENDDO
1434
1435!
1436!--       The near surface heat flux is derived from the heat flux
1437!--       distribution within the canopy
1438          shf(:,:) = canopy_heat_flux(0,:,:)
1439
1440          IF ( ASSOCIATED( shf_m ) )  shf_m = shf
1441
1442       ENDIF
1443
1444    ENDIF
1445
1446!
1447!-- If required, initialize dvrp-software
1448    IF ( dt_dvrp /= 9999999.9 )  CALL init_dvrp
1449
1450    IF ( ocean )  THEN
1451!
1452!--    Initialize quantities needed for the ocean model
1453       CALL init_ocean
1454
1455    ELSE
1456!
1457!--    Initialize quantities for handling cloud physics
1458!--    This routine must be called before lpm_init, because
1459!--    otherwise, array pt_d_t, needed in data_output_dvrp (called by
1460!--    lpm_init) is not defined.
1461       CALL init_cloud_physics
1462    ENDIF
1463
1464!
1465!-- If required, initialize particles
1466    IF ( particle_advection )  CALL lpm_init
1467
1468!
1469!-- Initialize the ws-scheme.   
1470    IF ( ws_scheme_sca .OR. ws_scheme_mom )  CALL ws_init       
1471
1472!
1473!-- Setting weighting factors for calculation of perturbation pressure
1474!-- and turbulent quantities from the RK substeps               
1475    IF ( TRIM(timestep_scheme) == 'runge-kutta-3' )  THEN      ! for RK3-method
1476
1477       weight_substep(1) = 1./6.
1478       weight_substep(2) = 3./10.
1479       weight_substep(3) = 8./15.
1480
1481       weight_pres(1)    = 1./3.
1482       weight_pres(2)    = 5./12.
1483       weight_pres(3)    = 1./4.
1484
1485    ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' )  THEN  ! for RK2-method
1486
1487       weight_substep(1) = 1./2.
1488       weight_substep(2) = 1./2.
1489         
1490       weight_pres(1)    = 1./2.
1491       weight_pres(2)    = 1./2.       
1492
1493    ELSE                                     ! for Euler- and leapfrog-method
1494
1495       weight_substep(1) = 1.0     
1496       weight_pres(1)    = 1.0                   
1497
1498    ENDIF
1499
1500!
1501!-- Initialize Rayleigh damping factors
1502    rdf    = 0.0
1503    rdf_sc = 0.0
1504    IF ( rayleigh_damping_factor /= 0.0 )  THEN
1505       IF ( .NOT. ocean )  THEN
1506          DO  k = nzb+1, nzt
1507             IF ( zu(k) >= rayleigh_damping_height )  THEN
1508                rdf(k) = rayleigh_damping_factor * &
1509                      ( SIN( pi * 0.5 * ( zu(k) - rayleigh_damping_height )    &
1510                                      / ( zu(nzt) - rayleigh_damping_height ) )&
1511                      )**2
1512             ENDIF
1513          ENDDO
1514       ELSE
1515          DO  k = nzt, nzb+1, -1
1516             IF ( zu(k) <= rayleigh_damping_height )  THEN
1517                rdf(k) = rayleigh_damping_factor * &
1518                      ( SIN( pi * 0.5 * ( rayleigh_damping_height - zu(k) )    &
1519                                      / ( rayleigh_damping_height - zu(nzb+1)))&
1520                      )**2
1521             ENDIF
1522          ENDDO
1523       ENDIF
1524    ENDIF
1525    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
1526
1527!
1528!-- Initialize the starting level and the vertical smoothing factor used for
1529!-- the external pressure gradient
1530    dp_smooth_factor = 1.0
1531    IF ( dp_external )  THEN
1532!
1533!--    Set the starting level dp_level_ind_b only if it has not been set before
1534!--    (e.g. in init_grid).
1535       IF ( dp_level_ind_b == 0 )  THEN
1536          ind_array = MINLOC( ABS( dp_level_b - zu ) )
1537          dp_level_ind_b = ind_array(1) - 1 + nzb 
1538                                        ! MINLOC uses lower array bound 1
1539       ENDIF
1540       IF ( dp_smooth )  THEN
1541          dp_smooth_factor(:dp_level_ind_b) = 0.0
1542          DO  k = dp_level_ind_b+1, nzt
1543             dp_smooth_factor(k) = 0.5 * ( 1.0 + SIN( pi * &
1544                  ( REAL( k - dp_level_ind_b ) /  &
1545                    REAL( nzt - dp_level_ind_b ) - 0.5 ) ) )
1546          ENDDO
1547       ENDIF
1548    ENDIF
1549
1550!
1551!-- Initialize damping zone for the potential temperature in case of
1552!-- non-cyclic lateral boundaries. The damping zone has the maximum value
1553!-- at the inflow boundary and decreases to zero at pt_damping_width.
1554    ptdf_x = 0.0
1555    ptdf_y = 0.0
1556    IF ( bc_lr_dirrad .OR. bc_lr_dirneu )  THEN
1557       DO i = nxl, nxr
1558          IF ( ( i * dx ) < pt_damping_width )  THEN
1559             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5 *        &
1560                         REAL( pt_damping_width - i * dx ) / (        &
1561                         REAL( pt_damping_width )            ) ) )**2 
1562          ENDIF
1563       ENDDO
1564    ELSEIF ( bc_lr_raddir .OR. bc_lr_neudir )  THEN
1565       DO i = nxl, nxr
1566          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
1567             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5 *                 &
1568                         REAL( ( i - nx ) * dx + pt_damping_width ) / (        &
1569                         REAL( pt_damping_width )                     ) ) )**2       
1570          ENDIF
1571       ENDDO 
1572    ELSEIF ( bc_ns_dirrad .OR. bc_ns_dirneu )  THEN
1573       DO j = nys, nyn
1574          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
1575             ptdf_y(j) = pt_damping_factor * ( SIN( pi * 0.5 *                 &
1576                         REAL( ( j - ny ) * dy + pt_damping_width ) / (        &
1577                         REAL( pt_damping_width )                     ) ) )**2       
1578          ENDIF
1579       ENDDO 
1580    ELSEIF ( bc_ns_raddir .OR. bc_ns_neudir )  THEN
1581       DO j = nys, nyn
1582          IF ( ( j * dy ) < pt_damping_width )  THEN
1583             ptdf_y(j) = pt_damping_factor * ( SIN( pi * 0.5 *        &
1584                         REAL( pt_damping_width - j * dy ) / (        &
1585                         REAL( pt_damping_width )            ) ) )**2       
1586          ENDIF
1587       ENDDO
1588    ENDIF
1589
1590!
1591!-- Initialize local summation arrays for routine flow_statistics.
1592!-- This is necessary because they may not yet have been initialized when they
1593!-- are called from flow_statistics (or - depending on the chosen model run -
1594!-- are never initialized)
1595    sums_divnew_l      = 0.0
1596    sums_divold_l      = 0.0
1597    sums_l_l           = 0.0
1598    sums_up_fraction_l = 0.0
1599    sums_wsts_bc_l     = 0.0
1600
1601!
1602!-- Pre-set masks for regional statistics. Default is the total model domain.
1603    rmask = 1.0
1604
1605!
1606!-- User-defined initializing actions. Check afterwards, if maximum number
1607!-- of allowed timeseries is exceeded
1608    CALL user_init
1609
1610    IF ( dots_num > dots_max )  THEN
1611       WRITE( message_string, * ) 'number of time series quantities exceeds', &
1612                                  ' its maximum of dots_max = ', dots_max,    &
1613                                  ' &Please increase dots_max in modules.f90.'
1614       CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 )   
1615    ENDIF
1616
1617!
1618!-- Input binary data file is not needed anymore. This line must be placed
1619!-- after call of user_init!
1620    CALL close_file( 13 )
1621
1622!
1623!-- Compute total sum of active mask grid points
1624!-- ngp_2dh: number of grid points of a horizontal cross section through the
1625!--          total domain
1626!-- ngp_3d:  number of grid points of the total domain
1627    ngp_2dh_outer_l   = 0
1628    ngp_2dh_outer     = 0
1629    ngp_2dh_s_inner_l = 0
1630    ngp_2dh_s_inner   = 0
1631    ngp_2dh_l         = 0
1632    ngp_2dh           = 0
1633    ngp_3d_inner_l    = 0.0
1634    ngp_3d_inner      = 0
1635    ngp_3d            = 0
1636    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
1637
1638    DO  sr = 0, statistic_regions
1639       DO  i = nxl, nxr
1640          DO  j = nys, nyn
1641             IF ( rmask(j,i,sr) == 1.0 )  THEN
1642!
1643!--             All xy-grid points
1644                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
1645!
1646!--             xy-grid points above topography
1647                DO  k = nzb_s_outer(j,i), nz + 1
1648                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1
1649                ENDDO
1650                DO  k = nzb_s_inner(j,i), nz + 1
1651                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) + 1
1652                ENDDO
1653!
1654!--             All grid points of the total domain above topography
1655                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + &
1656                                     ( nz - nzb_s_inner(j,i) + 2 )
1657             ENDIF
1658          ENDDO
1659       ENDDO
1660    ENDDO
1661
1662    sr = statistic_regions + 1
1663#if defined( __parallel )
1664    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1665    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,   &
1666                        comm2d, ierr )
1667    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1668    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr,  &
1669                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
1670    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1671    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),         &
1672                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
1673    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1674    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL, &
1675                        MPI_SUM, comm2d, ierr )
1676    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
1677#else
1678    ngp_2dh         = ngp_2dh_l
1679    ngp_2dh_outer   = ngp_2dh_outer_l
1680    ngp_2dh_s_inner = ngp_2dh_s_inner_l
1681    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
1682#endif
1683
1684    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * &
1685             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
1686
1687!
1688!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
1689!-- buoyancy, etc. A zero value will occur for cases where all grid points of
1690!-- the respective subdomain lie below the surface topography
1691    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   ) 
1692    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )),            &
1693                           ngp_3d_inner(:) )
1694    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) 
1695
1696    DEALLOCATE( ngp_2dh_l, ngp_2dh_outer_l, ngp_3d_inner_l, ngp_3d_inner_tmp )
1697
1698
1699 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.