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

Last change on this file since 1683 was 1683, checked in by knoop, 9 years ago

last commit documented

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