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

Last change on this file since 1026 was 1026, checked in by letzel, 12 years ago

last commit documented

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