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

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

routines for nudging and large scale forcing from external file added

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