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

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

last commit documented

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