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

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

cpp-switches removed + cpp-bugfixes + zero-settings for velocities inside topography re-activated

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