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

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

last commit documented

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