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

Last change on this file since 1507 was 1507, checked in by suehring, 9 years ago

Bugfix: set horizontal velocity components to zero inside topography

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