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

Last change on this file since 1734 was 1734, checked in by raasch, 9 years ago

no initial disturbances in case that the disturbance energy limit has been set zero

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