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

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

last commit documented

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