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

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

last commit documented

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