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

Last change on this file since 1430 was 1430, checked in by knoop, 10 years ago

last commit documented

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