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

Last change on this file since 673 was 673, checked in by suehring, 13 years ago

Right computation of the pressure using Runge-Kutta weighting coefficients. Consideration of the pressure gradient during the time integration removed. Removed bugfix concerning velocity variances.

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