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

Last change on this file since 766 was 760, checked in by raasch, 13 years ago

last commit documented

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