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

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