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

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

enable usage of large_scale subsidence in combination with large_scale_forcing

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