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

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