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

Last change on this file since 1311 was 1310, checked in by raasch, 11 years ago

update of GPL copyright

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