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

Last change on this file since 680 was 680, checked in by gryschka, 13 years ago

message string

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