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

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

last commit documented

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