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

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

last commit documented

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