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

Last change on this file since 1517 was 1508, checked in by suehring, 10 years ago

last commit documented

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