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

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