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

Last change on this file since 708 was 708, checked in by raasch, 13 years ago

last commit documented

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