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

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

unused variables remove from several routines

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