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

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

pointer free version can be generated with cpp switch nopointer

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