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

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

last commit documented

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