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

Last change on this file since 1407 was 1407, checked in by raasch, 11 years ago

last commit documented

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