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

Last change on this file since 1004 was 1004, checked in by raasch, 12 years ago

last commit documented

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