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

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

last commit documented

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