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

Last change on this file since 974 was 850, checked in by raasch, 13 years ago

last commit documented

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