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

Last change on this file since 1001 was 1001, checked in by raasch, 12 years ago

leapfrog timestep scheme and upstream-spline advection scheme completely removed from the code,
reading of dt_fixed from restart file removed

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