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

Last change on this file since 689 was 681, checked in by gryschka, 14 years ago

last commit documented

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