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

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

last commit documented

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