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

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

New:
---

Flow field initialization with given (e.g. measured) profiles. Profile data
for u-,v-velocity components + respective heights are given with new
inipar-parameters u_profile, v_profile, and uv_heights. Final profiles are
calculated from these given profiles by linear interpolation.
(check_parameters, header, init_3d_model, modules, parin)

Changed:


ug,vg replaced by u_init,v_init as the Dirichlet top boundary condition
(boundary_conds)

dirichlet_0 conditions moved from init_3d_model to
check_parameters (check_parameters, init_3d_model)

Errors:


bugfix: dirichlet_0 conditions moved from init_3d_model to
check_parameters (check_parameters, init_3d_model)

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