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

Last change on this file since 1576 was 1576, checked in by raasch, 9 years ago

last commit documented

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