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

Last change on this file since 1266 was 1242, checked in by heinze, 11 years ago

Last commmit documented

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