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

Last change on this file since 1817 was 1817, checked in by maronga, 9 years ago

further modularization of land surface model

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