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

Last change on this file since 1832 was 1832, checked in by hoffmann, 8 years ago

last commit documented

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