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

Last change on this file since 1771 was 1765, checked in by raasch, 8 years ago

last commit documented

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