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

Last change on this file since 1723 was 1708, checked in by maronga, 8 years ago

last commit documented

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