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

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

last commit documented

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