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

Last change on this file since 1318 was 1318, checked in by raasch, 10 years ago

former files/routines cpu_log and cpu_statistics combined to one module,
which also includes the former data module cpulog from the modules-file,
module interfaces removed

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