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

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

last commit documented

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