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

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

output of location messages to stdout

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