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

Last change on this file since 1609 was 1586, checked in by maronga, 10 years ago

last commit documented

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