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

Last change on this file since 1493 was 1485, checked in by kanani, 10 years ago

last commit documented

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