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

Last change on this file since 1317 was 1317, checked in by heinze, 10 years ago

last commit documented

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