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

Last change on this file since 1369 was 1362, checked in by hoffmann, 11 years ago

last commit documented

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