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

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

code has been put under the GNU General Public License (v3)

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