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

Last change on this file since 1788 was 1788, checked in by maronga, 9 years ago

added support for water and paved surfaced in land surface model / minor changes

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