source: palm/tags/release-3.8/SOURCE/init_3d_model.f90 @ 2704

Last change on this file since 2704 was 714, checked in by suehring, 13 years ago

last commit documented

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