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

Last change on this file since 996 was 996, checked in by raasch, 12 years ago

parameter use_prior_plot1d_parameters removed; little reformatting

  • Property svn:keywords set to Id
File size: 59.0 KB
Line 
1#if defined( __ibmy_special )
2@PROCESS NOOPTimize
3#endif
4 SUBROUTINE init_3d_model
5
6!------------------------------------------------------------------------------!
7! Current revisions:
8! ------------------
9! little reformatting
10!
11! Former revisions:
12! -----------------
13! $Id: init_3d_model.f90 996 2012-09-07 10:41:47Z raasch $
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             IF ( ibc_uv_b == 1 )  THEN
604!
605!--             Neumann condition
606                DO  i = nxl-1, nxr+1
607                   DO  j = nys-1, nyn+1
608                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i)
609                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i)
610                   ENDDO
611                ENDDO
612
613             ENDIF
614
615          ENDIF
616
617       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 ) &
618       THEN
619!
620!--       Use constructed initial profiles (velocity constant with height,
621!--       temperature profile with constant gradient)
622          DO  i = nxlg, nxrg
623             DO  j = nysg, nyng
624                pt(:,j,i) = pt_init
625                u(:,j,i)  = u_init
626                v(:,j,i)  = v_init
627             ENDDO
628          ENDDO
629
630!
631!--       Set initial horizontal velocities at the lowest computational grid
632!--       levels to zero in order to avoid too small time steps caused by the
633!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
634!--       in the limiting formula!). The original values are stored to be later
635!--       used for volume flow control.
636          DO  i = nxlg, nxrg
637             DO  j = nysg, nyng
638                u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0
639                v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0
640             ENDDO
641          ENDDO
642
643          IF ( humidity  .OR.  passive_scalar )  THEN
644             DO  i = nxlg, nxrg
645                DO  j = nysg, nyng
646                   q(:,j,i) = q_init
647                ENDDO
648             ENDDO
649          ENDIF
650
651          IF ( ocean )  THEN
652             DO  i = nxlg, nxrg
653                DO  j = nysg, nyng
654                   sa(:,j,i) = sa_init
655                ENDDO
656             ENDDO
657          ENDIF
658         
659          IF ( constant_diffusion )  THEN
660             km   = km_constant
661             kh   = km / prandtl_number
662             e    = 0.0
663          ELSEIF ( e_init > 0.0 )  THEN
664             DO  k = nzb+1, nzt
665                km(k,:,:) = 0.1 * l_grid(k) * SQRT( e_init )
666             ENDDO
667             km(nzb,:,:)   = km(nzb+1,:,:)
668             km(nzt+1,:,:) = km(nzt,:,:)
669             kh   = km / prandtl_number
670             e    = e_init
671          ELSE
672             IF ( .NOT. ocean )  THEN
673                kh   = 0.01   ! there must exist an initial diffusion, because
674                km   = 0.01   ! otherwise no TKE would be produced by the
675                              ! production terms, as long as not yet
676                              ! e = (u*/cm)**2 at k=nzb+1
677             ELSE
678                kh   = 0.00001
679                km   = 0.00001
680             ENDIF
681             e    = 0.0
682          ENDIF
683          rif   = 0.0
684          ts    = 0.0
685          us    = 0.0
686          usws  = 0.0
687          uswst = top_momentumflux_u
688          vsws  = 0.0
689          vswst = top_momentumflux_v
690          IF ( humidity  .OR.  passive_scalar ) qs = 0.0
691
692!
693!--       Compute initial temperature field and other constants used in case
694!--       of a sloping surface
695          IF ( sloping_surface )  CALL init_slope
696
697       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 ) &
698       THEN
699!
700!--       Initialization will completely be done by the user
701          CALL user_init_3d_model
702
703       ENDIF
704!
705!--    Bottom boundary
706       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
707          u(nzb,:,:) = 0.0
708          v(nzb,:,:) = 0.0
709       ENDIF
710
711!
712!--    Apply channel flow boundary condition
713       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
714          u(nzt+1,:,:) = 0.0
715          v(nzt+1,:,:) = 0.0
716       ENDIF
717
718!
719!--    Calculate virtual potential temperature
720       IF ( humidity ) vpt = pt * ( 1.0 + 0.61 * q )
721
722!
723!--    Store initial profiles for output purposes etc.
724       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
725       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
726       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
727          hom(nzb,1,5,:) = 0.0   
728          hom(nzb,1,6,:) = 0.0 
729       ENDIF
730       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
731       hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 )
732       hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 )
733
734       IF ( ocean )  THEN
735!
736!--       Store initial salinity profile
737          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
738       ENDIF
739
740       IF ( humidity )  THEN
741!
742!--       Store initial profile of total water content, virtual potential
743!--       temperature
744          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
745          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
746          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
747!
748!--          Store initial profile of specific humidity and potential
749!--          temperature
750             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
751             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
752          ENDIF
753       ENDIF
754
755       IF ( passive_scalar )  THEN
756!
757!--       Store initial scalar profile
758          hom(:,1,26,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
759       ENDIF
760
761!
762!--    Initialize fluxes at bottom surface
763       IF ( use_surface_fluxes )  THEN
764
765          IF ( constant_heatflux )  THEN
766!
767!--          Heat flux is prescribed
768             IF ( random_heatflux )  THEN
769                CALL disturb_heatflux
770             ELSE
771                shf = surface_heatflux
772!
773!--             Over topography surface_heatflux is replaced by wall_heatflux(0)
774                IF ( TRIM( topography ) /= 'flat' )  THEN
775                   DO  i = nxlg, nxrg
776                      DO  j = nysg, nyng
777                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
778                            shf(j,i) = wall_heatflux(0)
779                         ENDIF
780                      ENDDO
781                   ENDDO
782                ENDIF
783             ENDIF
784             IF ( ASSOCIATED( shf_m ) )  shf_m = shf
785          ENDIF
786
787!
788!--       Determine the near-surface water flux
789          IF ( humidity  .OR.  passive_scalar )  THEN
790             IF ( constant_waterflux )  THEN
791                qsws   = surface_waterflux
792!
793!--             Over topography surface_waterflux is replaced by
794!--             wall_humidityflux(0)
795                IF ( TRIM( topography ) /= 'flat' )  THEN
796                   wall_qflux = wall_humidityflux
797                   DO  i = nxlg, nxrg
798                      DO  j = nysg, nyng
799                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
800                            qsws(j,i) = wall_qflux(0)
801                         ENDIF
802                      ENDDO
803                   ENDDO
804                ENDIF
805                IF ( ASSOCIATED( qsws_m ) )  qsws_m = qsws
806             ENDIF
807          ENDIF
808
809       ENDIF
810
811!
812!--    Initialize fluxes at top surface
813!--    Currently, only the heatflux and salinity flux can be prescribed.
814!--    The latent flux is zero in this case!
815       IF ( use_top_fluxes )  THEN
816
817          IF ( constant_top_heatflux )  THEN
818!
819!--          Heat flux is prescribed
820             tswst = top_heatflux
821             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
822
823             IF ( humidity  .OR.  passive_scalar )  THEN
824                qswst = 0.0
825                IF ( ASSOCIATED( qswst_m ) )  qswst_m = qswst
826             ENDIF
827
828             IF ( ocean )  THEN
829                saswsb = bottom_salinityflux
830                saswst = top_salinityflux
831             ENDIF
832          ENDIF
833
834!
835!--       Initialization in case of a coupled model run
836          IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
837             tswst = 0.0
838             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
839          ENDIF
840
841       ENDIF
842
843!
844!--    Initialize Prandtl layer quantities
845       IF ( prandtl_layer )  THEN
846
847          z0 = roughness_length
848          z0h = z0h_factor * z0
849
850          IF ( .NOT. constant_heatflux )  THEN 
851!
852!--          Surface temperature is prescribed. Here the heat flux cannot be
853!--          simply estimated, because therefore rif, u* and theta* would have
854!--          to be computed by iteration. This is why the heat flux is assumed
855!--          to be zero before the first time step. It approaches its correct
856!--          value in the course of the first few time steps.
857             shf   = 0.0
858             IF ( ASSOCIATED( shf_m ) )  shf_m = 0.0
859          ENDIF
860
861          IF ( humidity  .OR.  passive_scalar )  THEN
862             IF ( .NOT. constant_waterflux )  THEN
863                qsws   = 0.0
864                IF ( ASSOCIATED( qsws_m ) )   qsws_m = 0.0
865             ENDIF
866          ENDIF
867
868       ENDIF
869
870
871!
872!--    For the moment, vertical velocity is zero
873       w = 0.0
874
875!
876!--    Initialize array sums (must be defined in first call of pres)
877       sums = 0.0
878
879!
880!--    In case of iterative solvers, p must get an initial value
881       IF ( psolver == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0
882
883!
884!--    Treating cloud physics, liquid water content and precipitation amount
885!--    are zero at beginning of the simulation
886       IF ( cloud_physics )  THEN
887          ql = 0.0
888          IF ( precipitation )  precipitation_amount = 0.0
889       ENDIF
890!
891!--    Impose vortex with vertical axis on the initial velocity profile
892       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
893          CALL init_rankine
894       ENDIF
895
896!
897!--    Impose temperature anomaly (advection test only)
898       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
899          CALL init_pt_anomaly
900       ENDIF
901
902!
903!--    If required, change the surface temperature at the start of the 3D run
904       IF ( pt_surface_initial_change /= 0.0 )  THEN
905          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
906       ENDIF
907
908!
909!--    If required, change the surface humidity/scalar at the start of the 3D
910!--    run
911       IF ( ( humidity .OR. passive_scalar ) .AND. &
912            q_surface_initial_change /= 0.0 )  THEN
913          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
914       ENDIF
915
916!
917!--    Initialize the random number generator (from numerical recipes)
918       CALL random_function_ini
919
920!
921!--    Initialize old and new time levels.
922       IF ( timestep_scheme(1:5) /= 'runge' )  THEN
923          e_m = e; pt_m = pt; u_m = u; v_m = v; w_m = w; kh_m = kh; km_m = km
924       ELSE
925          te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0
926       ENDIF
927       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
928
929       IF ( humidity  .OR.  passive_scalar )  THEN
930          IF ( ASSOCIATED( q_m ) )  q_m = q
931          IF ( timestep_scheme(1:5) == 'runge' )  tq_m = 0.0
932          q_p = q
933          IF ( humidity  .AND.  ASSOCIATED( vpt_m ) )  vpt_m = vpt
934       ENDIF
935
936       IF ( ocean )  THEN
937          tsa_m = 0.0
938          sa_p  = sa
939       ENDIF
940       
941
942    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.    &
943         TRIM( initializing_actions ) == 'cyclic_fill' )  &
944    THEN
945!
946!--    When reading data for cyclic fill of 3D prerun data files, read
947!--    some of the global variables from the restart file which are required
948!--    for initializing the inflow
949       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
950
951          DO  i = 0, io_blocks-1
952             IF ( i == io_group )  THEN
953                CALL read_parts_of_var_list
954                CALL close_file( 13 )
955             ENDIF
956#if defined( __parallel )
957             CALL MPI_BARRIER( comm2d, ierr )
958#endif
959          ENDDO
960
961       ENDIF
962
963!
964!--    Read binary data from restart file
965       DO  i = 0, io_blocks-1
966          IF ( i == io_group )  THEN
967             CALL read_3d_binary
968          ENDIF
969#if defined( __parallel )
970          CALL MPI_BARRIER( comm2d, ierr )
971#endif
972       ENDDO
973
974!
975!--    Initialization of the turbulence recycling method
976       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.  &
977            turbulent_inflow )  THEN
978!
979!--       First store the profiles to be used at the inflow.
980!--       These profiles are the (temporally) and horizontally averaged vertical
981!--       profiles from the prerun. Alternatively, prescribed profiles
982!--       for u,v-components can be used.
983          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,5) )
984
985          IF ( use_prescribed_profile_data )  THEN
986             mean_inflow_profiles(:,1) = u_init            ! u
987             mean_inflow_profiles(:,2) = v_init            ! v
988          ELSE
989             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
990             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
991          ENDIF
992          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
993          mean_inflow_profiles(:,5) = hom_sum(:,8,0)       ! e
994
995!
996!--       If necessary, adjust the horizontal flow field to the prescribed
997!--       profiles
998          IF ( use_prescribed_profile_data )  THEN
999             DO  i = nxlg, nxrg
1000                DO  j = nysg, nyng
1001                   DO  k = nzb, nzt+1
1002                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1003                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
1004                   ENDDO
1005                ENDDO
1006             ENDDO
1007          ENDIF
1008
1009!
1010!--       Use these mean profiles at the inflow (provided that Dirichlet
1011!--       conditions are used)
1012          IF ( inflow_l )  THEN
1013             DO  j = nysg, nyng
1014                DO  k = nzb, nzt+1
1015                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1016                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
1017                   w(k,j,nxlg:-1)  = 0.0
1018                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
1019                   e(k,j,nxlg:-1)  = mean_inflow_profiles(k,5)
1020                ENDDO
1021             ENDDO
1022          ENDIF
1023
1024!
1025!--       Calculate the damping factors to be used at the inflow. For a
1026!--       turbulent inflow the turbulent fluctuations have to be limited
1027!--       vertically because otherwise the turbulent inflow layer will grow
1028!--       in time.
1029          IF ( inflow_damping_height == 9999999.9 )  THEN
1030!
1031!--          Default: use the inversion height calculated by the prerun; if
1032!--          this is zero, inflow_damping_height must be explicitly
1033!--          specified.
1034             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0 )  THEN
1035                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1036             ELSE
1037                WRITE( message_string, * ) 'inflow_damping_height must be ',&
1038                     'explicitly specified because&the inversion height ', &
1039                     'calculated by the prerun is zero.'
1040                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
1041             ENDIF
1042
1043          ENDIF
1044
1045          IF ( inflow_damping_width == 9999999.9 )  THEN
1046!
1047!--          Default for the transition range: one tenth of the undamped
1048!--          layer
1049             inflow_damping_width = 0.1 * inflow_damping_height
1050
1051          ENDIF
1052
1053          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
1054
1055          DO  k = nzb, nzt+1
1056
1057             IF ( zu(k) <= inflow_damping_height )  THEN
1058                inflow_damping_factor(k) = 1.0
1059             ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) )  THEN
1060                inflow_damping_factor(k) = 1.0 -                               &
1061                                           ( zu(k) - inflow_damping_height ) / &
1062                                           inflow_damping_width
1063             ELSE
1064                inflow_damping_factor(k) = 0.0
1065             ENDIF
1066
1067          ENDDO
1068
1069       ENDIF
1070
1071!
1072!--    Inside buildings set velocities and TKE back to zero
1073       IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND.  &
1074            topography /= 'flat' )  THEN
1075!
1076!--       Inside buildings set velocities and TKE back to zero.
1077!--       Other scalars (pt, q, s, km, kh, p, sa, ...) are ignored at present,
1078!--       maybe revise later.
1079          IF ( timestep_scheme(1:5) == 'runge' )  THEN
1080             DO  i = nxlg, nxrg
1081                DO  j = nysg, nyng
1082                   u  (nzb:nzb_u_inner(j,i),j,i)   = 0.0
1083                   v  (nzb:nzb_v_inner(j,i),j,i)   = 0.0
1084                   w  (nzb:nzb_w_inner(j,i),j,i)   = 0.0
1085                   e  (nzb:nzb_w_inner(j,i),j,i)   = 0.0
1086                   u_m(nzb:nzb_u_inner(j,i),j,i)   = 0.0
1087                   v_m(nzb:nzb_v_inner(j,i),j,i)   = 0.0
1088                   w_m(nzb:nzb_w_inner(j,i),j,i)   = 0.0
1089                   e_m(nzb:nzb_w_inner(j,i),j,i)   = 0.0
1090                   tu_m(nzb:nzb_u_inner(j,i),j,i)  = 0.0
1091                   tv_m(nzb:nzb_v_inner(j,i),j,i)  = 0.0
1092                   tw_m(nzb:nzb_w_inner(j,i),j,i)  = 0.0
1093                   te_m(nzb:nzb_w_inner(j,i),j,i)  = 0.0
1094                   tpt_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1095                ENDDO
1096             ENDDO
1097          ELSE
1098             DO  i = nxlg, nxrg
1099                DO  j = nysg, nyng
1100                   u  (nzb:nzb_u_inner(j,i),j,i) = 0.0
1101                   v  (nzb:nzb_v_inner(j,i),j,i) = 0.0
1102                   w  (nzb:nzb_w_inner(j,i),j,i) = 0.0
1103                   e  (nzb:nzb_w_inner(j,i),j,i) = 0.0
1104                   u_m(nzb:nzb_u_inner(j,i),j,i) = 0.0
1105                   v_m(nzb:nzb_v_inner(j,i),j,i) = 0.0
1106                   w_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1107                   e_m(nzb:nzb_w_inner(j,i),j,i) = 0.0
1108                   u_p(nzb:nzb_u_inner(j,i),j,i) = 0.0
1109                   v_p(nzb:nzb_v_inner(j,i),j,i) = 0.0
1110                   w_p(nzb:nzb_w_inner(j,i),j,i) = 0.0
1111                   e_p(nzb:nzb_w_inner(j,i),j,i) = 0.0
1112                ENDDO
1113             ENDDO
1114          ENDIF
1115
1116       ENDIF
1117
1118!
1119!--    Calculate initial temperature field and other constants used in case
1120!--    of a sloping surface
1121       IF ( sloping_surface )  CALL init_slope
1122
1123!
1124!--    Initialize new time levels (only done in order to set boundary values
1125!--    including ghost points)
1126       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
1127       IF ( humidity  .OR.  passive_scalar )  q_p = q
1128       IF ( ocean )  sa_p = sa
1129
1130!
1131!--    Allthough tendency arrays are set in prognostic_equations, they have
1132!--    have to be predefined here because they are used (but multiplied with 0)
1133!--    there before they are set.
1134       IF ( timestep_scheme(1:5) == 'runge' )  THEN
1135          te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0
1136          IF ( humidity  .OR.  passive_scalar )  tq_m = 0.0
1137          IF ( ocean )  tsa_m = 0.0
1138       ENDIF
1139
1140    ELSE
1141!
1142!--    Actually this part of the programm should not be reached
1143       message_string = 'unknown initializing problem'
1144       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
1145    ENDIF
1146
1147
1148    IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1149!
1150!--    Initialize old timelevels needed for radiation boundary conditions
1151       IF ( outflow_l )  THEN
1152          u_m_l(:,:,:) = u(:,:,1:2)
1153          v_m_l(:,:,:) = v(:,:,0:1)
1154          w_m_l(:,:,:) = w(:,:,0:1)
1155       ENDIF
1156       IF ( outflow_r )  THEN
1157          u_m_r(:,:,:) = u(:,:,nx-1:nx)
1158          v_m_r(:,:,:) = v(:,:,nx-1:nx)
1159          w_m_r(:,:,:) = w(:,:,nx-1:nx)
1160       ENDIF
1161       IF ( outflow_s )  THEN
1162          u_m_s(:,:,:) = u(:,0:1,:)
1163          v_m_s(:,:,:) = v(:,1:2,:)
1164          w_m_s(:,:,:) = w(:,0:1,:)
1165       ENDIF
1166       IF ( outflow_n )  THEN
1167          u_m_n(:,:,:) = u(:,ny-1:ny,:)
1168          v_m_n(:,:,:) = v(:,ny-1:ny,:)
1169          w_m_n(:,:,:) = w(:,ny-1:ny,:)
1170       ENDIF
1171       
1172    ENDIF
1173
1174!
1175!-- Calculate the initial volume flow at the right and north boundary
1176    IF ( conserve_volume_flow )  THEN
1177
1178       IF ( use_prescribed_profile_data )  THEN
1179
1180          volume_flow_initial_l = 0.0
1181          volume_flow_area_l    = 0.0
1182
1183          IF ( nxr == nx )  THEN
1184             DO  j = nys, nyn
1185                DO  k = nzb_2d(j,nx)+1, nzt
1186                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1187                                              u_init(k) * dzw(k)
1188                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1189                ENDDO
1190             ENDDO
1191          ENDIF
1192         
1193          IF ( nyn == ny )  THEN
1194             DO  i = nxl, nxr
1195                DO  k = nzb_2d(ny,i)+1, nzt 
1196                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1197                                              v_init(k) * dzw(k)
1198                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1199                ENDDO
1200             ENDDO
1201          ENDIF
1202
1203#if defined( __parallel )
1204          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1205                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1206          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1207                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1208
1209#else
1210          volume_flow_initial = volume_flow_initial_l
1211          volume_flow_area    = volume_flow_area_l
1212#endif 
1213
1214       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1215
1216          volume_flow_initial_l = 0.0
1217          volume_flow_area_l    = 0.0
1218
1219          IF ( nxr == nx )  THEN
1220             DO  j = nys, nyn
1221                DO  k = nzb_2d(j,nx)+1, nzt
1222                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1223                                              hom_sum(k,1,0) * dzw(k)
1224                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1225                ENDDO
1226             ENDDO
1227          ENDIF
1228         
1229          IF ( nyn == ny )  THEN
1230             DO  i = nxl, nxr
1231                DO  k = nzb_2d(ny,i)+1, nzt 
1232                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1233                                              hom_sum(k,2,0) * dzw(k)
1234                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1235                ENDDO
1236             ENDDO
1237          ENDIF
1238
1239#if defined( __parallel )
1240          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1241                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1242          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1243                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1244
1245#else
1246          volume_flow_initial = volume_flow_initial_l
1247          volume_flow_area    = volume_flow_area_l
1248#endif 
1249
1250       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1251
1252          volume_flow_initial_l = 0.0
1253          volume_flow_area_l    = 0.0
1254
1255          IF ( nxr == nx )  THEN
1256             DO  j = nys, nyn
1257                DO  k = nzb_2d(j,nx)+1, nzt
1258                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1259                                              u(k,j,nx) * dzw(k)
1260                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1261                ENDDO
1262             ENDDO
1263          ENDIF
1264         
1265          IF ( nyn == ny )  THEN
1266             DO  i = nxl, nxr
1267                DO  k = nzb_2d(ny,i)+1, nzt 
1268                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1269                                              v(k,ny,i) * dzw(k)
1270                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1271                ENDDO
1272             ENDDO
1273          ENDIF
1274
1275#if defined( __parallel )
1276          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1277                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1278          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1279                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1280
1281#else
1282          volume_flow_initial = volume_flow_initial_l
1283          volume_flow_area    = volume_flow_area_l
1284#endif 
1285
1286       ENDIF
1287
1288!
1289!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated
1290!--    from u|v_bulk instead
1291       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
1292          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
1293          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
1294       ENDIF
1295
1296    ENDIF
1297
1298!
1299!-- Initialize quantities for special advections schemes
1300    CALL init_advec
1301
1302!
1303!-- Impose random perturbation on the horizontal velocity field and then
1304!-- remove the divergences from the velocity field at the initial stage
1305    IF ( create_disturbances .AND. &
1306         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.  &
1307         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
1308
1309       CALL disturb_field( nzb_u_inner, tend, u )
1310       CALL disturb_field( nzb_v_inner, tend, v )
1311       n_sor = nsor_ini
1312       CALL pres
1313       n_sor = nsor
1314    ENDIF
1315
1316!
1317!-- Initialization of the leaf area density
1318    IF ( plant_canopy )  THEN
1319 
1320       SELECT CASE ( TRIM( canopy_mode ) )
1321
1322          CASE( 'block' )
1323
1324             DO  i = nxlg, nxrg
1325                DO  j = nysg, nyng
1326                   lad_s(:,j,i) = lad(:)
1327                   cdc(:,j,i)   = drag_coefficient
1328                   IF ( passive_scalar )  THEN
1329                      sls(:,j,i) = leaf_surface_concentration
1330                      sec(:,j,i) = scalar_exchange_coefficient
1331                   ENDIF
1332                ENDDO
1333             ENDDO
1334
1335          CASE DEFAULT
1336
1337!
1338!--          The DEFAULT case is reached either if the parameter
1339!--          canopy mode contains a wrong character string or if the
1340!--          user has coded a special case in the user interface.
1341!--          There, the subroutine user_init_plant_canopy checks
1342!--          which of these two conditions applies.
1343             CALL user_init_plant_canopy
1344 
1345          END SELECT
1346
1347       CALL exchange_horiz( lad_s, nbgp )
1348       CALL exchange_horiz( cdc, nbgp )
1349
1350       IF ( passive_scalar )  THEN
1351          CALL exchange_horiz( sls, nbgp )
1352          CALL exchange_horiz( sec, nbgp )
1353       ENDIF
1354
1355!
1356!--    Sharp boundaries of the plant canopy in horizontal directions
1357!--    In vertical direction the interpolation is retained, as the leaf
1358!--    area density is initialised by prescribing a vertical profile
1359!--    consisting of piecewise linear segments. The upper boundary
1360!--    of the plant canopy is now defined by lad_w(pch_index,:,:) = 0.0.
1361
1362       DO  i = nxl, nxr
1363          DO  j = nys, nyn
1364             DO  k = nzb, nzt+1 
1365                IF ( lad_s(k,j,i) > 0.0 )  THEN
1366                   lad_u(k,j,i)   = lad_s(k,j,i) 
1367                   lad_u(k,j,i+1) = lad_s(k,j,i)
1368                   lad_v(k,j,i)   = lad_s(k,j,i)
1369                   lad_v(k,j+1,i) = lad_s(k,j,i)
1370                ENDIF
1371             ENDDO
1372             DO  k = nzb, nzt
1373                lad_w(k,j,i) = 0.5 * ( lad_s(k+1,j,i) + lad_s(k,j,i) )
1374             ENDDO
1375          ENDDO
1376       ENDDO
1377
1378       lad_w(pch_index,:,:) = 0.0
1379       lad_w(nzt+1,:,:)     = lad_w(nzt,:,:)
1380
1381       CALL exchange_horiz( lad_u, nbgp )
1382       CALL exchange_horiz( lad_v, nbgp )
1383       CALL exchange_horiz( lad_w, nbgp )
1384
1385!
1386!--    Initialisation of the canopy heat source distribution
1387       IF ( cthf /= 0.0 )  THEN
1388!
1389!--       Piecewise evaluation of the leaf area index by
1390!--       integration of the leaf area density
1391          lai(:,:,:) = 0.0
1392          DO  i = nxlg, nxrg
1393             DO  j = nysg, nyng
1394                DO  k = pch_index-1, 0, -1
1395                   lai(k,j,i) = lai(k+1,j,i) +                   &
1396                                ( 0.5 * ( lad_w(k+1,j,i) +       &
1397                                          lad_s(k+1,j,i) ) *     &
1398                                  ( zw(k+1) - zu(k+1) ) )  +     &
1399                                ( 0.5 * ( lad_w(k,j,i)   +       &
1400                                          lad_s(k+1,j,i) ) *     &
1401                                  ( zu(k+1) - zw(k) ) )
1402                ENDDO
1403             ENDDO
1404          ENDDO
1405
1406!
1407!--       Evaluation of the upward kinematic vertical heat flux within the
1408!--       canopy
1409          DO  i = nxlg, nxrg
1410             DO  j = nysg, nyng
1411                DO  k = 0, pch_index
1412                   canopy_heat_flux(k,j,i) = cthf *                    &
1413                                             exp( -0.6 * lai(k,j,i) )
1414                ENDDO
1415             ENDDO
1416          ENDDO
1417
1418!
1419!--       The near surface heat flux is derived from the heat flux
1420!--       distribution within the canopy
1421          shf(:,:) = canopy_heat_flux(0,:,:)
1422
1423          IF ( ASSOCIATED( shf_m ) )  shf_m = shf
1424
1425       ENDIF
1426
1427    ENDIF
1428
1429!
1430!-- If required, initialize dvrp-software
1431    IF ( dt_dvrp /= 9999999.9 )  CALL init_dvrp
1432
1433    IF ( ocean )  THEN
1434!
1435!--    Initialize quantities needed for the ocean model
1436       CALL init_ocean
1437
1438    ELSE
1439!
1440!--    Initialize quantities for handling cloud physics
1441!--    This routine must be called before lpm_init, because
1442!--    otherwise, array pt_d_t, needed in data_output_dvrp (called by
1443!--    lpm_init) is not defined.
1444       CALL init_cloud_physics
1445    ENDIF
1446
1447!
1448!-- If required, initialize particles
1449    IF ( particle_advection )  CALL lpm_init
1450
1451!
1452!-- Initialize the ws-scheme.   
1453    IF ( ws_scheme_sca .OR. ws_scheme_mom )  CALL ws_init       
1454
1455!
1456!-- Setting weighting factors for calculation of perturbation pressure
1457!-- and turbulent quantities from the RK substeps               
1458    IF ( TRIM(timestep_scheme) == 'runge-kutta-3' )  THEN      ! for RK3-method
1459
1460       weight_substep(1) = 1./6.
1461       weight_substep(2) = 3./10.
1462       weight_substep(3) = 8./15.
1463
1464       weight_pres(1)    = 1./3.
1465       weight_pres(2)    = 5./12.
1466       weight_pres(3)    = 1./4.
1467
1468    ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' )  THEN  ! for RK2-method
1469
1470       weight_substep(1) = 1./2.
1471       weight_substep(2) = 1./2.
1472         
1473       weight_pres(1)    = 1./2.
1474       weight_pres(2)    = 1./2.       
1475
1476    ELSE                                     ! for Euler- and leapfrog-method
1477
1478       weight_substep(1) = 1.0     
1479       weight_pres(1)    = 1.0                   
1480
1481    ENDIF
1482
1483!
1484!-- Initialize Rayleigh damping factors
1485    rdf    = 0.0
1486    rdf_sc = 0.0
1487    IF ( rayleigh_damping_factor /= 0.0 )  THEN
1488       IF ( .NOT. ocean )  THEN
1489          DO  k = nzb+1, nzt
1490             IF ( zu(k) >= rayleigh_damping_height )  THEN
1491                rdf(k) = rayleigh_damping_factor * &
1492                      ( SIN( pi * 0.5 * ( zu(k) - rayleigh_damping_height )    &
1493                                      / ( zu(nzt) - rayleigh_damping_height ) )&
1494                      )**2
1495             ENDIF
1496          ENDDO
1497       ELSE
1498          DO  k = nzt, nzb+1, -1
1499             IF ( zu(k) <= rayleigh_damping_height )  THEN
1500                rdf(k) = rayleigh_damping_factor * &
1501                      ( SIN( pi * 0.5 * ( rayleigh_damping_height - zu(k) )    &
1502                                      / ( rayleigh_damping_height - zu(nzb+1)))&
1503                      )**2
1504             ENDIF
1505          ENDDO
1506       ENDIF
1507    ENDIF
1508    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
1509
1510!
1511!-- Initialize the starting level and the vertical smoothing factor used for
1512!-- the external pressure gradient
1513    dp_smooth_factor = 1.0
1514    IF ( dp_external )  THEN
1515!
1516!--    Set the starting level dp_level_ind_b only if it has not been set before
1517!--    (e.g. in init_grid).
1518       IF ( dp_level_ind_b == 0 )  THEN
1519          ind_array = MINLOC( ABS( dp_level_b - zu ) )
1520          dp_level_ind_b = ind_array(1) - 1 + nzb 
1521                                        ! MINLOC uses lower array bound 1
1522       ENDIF
1523       IF ( dp_smooth )  THEN
1524          dp_smooth_factor(:dp_level_ind_b) = 0.0
1525          DO  k = dp_level_ind_b+1, nzt
1526             dp_smooth_factor(k) = 0.5 * ( 1.0 + SIN( pi * &
1527                  ( REAL( k - dp_level_ind_b ) /  &
1528                    REAL( nzt - dp_level_ind_b ) - 0.5 ) ) )
1529          ENDDO
1530       ENDIF
1531    ENDIF
1532
1533!
1534!-- Initialize damping zone for the potential temperature in case of
1535!-- non-cyclic lateral boundaries. The damping zone has the maximum value
1536!-- at the inflow boundary and decreases to zero at pt_damping_width.
1537    ptdf_x = 0.0
1538    ptdf_y = 0.0
1539    IF ( bc_lr_dirrad  .OR.  bc_lr_dirneu )  THEN
1540       DO  i = nxl, nxr
1541          IF ( ( i * dx ) < pt_damping_width )  THEN
1542             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5 *        &
1543                         REAL( pt_damping_width - i * dx ) / (        &
1544                         REAL( pt_damping_width )            ) ) )**2 
1545          ENDIF
1546       ENDDO
1547    ELSEIF ( bc_lr_raddir  .OR.  bc_lr_neudir )  THEN
1548       DO  i = nxl, nxr
1549          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
1550             ptdf_x(i) = pt_damping_factor *                                      &
1551                         SIN( pi * 0.5 * ( ( i - nx ) * dx + pt_damping_width ) / &
1552                                         REAL( pt_damping_width ) )**2
1553          ENDIF
1554       ENDDO 
1555    ELSEIF ( bc_ns_dirrad  .OR.  bc_ns_dirneu )  THEN
1556       DO  j = nys, nyn
1557          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
1558             ptdf_y(j) = pt_damping_factor *                                      &
1559                         SIN( pi * 0.5 * ( ( j - ny ) * dy + pt_damping_width ) / &
1560                                         REAL( pt_damping_width ) )**2
1561          ENDIF
1562       ENDDO 
1563    ELSEIF ( bc_ns_raddir  .OR.  bc_ns_neudir )  THEN
1564       DO  j = nys, nyn
1565          IF ( ( j * dy ) < pt_damping_width )  THEN
1566             ptdf_y(j) = pt_damping_factor *                             &
1567                         SIN( pi * 0.5 * ( pt_damping_width - j * dy ) / &
1568                                         REAL( pt_damping_width ) )**2
1569          ENDIF
1570       ENDDO
1571    ENDIF
1572
1573!
1574!-- Initialize local summation arrays for routine flow_statistics.
1575!-- This is necessary because they may not yet have been initialized when they
1576!-- are called from flow_statistics (or - depending on the chosen model run -
1577!-- are never initialized)
1578    sums_divnew_l      = 0.0
1579    sums_divold_l      = 0.0
1580    sums_l_l           = 0.0
1581    sums_up_fraction_l = 0.0
1582    sums_wsts_bc_l     = 0.0
1583
1584!
1585!-- Pre-set masks for regional statistics. Default is the total model domain.
1586    rmask = 1.0
1587
1588!
1589!-- User-defined initializing actions. Check afterwards, if maximum number
1590!-- of allowed timeseries is exceeded
1591    CALL user_init
1592
1593    IF ( dots_num > dots_max )  THEN
1594       WRITE( message_string, * ) 'number of time series quantities exceeds', &
1595                                  ' its maximum of dots_max = ', dots_max,    &
1596                                  ' &Please increase dots_max in modules.f90.'
1597       CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 )   
1598    ENDIF
1599
1600!
1601!-- Input binary data file is not needed anymore. This line must be placed
1602!-- after call of user_init!
1603    CALL close_file( 13 )
1604
1605!
1606!-- Compute total sum of active mask grid points
1607!-- ngp_2dh: number of grid points of a horizontal cross section through the
1608!--          total domain
1609!-- ngp_3d:  number of grid points of the total domain
1610    ngp_2dh_outer_l   = 0
1611    ngp_2dh_outer     = 0
1612    ngp_2dh_s_inner_l = 0
1613    ngp_2dh_s_inner   = 0
1614    ngp_2dh_l         = 0
1615    ngp_2dh           = 0
1616    ngp_3d_inner_l    = 0.0
1617    ngp_3d_inner      = 0
1618    ngp_3d            = 0
1619    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
1620
1621    DO  sr = 0, statistic_regions
1622       DO  i = nxl, nxr
1623          DO  j = nys, nyn
1624             IF ( rmask(j,i,sr) == 1.0 )  THEN
1625!
1626!--             All xy-grid points
1627                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
1628!
1629!--             xy-grid points above topography
1630                DO  k = nzb_s_outer(j,i), nz + 1
1631                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1
1632                ENDDO
1633                DO  k = nzb_s_inner(j,i), nz + 1
1634                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) + 1
1635                ENDDO
1636!
1637!--             All grid points of the total domain above topography
1638                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + &
1639                                     ( nz - nzb_s_inner(j,i) + 2 )
1640             ENDIF
1641          ENDDO
1642       ENDDO
1643    ENDDO
1644
1645    sr = statistic_regions + 1
1646#if defined( __parallel )
1647    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1648    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,   &
1649                        comm2d, ierr )
1650    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1651    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr,  &
1652                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
1653    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1654    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),         &
1655                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
1656    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1657    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL, &
1658                        MPI_SUM, comm2d, ierr )
1659    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
1660#else
1661    ngp_2dh         = ngp_2dh_l
1662    ngp_2dh_outer   = ngp_2dh_outer_l
1663    ngp_2dh_s_inner = ngp_2dh_s_inner_l
1664    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
1665#endif
1666
1667    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * &
1668             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
1669
1670!
1671!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
1672!-- buoyancy, etc. A zero value will occur for cases where all grid points of
1673!-- the respective subdomain lie below the surface topography
1674    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   ) 
1675    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )),            &
1676                           ngp_3d_inner(:) )
1677    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) 
1678
1679    DEALLOCATE( ngp_2dh_l, ngp_2dh_outer_l, ngp_3d_inner_l, ngp_3d_inner_tmp )
1680
1681
1682 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.