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

Last change on this file since 1827 was 1827, checked in by maronga, 8 years ago

last commit documented

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