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

Last change on this file since 1736 was 1735, checked in by raasch, 8 years ago

last commit documented

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