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

Last change on this file since 1849 was 1849, checked in by hoffmann, 9 years ago

lpm_droplet_condensation improved, microphysics partially modularized

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