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

Last change on this file since 147 was 147, checked in by raasch, 16 years ago

further updates for turbulent inflow: reading input data of a precursor run using a smaller total domain is working

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