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

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

Parallel random number generator added (preliminary version).

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