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

Last change on this file since 1403 was 1403, checked in by raasch, 10 years ago

last commit documented

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