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

Last change on this file since 1319 was 1319, checked in by raasch, 10 years ago

last commit documented

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