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

Last change on this file since 1552 was 1552, checked in by maronga, 9 years ago

last commit documented

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