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

Last change on this file since 826 was 826, checked in by raasch, 12 years ago

last commit documented

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