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

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

last commit documented

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