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

Last change on this file since 787 was 787, checked in by heinze, 12 years ago

bugfix: call init_advec in every case - not only for inital runs

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