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

Last change on this file since 1079 was 1073, checked in by raasch, 12 years ago

FORTRAN bugfix for r1072

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