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

Last change on this file since 1041 was 1037, checked in by raasch, 12 years ago

last commit documented

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