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

Last change on this file since 1316 was 1316, checked in by heinze, 10 years ago

bugfix: allocation of w_subs

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