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

Last change on this file since 1188 was 1182, checked in by raasch, 11 years ago

last commit documented, rc-file for example run updated

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