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

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

land surface model released

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