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

Last change on this file since 1658 was 1616, checked in by suehring, 9 years ago

last commit documented

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