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

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

preliminary bugfix for volume flow control

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