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

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

various bugfixes and modifications of the atmosphere-land-surface-radiation interaction. Completely re-written routine to calculate surface fluxes (surface_layer_fluxes.f90) that replaces prandtl_fluxes. Minor formatting corrections and renamings

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