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

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

New:


openACC porting of reduction operations
additional 3D-flag arrays for replacing the 2D-index arrays nzb_s_inner and nzb_diff_s_inner
(flow_statistics, init_grid, init_3d_model, modules, palm, pres, time_integration)

Changed:


for PGI/openACC performance reasons set default compile options for openACC to "-ta=nocache",
and set environment variable PGI_ACC_SYNCHRONOUS=1
(MAKE.inc.pgi.openacc, palm_simple_run)

wall_flags_0 changed to 32bit INTEGER, additional array wall_flags_00 introduced to hold
bits 32-63
(advec_ws, init_grid, modules, palm)

Errors:


dummy argument tri in 1d-routines replaced by tri_for_1d because of name
conflict with arry tri in module arrays_3d
(tridia_solver)

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