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

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

output of location messages complemented, output of location bar added
(Makefile, check_parameters, cpulog, init_pegrid, init_3d_model, message, palm, parin, time_integration, new: progress_bar)
preprocessor switch intel_compiler added, -r8 compiler options removed
(.mrun.config.default, .mrun.config.imuk, .mrun.config.kiaps)

batch_job added to envpar-NAMELIST
(mrun, parin)

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