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

Last change on this file since 1428 was 1412, checked in by suehring, 10 years ago

last commit documented

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