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

Last change on this file since 134 was 132, checked in by letzel, 17 years ago

Vertical profiles now based on nzb_s_inner; they are divided by
ngp_2dh_s_inner (scalars, procucts of scalars and velocity components) and
ngp_2dh (staggered velocity components and their products), respectively.

Allow new case bc_uv_t = 'dirichlet_0' for channel flow.

  • Property svn:keywords set to Id
File size: 39.0 KB
Line 
1#if defined( __ibmy_special )
2@PROCESS NOOPTimize
3#endif
4 SUBROUTINE init_3d_model
5
6!------------------------------------------------------------------------------!
7! Actual revisions:
8! -----------------
9! New counter ngp_2dh_s_inner.
10! Allow new case bc_uv_t = 'dirichlet_0' for channel flow.
11! Corrected calculation of initial volume flow for 'set_1d-model_profiles' and
12! 'set_constant_profiles' in case of buildings in the reference cross-sections.
13!
14! Former revisions:
15! -----------------
16! $Id: init_3d_model.f90 132 2007-11-20 09:46:11Z letzel $
17!
18! 108 2007-08-24 15:10:38Z letzel
19! Flux initialization in case of coupled runs, +momentum fluxes at top boundary,
20! +arrays for phase speed c_u, c_v, c_w, indices for u|v|w_m_l|r changed
21! +qswst_remote in case of atmosphere model with humidity coupled to ocean
22! Rayleigh damping for ocean, optionally calculate km and kh from initial
23! TKE e_init
24!
25! 97 2007-06-21 08:23:15Z raasch
26! Initialization of salinity, call of init_ocean
27!
28! 87 2007-05-22 15:46:47Z raasch
29! var_hom and var_sum renamed pr_palm
30!
31! 75 2007-03-22 09:54:05Z raasch
32! Arrays for radiation boundary conditions are allocated (u_m_l, u_m_r, etc.),
33! bugfix for cases with the outflow damping layer extending over more than one
34! subdomain, moisture renamed humidity,
35! new initializing action "by_user" calls user_init_3d_model,
36! precipitation_amount/rate, ts_value are allocated, +module netcdf_control,
37! initial velocities at nzb+1 are regarded for volume
38! flow control in case they have been set zero before (to avoid small timesteps)
39! -uvmean_outflow, uxrp, vynp eliminated
40!
41! 19 2007-02-23 04:53:48Z raasch
42! +handling of top fluxes
43!
44! RCS Log replace by Id keyword, revision history cleaned up
45!
46! Revision 1.49  2006/08/22 15:59:07  raasch
47! No optimization of this file on the ibmy (Yonsei Univ.)
48!
49! Revision 1.1  1998/03/09 16:22:22  raasch
50! Initial revision
51!
52!
53! Description:
54! ------------
55! Allocation of arrays and initialization of the 3D model via
56! a) pre-run the 1D model
57! or
58! b) pre-set constant linear profiles
59! or
60! c) read values of a previous run
61!------------------------------------------------------------------------------!
62
63    USE arrays_3d
64    USE averaging
65    USE cloud_parameters
66    USE constants
67    USE control_parameters
68    USE cpulog
69    USE indices
70    USE interfaces
71    USE model_1d
72    USE netcdf_control
73    USE particle_attributes
74    USE pegrid
75    USE profil_parameter
76    USE random_function_mod
77    USE statistics
78
79    IMPLICIT NONE
80
81    INTEGER ::  i, j, k, sr
82
83    INTEGER, DIMENSION(:), ALLOCATABLE ::  ngp_2dh_l, ngp_3d_inner_l
84
85    INTEGER, DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l,  &
86         ngp_2dh_s_inner_l
87
88    REAL, DIMENSION(1:2) ::  volume_flow_area_l, volume_flow_initial_l
89
90
91!
92!-- Allocate arrays
93    ALLOCATE( ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions), &
94              ngp_3d(0:statistic_regions),                                  &
95              ngp_3d_inner(0:statistic_regions),                            &
96              ngp_3d_inner_l(0:statistic_regions),                          &
97              sums_divnew_l(0:statistic_regions),                           &
98              sums_divold_l(0:statistic_regions) )
99    ALLOCATE( rdf(nzb+1:nzt) )
100    ALLOCATE( hom_sum(nzb:nzt+1,pr_palm+max_pr_user,0:statistic_regions),   &
101              ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                 &
102              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),               &
103              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),               &
104              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),             &
105              rmask(nys-1:nyn+1,nxl-1:nxr+1,0:statistic_regions),           &
106              sums(nzb:nzt+1,pr_palm+max_pr_user),                          &
107              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),   &
108              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1), &
109              sums_up_fraction_l(10,3,0:statistic_regions),                 &
110              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                &
111              ts_value(var_ts,0:statistic_regions) )
112    ALLOCATE( km_damp_x(nxl-1:nxr+1), km_damp_y(nys-1:nyn+1) )
113
114    ALLOCATE( rif_1(nys-1:nyn+1,nxl-1:nxr+1), shf_1(nys-1:nyn+1,nxl-1:nxr+1), &
115              ts(nys-1:nyn+1,nxl-1:nxr+1), tswst_1(nys-1:nyn+1,nxl-1:nxr+1),  &
116              us(nys-1:nyn+1,nxl-1:nxr+1), usws_1(nys-1:nyn+1,nxl-1:nxr+1),   &
117              uswst_1(nys-1:nyn+1,nxl-1:nxr+1),                               &
118              vsws_1(nys-1:nyn+1,nxl-1:nxr+1),                                &
119              vswst_1(nys-1:nyn+1,nxl-1:nxr+1), z0(nys-1:nyn+1,nxl-1:nxr+1) )
120
121    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
122!
123!--    Leapfrog scheme needs two timelevels of diffusion quantities
124       ALLOCATE( rif_2(nys-1:nyn+1,nxl-1:nxr+1),   &
125                 shf_2(nys-1:nyn+1,nxl-1:nxr+1),   &
126                 tswst_2(nys-1:nyn+1,nxl-1:nxr+1), &
127                 usws_2(nys-1:nyn+1,nxl-1:nxr+1),  &
128                 uswst_2(nys-1:nyn+1,nxl-1:nxr+1), &
129                 vswst_2(nys-1:nyn+1,nxl-1:nxr+1), &
130                 vsws_2(nys-1:nyn+1,nxl-1:nxr+1) )
131    ENDIF
132
133    ALLOCATE( d(nzb+1:nzta,nys:nyna,nxl:nxra),         &
134              e_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
135              e_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
136              e_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
137              kh_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
138              km_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
139              p(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),    &
140              pt_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
141              pt_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
142              pt_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
143              tend(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
144              u_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
145              u_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
146              u_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
147              v_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
148              v_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
149              v_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
150              w_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
151              w_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
152              w_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
153
154    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
155       ALLOCATE( kh_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
156                 km_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
157    ENDIF
158
159    IF ( humidity  .OR.  passive_scalar ) THEN
160!
161!--    2D-humidity/scalar arrays
162       ALLOCATE ( qs(nys-1:nyn+1,nxl-1:nxr+1),     &
163                  qsws_1(nys-1:nyn+1,nxl-1:nxr+1), &
164                  qswst_1(nys-1:nyn+1,nxl-1:nxr+1) )
165
166       IF ( timestep_scheme(1:5) /= 'runge' )  THEN
167          ALLOCATE( qsws_2(nys-1:nyn+1,nxl-1:nxr+1), &
168                    qswst_2(nys-1:nyn+1,nxl-1:nxr+1) )
169       ENDIF
170!
171!--    3D-humidity/scalar arrays
172       ALLOCATE( q_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
173                 q_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
174                 q_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
175
176!
177!--    3D-arrays needed for humidity only
178       IF ( humidity )  THEN
179          ALLOCATE( vpt_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
180
181          IF ( timestep_scheme(1:5) /= 'runge' )  THEN
182             ALLOCATE( vpt_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
183          ENDIF
184
185          IF ( cloud_physics ) THEN
186!
187!--          Liquid water content
188             ALLOCATE ( ql_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
189!
190!--          Precipitation amount and rate (only needed if output is switched)
191             ALLOCATE( precipitation_amount(nys-1:nyn+1,nxl-1:nxr+1), &
192                       precipitation_rate(nys-1:nyn+1,nxl-1:nxr+1) )
193          ENDIF
194
195          IF ( cloud_droplets )  THEN
196!
197!--          Liquid water content, change in liquid water content,
198!--          real volume of particles (with weighting), volume of particles
199             ALLOCATE ( ql_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
200                        ql_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
201                        ql_v(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
202                        ql_vp(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
203          ENDIF
204
205       ENDIF
206
207    ENDIF
208
209    IF ( ocean )  THEN
210       ALLOCATE( saswsb_1(nys-1:nyn+1,nxl-1:nxr+1), &
211                 saswst_1(nys-1:nyn+1,nxl-1:nxr+1) )
212       ALLOCATE( rho_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1), &
213                 sa_1(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
214                 sa_2(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1),  &
215                 sa_3(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
216       rho => rho_1  ! routine calc_mean_profile requires density to be a
217                     ! pointer
218       IF ( humidity_remote )  THEN
219          ALLOCATE( qswst_remote(nys-1:nyn+1,nxl-1:nxr+1) )
220          qswst_remote = 0.0
221       ENDIF
222    ENDIF
223
224!
225!-- 3D-array for storing the dissipation, needed for calculating the sgs
226!-- particle velocities
227    IF ( use_sgs_for_particles )  THEN
228       ALLOCATE ( diss(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
229    ENDIF
230
231    IF ( dt_dosp /= 9999999.9 )  THEN
232       ALLOCATE( spectrum_x( 1:nx/2, 1:10, 1:10 ), &
233                 spectrum_y( 1:ny/2, 1:10, 1:10 ) )
234    ENDIF
235
236!
237!-- 4D-array for storing the Rif-values at vertical walls
238    IF ( topography /= 'flat' )  THEN
239       ALLOCATE( rif_wall(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1,1:4) )
240       rif_wall = 0.0
241    ENDIF
242
243!
244!-- Velocities at nzb+1 needed for volume flow control
245    IF ( conserve_volume_flow )  THEN
246       ALLOCATE( u_nzb_p1_for_vfc(nys:nyn), v_nzb_p1_for_vfc(nxl:nxr) )
247       u_nzb_p1_for_vfc = 0.0
248       v_nzb_p1_for_vfc = 0.0
249    ENDIF
250
251!
252!-- Arrays to store velocity data from t-dt and the phase speeds which
253!-- are needed for radiation boundary conditions
254    IF ( outflow_l )  THEN
255       ALLOCATE( u_m_l(nzb:nzt+1,nys-1:nyn+1,1:2), &
256                 v_m_l(nzb:nzt+1,nys-1:nyn+1,0:1), &
257                 w_m_l(nzb:nzt+1,nys-1:nyn+1,0:1) )
258    ENDIF
259    IF ( outflow_r )  THEN
260       ALLOCATE( u_m_r(nzb:nzt+1,nys-1:nyn+1,nx-1:nx), &
261                 v_m_r(nzb:nzt+1,nys-1:nyn+1,nx-1:nx), &
262                 w_m_r(nzb:nzt+1,nys-1:nyn+1,nx-1:nx) )
263    ENDIF
264    IF ( outflow_l  .OR.  outflow_r )  THEN
265       ALLOCATE( c_u(nzb:nzt+1,nys-1:nyn+1), c_v(nzb:nzt+1,nys-1:nyn+1), &
266                 c_w(nzb:nzt+1,nys-1:nyn+1) )
267    ENDIF
268    IF ( outflow_s )  THEN
269       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxl-1:nxr+1), &
270                 v_m_s(nzb:nzt+1,1:2,nxl-1:nxr+1), &
271                 w_m_s(nzb:nzt+1,0:1,nxl-1:nxr+1) )
272    ENDIF
273    IF ( outflow_n )  THEN
274       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxl-1:nxr+1), &
275                 v_m_n(nzb:nzt+1,ny-1:ny,nxl-1:nxr+1), &
276                 w_m_n(nzb:nzt+1,ny-1:ny,nxl-1:nxr+1) )
277    ENDIF
278    IF ( outflow_s  .OR.  outflow_n )  THEN
279       ALLOCATE( c_u(nzb:nzt+1,nxl-1:nxr+1), c_v(nzb:nzt+1,nxl-1:nxr+1), &
280                 c_w(nzb:nzt+1,nxl-1:nxr+1) )
281    ENDIF
282
283!
284!-- Initial assignment of the pointers
285    IF ( timestep_scheme(1:5) /= 'runge' )  THEN
286
287       rif_m   => rif_1;    rif   => rif_2
288       shf_m   => shf_1;    shf   => shf_2
289       tswst_m => tswst_1;  tswst => tswst_2
290       usws_m  => usws_1;   usws  => usws_2
291       uswst_m => uswst_1;  uswst => uswst_2
292       vsws_m  => vsws_1;   vsws  => vsws_2
293       vswst_m => vswst_1;  vswst => vswst_2
294       e_m  => e_1;   e  => e_2;   e_p  => e_3;   te_m  => e_3
295       kh_m => kh_1;  kh => kh_2
296       km_m => km_1;  km => km_2
297       pt_m => pt_1;  pt => pt_2;  pt_p => pt_3;  tpt_m => pt_3
298       u_m  => u_1;   u  => u_2;   u_p  => u_3;   tu_m  => u_3
299       v_m  => v_1;   v  => v_2;   v_p  => v_3;   tv_m  => v_3
300       w_m  => w_1;   w  => w_2;   w_p  => w_3;   tw_m  => w_3
301
302       IF ( humidity  .OR.  passive_scalar )  THEN
303          qsws_m  => qsws_1;   qsws  => qsws_2
304          qswst_m => qswst_1;  qswst => qswst_2
305          q_m    => q_1;     q    => q_2;     q_p => q_3;     tq_m => q_3
306          IF ( humidity )        vpt_m  => vpt_1;   vpt  => vpt_2
307          IF ( cloud_physics )   ql   => ql_1
308          IF ( cloud_droplets )  THEN
309             ql   => ql_1
310             ql_c => ql_2
311          ENDIF
312       ENDIF
313
314    ELSE
315
316       rif   => rif_1
317       shf   => shf_1
318       tswst => tswst_1
319       usws  => usws_1
320       uswst => uswst_1
321       vsws  => vsws_1
322       vswst => vswst_1
323       e     => e_1;   e_p  => e_2;   te_m  => e_3;   e_m  => e_3
324       kh    => kh_1
325       km    => km_1
326       pt    => pt_1;  pt_p => pt_2;  tpt_m => pt_3;  pt_m => pt_3
327       u     => u_1;   u_p  => u_2;   tu_m  => u_3;   u_m  => u_3
328       v     => v_1;   v_p  => v_2;   tv_m  => v_3;   v_m  => v_3
329       w     => w_1;   w_p  => w_2;   tw_m  => w_3;   w_m  => w_3
330
331       IF ( humidity  .OR.  passive_scalar )  THEN
332          qsws   => qsws_1
333          qswst  => qswst_1
334          q      => q_1;     q_p  => q_2;     tq_m  => q_3;    q_m => q_3
335          IF ( humidity )        vpt  => vpt_1
336          IF ( cloud_physics )   ql   => ql_1
337          IF ( cloud_droplets )  THEN
338             ql   => ql_1
339             ql_c => ql_2
340          ENDIF
341       ENDIF
342
343       IF ( ocean )  THEN
344          saswsb => saswsb_1
345          saswst => saswst_1
346          sa     => sa_1;    sa_p => sa_2;    tsa_m => sa_3
347       ENDIF
348
349    ENDIF
350
351!
352!-- Initialize model variables
353    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
354!
355!--    First model run of a possible job queue.
356!--    Initial profiles of the variables must be computes.
357       IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
358!
359!--       Use solutions of the 1D model as initial profiles,
360!--       start 1D model
361          CALL init_1d_model
362!
363!--       Transfer initial profiles to the arrays of the 3D model
364          DO  i = nxl-1, nxr+1
365             DO  j = nys-1, nyn+1
366                e(:,j,i)  = e1d
367                kh(:,j,i) = kh1d
368                km(:,j,i) = km1d
369                pt(:,j,i) = pt_init
370                u(:,j,i)  = u1d
371                v(:,j,i)  = v1d
372             ENDDO
373          ENDDO
374
375          IF ( humidity  .OR.  passive_scalar )  THEN
376             DO  i = nxl-1, nxr+1
377                DO  j = nys-1, nyn+1
378                   q(:,j,i) = q_init
379                ENDDO
380             ENDDO
381          ENDIF
382
383          IF ( .NOT. constant_diffusion )  THEN
384             DO  i = nxl-1, nxr+1
385                DO  j = nys-1, nyn+1
386                   e(:,j,i)  = e1d
387                ENDDO
388             ENDDO
389!
390!--          Store initial profiles for output purposes etc.
391             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
392
393             IF ( prandtl_layer )  THEN
394                rif  = rif1d(nzb+1)
395                ts   = 0.0  ! could actually be computed more accurately in the
396                            ! 1D model. Update when opportunity arises.
397                us   = us1d
398                usws = usws1d
399                vsws = vsws1d
400             ELSE
401                ts   = 0.0  ! must be set, because used in
402                rif  = 0.0  ! flowste
403                us   = 0.0
404                usws = 0.0
405                vsws = 0.0
406             ENDIF
407
408          ELSE
409             e    = 0.0  ! must be set, because used in
410             rif  = 0.0  ! flowste
411             ts   = 0.0
412             us   = 0.0
413             usws = 0.0
414             vsws = 0.0
415          ENDIF
416          uswst = top_momentumflux_u
417          vswst = top_momentumflux_v
418
419!
420!--       In every case qs = 0.0 (see also pt)
421!--       This could actually be computed more accurately in the 1D model.
422!--       Update when opportunity arises!
423          IF ( humidity  .OR.  passive_scalar )  qs = 0.0
424
425!
426!--       inside buildings set velocities back to zero
427          IF ( topography /= 'flat' )  THEN
428             DO  i = nxl-1, nxr+1
429                DO  j = nys-1, nyn+1
430                   u(nzb:nzb_u_inner(j,i),j,i) = 0.0
431                   v(nzb:nzb_v_inner(j,i),j,i) = 0.0
432                ENDDO
433             ENDDO
434             IF ( conserve_volume_flow )  THEN
435                IF ( nxr == nx )  THEN
436                   DO  j = nys, nyn
437                      DO  k = nzb + 1, nzb_u_inner(j,nx)
438                         u_nzb_p1_for_vfc(j) = u1d(k) * dzu(k)
439                      ENDDO
440                   ENDDO
441                ENDIF
442                IF ( nyn == ny )  THEN
443                   DO  i = nxl, nxr
444                      DO  k = nzb + 1, nzb_v_inner(ny,i)
445                         v_nzb_p1_for_vfc(i) = v1d(k) * dzu(k)
446                      ENDDO
447                   ENDDO
448                ENDIF
449             ENDIF
450!
451!--          WARNING: The extra boundary conditions set after running the
452!--          -------  1D model impose an error on the divergence one layer
453!--                   below the topography; need to correct later
454!--          ATTENTION: Provisional correction for Piacsek & Williams
455!--          ---------  advection scheme: keep u and v zero one layer below
456!--                     the topography.
457             IF ( ibc_uv_b == 0 )  THEN
458!
459!--             Satisfying the Dirichlet condition with an extra layer below
460!--             the surface where the u and v component change their sign.
461                DO  i = nxl-1, nxr+1
462                   DO  j = nys-1, nyn+1
463                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = -u(1,j,i)
464                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = -v(1,j,i)
465                   ENDDO
466                ENDDO
467
468             ELSE
469!
470!--             Neumann condition
471                DO  i = nxl-1, nxr+1
472                   DO  j = nys-1, nyn+1
473                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i)
474                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i)
475                   ENDDO
476                ENDDO
477
478             ENDIF
479
480          ENDIF
481
482       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 ) &
483       THEN
484!
485!--       Use constructed initial profiles (velocity constant with height,
486!--       temperature profile with constant gradient)
487          DO  i = nxl-1, nxr+1
488             DO  j = nys-1, nyn+1
489                pt(:,j,i) = pt_init
490                u(:,j,i)  = u_init
491                v(:,j,i)  = v_init
492             ENDDO
493          ENDDO
494
495!
496!--       Set initial horizontal velocities at the lowest computational grid levels
497!--       to zero in order to avoid too small time steps caused by the diffusion
498!--       limit in the initial phase of a run (at k=1, dz/2 occurs in the
499!--       limiting formula!). The original values are stored to be later used for
500!--       volume flow control.
501          DO  i = nxl-1, nxr+1
502             DO  j = nys-1, nyn+1
503                u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0
504                v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0
505             ENDDO
506          ENDDO
507          IF ( conserve_volume_flow )  THEN
508             IF ( nxr == nx )  THEN
509                DO  j = nys, nyn
510                   DO  k = nzb + 1, nzb_u_inner(j,nx) + 1
511                      u_nzb_p1_for_vfc(j) = u_init(k) * dzu(k)
512                   ENDDO
513                ENDDO
514             ENDIF
515             IF ( nyn == ny )  THEN
516                DO  i = nxl, nxr
517                   DO  k = nzb + 1, nzb_v_inner(ny,i) + 1
518                      v_nzb_p1_for_vfc(i) = v_init(k) * dzu(k)
519                   ENDDO
520                ENDDO
521             ENDIF
522          ENDIF
523
524          IF ( humidity  .OR.  passive_scalar )  THEN
525             DO  i = nxl-1, nxr+1
526                DO  j = nys-1, nyn+1
527                   q(:,j,i) = q_init
528                ENDDO
529             ENDDO
530          ENDIF
531
532          IF ( ocean )  THEN
533             DO  i = nxl-1, nxr+1
534                DO  j = nys-1, nyn+1
535                   sa(:,j,i) = sa_init
536                ENDDO
537             ENDDO
538          ENDIF
539         
540          IF ( constant_diffusion )  THEN
541             km   = km_constant
542             kh   = km / prandtl_number
543             e    = 0.0
544          ELSEIF ( e_init > 0.0 )  THEN
545             DO  k = nzb+1, nzt
546                km(k,:,:) = 0.1 * l_grid(k) * SQRT( e_init )
547             ENDDO
548             km(nzb,:,:)   = km(nzb+1,:,:)
549             km(nzt+1,:,:) = km(nzt,:,:)
550             kh   = km / prandtl_number
551             e    = e_init
552          ELSE
553             IF ( .NOT. ocean )  THEN
554                kh   = 0.01   ! there must exist an initial diffusion, because
555                km   = 0.01   ! otherwise no TKE would be produced by the
556                              ! production terms, as long as not yet
557                              ! e = (u*/cm)**2 at k=nzb+1
558             ELSE
559                kh   = 0.00001
560                km   = 0.00001
561             ENDIF
562             e    = 0.0
563          ENDIF
564          rif   = 0.0
565          ts    = 0.0
566          us    = 0.0
567          usws  = 0.0
568          uswst = top_momentumflux_u
569          vsws  = 0.0
570          vswst = top_momentumflux_v
571          IF ( humidity  .OR.  passive_scalar ) qs = 0.0
572
573!
574!--       Compute initial temperature field and other constants used in case
575!--       of a sloping surface
576          IF ( sloping_surface )  CALL init_slope
577
578       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 ) &
579       THEN
580!
581!--       Initialization will completely be done by the user
582          CALL user_init_3d_model
583
584       ENDIF
585
586!
587!--    apply channel flow boundary condition
588       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
589
590          u(nzt+1,:,:) = 0.0
591          v(nzt+1,:,:) = 0.0
592
593!--       for the Dirichlet condition to be correctly applied at the top, set
594!--       ug and vg to zero there
595          ug(nzt+1)    = 0.0
596          vg(nzt+1)    = 0.0
597
598       ENDIF
599
600!
601!--    Calculate virtual potential temperature
602       IF ( humidity ) vpt = pt * ( 1.0 + 0.61 * q )
603
604!
605!--    Store initial profiles for output purposes etc.
606       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
607       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
608       IF ( ibc_uv_b == 0 )  THEN
609          hom(nzb,1,5,:) = -hom(nzb+1,1,5,:)  ! due to satisfying the Dirichlet
610          hom(nzb,1,6,:) = -hom(nzb+1,1,6,:)  ! condition with an extra layer
611              ! below the surface where the u and v component change their sign
612       ENDIF
613       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
614       hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 )
615       hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 )
616
617       IF ( ocean )  THEN
618!
619!--       Store initial salinity profile
620          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
621       ENDIF
622
623       IF ( humidity )  THEN
624!
625!--       Store initial profile of total water content, virtual potential
626!--       temperature
627          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
628          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
629          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
630!
631!--          Store initial profile of specific humidity and potential
632!--          temperature
633             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
634             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
635          ENDIF
636       ENDIF
637
638       IF ( passive_scalar )  THEN
639!
640!--       Store initial scalar profile
641          hom(:,1,26,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
642       ENDIF
643
644!
645!--    Initialize fluxes at bottom surface
646       IF ( use_surface_fluxes )  THEN
647
648          IF ( constant_heatflux )  THEN
649!
650!--          Heat flux is prescribed
651             IF ( random_heatflux )  THEN
652                CALL disturb_heatflux
653             ELSE
654                shf = surface_heatflux
655!
656!--             Over topography surface_heatflux is replaced by wall_heatflux(0)
657                IF ( TRIM( topography ) /= 'flat' )  THEN
658                   DO  i = nxl-1, nxr+1
659                      DO  j = nys-1, nyn+1
660                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
661                            shf(j,i) = wall_heatflux(0)
662                         ENDIF
663                      ENDDO
664                   ENDDO
665                ENDIF
666             ENDIF
667             IF ( ASSOCIATED( shf_m ) )  shf_m = shf
668          ENDIF
669
670!
671!--       Determine the near-surface water flux
672          IF ( humidity  .OR.  passive_scalar )  THEN
673             IF ( constant_waterflux )  THEN
674                qsws   = surface_waterflux
675                IF ( ASSOCIATED( qsws_m ) )  qsws_m = qsws
676             ENDIF
677          ENDIF
678
679       ENDIF
680
681!
682!--    Initialize fluxes at top surface
683!--    Currently, only the heatflux and salinity flux can be prescribed.
684!--    The latent flux is zero in this case!
685       IF ( use_top_fluxes )  THEN
686
687          IF ( constant_top_heatflux )  THEN
688!
689!--          Heat flux is prescribed
690             tswst = top_heatflux
691             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
692
693             IF ( humidity  .OR.  passive_scalar )  THEN
694                qswst = 0.0
695                IF ( ASSOCIATED( qswst_m ) )  qswst_m = qswst
696             ENDIF
697
698             IF ( ocean )  THEN
699                saswsb = bottom_salinityflux
700                saswst = top_salinityflux
701             ENDIF
702          ENDIF
703
704!
705!--       Initialization in case of a coupled model run
706          IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
707             tswst = 0.0
708             IF ( ASSOCIATED( tswst_m ) )  tswst_m = tswst
709          ENDIF
710
711       ENDIF
712
713!
714!--    Initialize Prandtl layer quantities
715       IF ( prandtl_layer )  THEN
716
717          z0 = roughness_length
718
719          IF ( .NOT. constant_heatflux )  THEN 
720!
721!--          Surface temperature is prescribed. Here the heat flux cannot be
722!--          simply estimated, because therefore rif, u* and theta* would have
723!--          to be computed by iteration. This is why the heat flux is assumed
724!--          to be zero before the first time step. It approaches its correct
725!--          value in the course of the first few time steps.
726             shf   = 0.0
727             IF ( ASSOCIATED( shf_m ) )  shf_m = 0.0
728          ENDIF
729
730          IF ( humidity  .OR.  passive_scalar )  THEN
731             IF ( .NOT. constant_waterflux )  THEN
732                qsws   = 0.0
733                IF ( ASSOCIATED( qsws_m ) )   qsws_m = 0.0
734             ENDIF
735          ENDIF
736
737       ENDIF
738
739!
740!--    Calculate the initial volume flow at the right and north boundary
741       IF ( conserve_volume_flow )  THEN
742
743          volume_flow_initial_l = 0.0
744          volume_flow_area_l    = 0.0
745 
746          IF ( nxr == nx )  THEN
747             DO  j = nys, nyn
748                DO  k = nzb_2d(j,nx) + 1, nzt
749                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
750                                              u(k,j,nx) * dzu(k)
751                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzu(k)
752                ENDDO
753!
754!--             Correction if velocity at nzb+1 has been set zero further above
755                volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
756                                           u_nzb_p1_for_vfc(j)
757             ENDDO
758          ENDIF
759
760          IF ( nyn == ny )  THEN
761             DO  i = nxl, nxr
762                DO  k = nzb_2d(ny,i) + 1, nzt 
763                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
764                                              v(k,ny,i) * dzu(k)
765                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzu(k)
766                ENDDO
767!
768!--             Correction if velocity at nzb+1 has been set zero further above
769                volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
770                                           v_nzb_p1_for_vfc(i)
771             ENDDO
772          ENDIF
773
774#if defined( __parallel )
775          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
776                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
777          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
778                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
779#else
780          volume_flow_initial = volume_flow_initial_l
781          volume_flow_area    = volume_flow_area_l
782#endif 
783       ENDIF
784
785!
786!--    For the moment, perturbation pressure and vertical velocity are zero
787       p = 0.0; w = 0.0
788
789!
790!--    Initialize array sums (must be defined in first call of pres)
791       sums = 0.0
792
793!
794!--    Treating cloud physics, liquid water content and precipitation amount
795!--    are zero at beginning of the simulation
796       IF ( cloud_physics )  THEN
797          ql = 0.0
798          IF ( precipitation )  precipitation_amount = 0.0
799       ENDIF
800
801!
802!--    Initialize spectra
803       IF ( dt_dosp /= 9999999.9 )  THEN
804          spectrum_x = 0.0
805          spectrum_y = 0.0
806       ENDIF
807
808!
809!--    Impose vortex with vertical axis on the initial velocity profile
810       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
811          CALL init_rankine
812       ENDIF
813
814!
815!--    Impose temperature anomaly (advection test only)
816       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
817          CALL init_pt_anomaly
818       ENDIF
819
820!
821!--    If required, change the surface temperature at the start of the 3D run
822       IF ( pt_surface_initial_change /= 0.0 )  THEN
823          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
824       ENDIF
825
826!
827!--    If required, change the surface humidity/scalar at the start of the 3D
828!--    run
829       IF ( ( humidity .OR. passive_scalar ) .AND. &
830            q_surface_initial_change /= 0.0 )  THEN
831          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
832       ENDIF
833
834!
835!--    Initialize the random number generator (from numerical recipes)
836       CALL random_function_ini
837
838!
839!--    Impose random perturbation on the horizontal velocity field and then
840!--    remove the divergences from the velocity field
841       IF ( create_disturbances )  THEN
842          CALL disturb_field( nzb_u_inner, tend, u )
843          CALL disturb_field( nzb_v_inner, tend, v )
844          n_sor = nsor_ini
845          CALL pres
846          n_sor = nsor
847       ENDIF
848
849!
850!--    Once again set the perturbation pressure explicitly to zero in order to
851!--    assure that it does not generate any divergences in the first time step.
852!--    At t=0 the velocity field is free of divergence (as constructed above).
853!--    Divergences being created during a time step are not yet known and thus
854!--    cannot be corrected during the time step yet.
855       p = 0.0
856
857!
858!--    Initialize old and new time levels.
859       IF ( timestep_scheme(1:5) /= 'runge' )  THEN
860          e_m = e; pt_m = pt; u_m = u; v_m = v; w_m = w; kh_m = kh; km_m = km
861       ELSE
862          te_m = 0.0; tpt_m = 0.0; tu_m = 0.0; tv_m = 0.0; tw_m = 0.0
863       ENDIF
864       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
865
866       IF ( humidity  .OR.  passive_scalar )  THEN
867          IF ( ASSOCIATED( q_m ) )  q_m = q
868          IF ( timestep_scheme(1:5) == 'runge' )  tq_m = 0.0
869          q_p = q
870          IF ( humidity  .AND.  ASSOCIATED( vpt_m ) )  vpt_m = vpt
871       ENDIF
872
873       IF ( ocean )  THEN
874          tsa_m = 0.0
875          sa_p  = sa
876       ENDIF
877
878!
879!--    Initialize old timelevels needed for radiation boundary conditions
880       IF ( outflow_l )  THEN
881          u_m_l(:,:,:) = u(:,:,1:2)
882          v_m_l(:,:,:) = v(:,:,0:1)
883          w_m_l(:,:,:) = w(:,:,0:1)
884       ENDIF
885       IF ( outflow_r )  THEN
886          u_m_r(:,:,:) = u(:,:,nx-1:nx)
887          v_m_r(:,:,:) = v(:,:,nx-1:nx)
888          w_m_r(:,:,:) = w(:,:,nx-1:nx)
889       ENDIF
890       IF ( outflow_s )  THEN
891          u_m_s(:,:,:) = u(:,0:1,:)
892          v_m_s(:,:,:) = v(:,1:2,:)
893          w_m_s(:,:,:) = w(:,0:1,:)
894       ENDIF
895       IF ( outflow_n )  THEN
896          u_m_n(:,:,:) = u(:,ny-1:ny,:)
897          v_m_n(:,:,:) = v(:,ny-1:ny,:)
898          w_m_n(:,:,:) = w(:,ny-1:ny,:)
899       ENDIF
900
901    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data' ) &
902    THEN
903!
904!--    Read binary data from restart file
905       CALL read_3d_binary
906
907!
908!--    Calculate initial temperature field and other constants used in case
909!--    of a sloping surface
910       IF ( sloping_surface )  CALL init_slope
911
912!
913!--    Initialize new time levels (only done in order to set boundary values
914!--    including ghost points)
915       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
916       IF ( humidity  .OR.  passive_scalar )  q_p = q
917       IF ( ocean )  sa_p = sa
918
919    ELSE
920!
921!--    Actually this part of the programm should not be reached
922       IF ( myid == 0 )  PRINT*,'+++ init_3d_model: unknown initializing ', &
923                                                    'problem'
924       CALL local_stop
925    ENDIF
926
927!
928!-- If required, initialize dvrp-software
929    IF ( dt_dvrp /= 9999999.9 )  CALL init_dvrp
930
931    IF ( ocean )  THEN
932!
933!--    Initialize quantities needed for the ocean model
934       CALL init_ocean
935    ELSE
936!
937!--    Initialize quantities for handling cloud physics
938!--    This routine must be called before init_particles, because
939!--    otherwise, array pt_d_t, needed in data_output_dvrp (called by
940!--    init_particles) is not defined.
941       CALL init_cloud_physics
942    ENDIF
943
944!
945!-- If required, initialize particles
946    IF ( particle_advection )  CALL init_particles
947
948!
949!-- Initialize quantities for special advections schemes
950    CALL init_advec
951
952!
953!-- Initialize Rayleigh damping factors
954    rdf = 0.0
955    IF ( rayleigh_damping_factor /= 0.0 )  THEN
956       IF ( .NOT. ocean )  THEN
957          DO  k = nzb+1, nzt
958             IF ( zu(k) >= rayleigh_damping_height )  THEN
959                rdf(k) = rayleigh_damping_factor * &
960                      ( SIN( pi * 0.5 * ( zu(k) - rayleigh_damping_height )    &
961                                      / ( zu(nzt) - rayleigh_damping_height ) )&
962                      )**2
963             ENDIF
964          ENDDO
965       ELSE
966          DO  k = nzt, nzb+1, -1
967             IF ( zu(k) <= rayleigh_damping_height )  THEN
968                rdf(k) = rayleigh_damping_factor * &
969                      ( SIN( pi * 0.5 * ( rayleigh_damping_height - zu(k) )    &
970                                      / ( rayleigh_damping_height - zu(nzb+1)))&
971                      )**2
972             ENDIF
973          ENDDO
974       ENDIF
975    ENDIF
976
977!
978!-- Initialize diffusivities used within the outflow damping layer in case of
979!-- non-cyclic lateral boundaries. A linear increase is assumed over the first
980!-- half of the width of the damping layer
981    IF ( bc_lr == 'dirichlet/radiation' )  THEN
982
983       DO  i = nxl-1, nxr+1
984          IF ( i >= nx - outflow_damping_width )  THEN
985             km_damp_x(i) = km_damp_max * MIN( 1.0,                    &
986                            ( i - ( nx - outflow_damping_width ) ) /   &
987                            REAL( outflow_damping_width/2 )            &
988                                             )
989          ELSE
990             km_damp_x(i) = 0.0
991          ENDIF
992       ENDDO
993
994    ELSEIF ( bc_lr == 'radiation/dirichlet' )  THEN
995
996       DO  i = nxl-1, nxr+1
997          IF ( i <= outflow_damping_width )  THEN
998             km_damp_x(i) = km_damp_max * MIN( 1.0,                    &
999                            ( outflow_damping_width - i ) /            &
1000                            REAL( outflow_damping_width/2 )            &
1001                                             )
1002          ELSE
1003             km_damp_x(i) = 0.0
1004          ENDIF
1005       ENDDO
1006
1007    ENDIF
1008
1009    IF ( bc_ns == 'dirichlet/radiation' )  THEN
1010
1011       DO  j = nys-1, nyn+1
1012          IF ( j >= ny - outflow_damping_width )  THEN
1013             km_damp_y(j) = km_damp_max * MIN( 1.0,                    &
1014                            ( j - ( ny - outflow_damping_width ) ) /   &
1015                            REAL( outflow_damping_width/2 )            &
1016                                             )
1017          ELSE
1018             km_damp_y(j) = 0.0
1019          ENDIF
1020       ENDDO
1021
1022    ELSEIF ( bc_ns == 'radiation/dirichlet' )  THEN
1023
1024       DO  j = nys-1, nyn+1
1025          IF ( j <= outflow_damping_width )  THEN
1026             km_damp_y(j) = km_damp_max * MIN( 1.0,                    &
1027                            ( outflow_damping_width - j ) /            &
1028                            REAL( outflow_damping_width/2 )            &
1029                                             )
1030          ELSE
1031             km_damp_y(j) = 0.0
1032          ENDIF
1033       ENDDO
1034
1035    ENDIF
1036
1037!
1038!-- Initialize local summation arrays for UP flow_statistics. This is necessary
1039!-- because they may not yet have been initialized when they are called from
1040!-- flow_statistics (or - depending on the chosen model run - are never
1041!-- initialized)
1042    sums_divnew_l      = 0.0
1043    sums_divold_l      = 0.0
1044    sums_l_l           = 0.0
1045    sums_up_fraction_l = 0.0
1046    sums_wsts_bc_l     = 0.0
1047
1048!
1049!-- Pre-set masks for regional statistics. Default is the total model domain.
1050    rmask = 1.0
1051
1052!
1053!-- User-defined initializing actions. Check afterwards, if maximum number
1054!-- of allowed timeseries is not exceeded
1055    CALL user_init
1056
1057    IF ( dots_num > dots_max )  THEN
1058       IF ( myid == 0 )  THEN
1059          PRINT*, '+++ user_init: number of time series quantities exceeds', &
1060                  ' its maximum of dots_max = ', dots_max
1061          PRINT*, '    Please increase dots_max in modules.f90.'
1062       ENDIF
1063       CALL local_stop
1064    ENDIF
1065
1066!
1067!-- Input binary data file is not needed anymore. This line must be placed
1068!-- after call of user_init!
1069    CALL close_file( 13 )
1070
1071!
1072!-- Compute total sum of active mask grid points
1073!-- ngp_2dh: number of grid points of a horizontal cross section through the
1074!--          total domain
1075!-- ngp_3d:  number of grid points of the total domain
1076    ngp_2dh_outer_l   = 0
1077    ngp_2dh_outer     = 0
1078    ngp_2dh_s_inner_l = 0
1079    ngp_2dh_s_inner   = 0
1080    ngp_2dh_l         = 0
1081    ngp_2dh           = 0
1082    ngp_3d_inner_l    = 0
1083    ngp_3d_inner      = 0
1084    ngp_3d            = 0
1085    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
1086
1087    DO  sr = 0, statistic_regions
1088       DO  i = nxl, nxr
1089          DO  j = nys, nyn
1090             IF ( rmask(j,i,sr) == 1.0 )  THEN
1091!
1092!--             All xy-grid points
1093                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
1094!
1095!--             xy-grid points above topography
1096                DO  k = nzb_s_outer(j,i), nz + 1
1097                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1
1098                ENDDO
1099                DO  k = nzb_s_inner(j,i), nz + 1
1100                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) + 1
1101                ENDDO
1102!
1103!--             All grid points of the total domain above topography
1104                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr) + &
1105                                     ( nz - nzb_s_inner(j,i) + 2 )
1106             ENDIF
1107          ENDDO
1108       ENDDO
1109    ENDDO
1110
1111    sr = statistic_regions + 1
1112#if defined( __parallel )
1113    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,  &
1114                        comm2d, ierr )
1115    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr, &
1116                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
1117    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),        &
1118                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
1119    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner(0), sr, MPI_INTEGER, &
1120                        MPI_SUM, comm2d, ierr )
1121#else
1122    ngp_2dh         = ngp_2dh_l
1123    ngp_2dh_outer   = ngp_2dh_outer_l
1124    ngp_2dh_s_inner = ngp_2dh_s_inner_l
1125    ngp_3d_inner    = ngp_3d_inner_l
1126#endif
1127
1128    ngp_3d = ngp_2dh * ( nz + 2 )
1129
1130!
1131!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
1132!-- buoyancy, etc. A zero value will occur for cases where all grid points of
1133!-- the respective subdomain lie below the surface topography
1134    ngp_2dh_outer = MAX( 1, ngp_2dh_outer(:,:) ) 
1135    ngp_3d_inner  = MAX( 1, ngp_3d_inner(:)    )
1136
1137    DEALLOCATE( ngp_2dh_l, ngp_2dh_outer_l, ngp_3d_inner_l )
1138
1139
1140 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.