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

Last change on this file since 1835 was 1834, checked in by raasch, 9 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 70.0 KB
Line 
1!> @file init_3d_model.f90
2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2016 Leibniz Universitaet Hannover
17!--------------------------------------------------------------------------------!
18!
19! Current revisions:
20! ------------------
21!
22!
23! Former revisions:
24! -----------------
25! $Id: init_3d_model.f90 1834 2016-04-07 14:34:20Z raasch $
26!
27! 1833 2016-04-07 14:23:03Z raasch
28! initialization of spectra quantities moved to spectra_mod
29!
30! 1831 2016-04-07 13:15:51Z hoffmann
31! turbulence renamed collision_turbulence
32!
33! 1826 2016-04-07 12:01:39Z maronga
34! Renamed radiation calls.
35! Renamed canopy model calls.
36!
37! 1822 2016-04-07 07:49:42Z hoffmann
38! icloud_scheme replaced by microphysics_*
39!
40! 1817 2016-04-06 15:44:20Z maronga
41! Renamed lsm calls.
42!
43! 1815 2016-04-06 13:49:59Z raasch
44! zero-settings for velocities inside topography re-activated (was deactivated
45! in r1762)
46!
47! 1788 2016-03-10 11:01:04Z maronga
48! Added z0q.
49! Syntax layout improved.
50!
51! 1783 2016-03-06 18:36:17Z raasch
52! netcdf module name changed + related changes
53!
54! 1764 2016-02-28 12:45:19Z raasch
55! bugfix: increase size of volume_flow_area_l and volume_flow_initial_l by 1
56!
57! 1762 2016-02-25 12:31:13Z hellstea
58! Introduction of nested domain feature
59!
60! 1738 2015-12-18 13:56:05Z raasch
61! calculate mean surface level height for each statistic region
62!
63! 1734 2015-12-02 12:17:12Z raasch
64! no initial disturbances in case that the disturbance energy limit has been
65! set zero
66!
67! 1707 2015-11-02 15:24:52Z maronga
68! Bugfix: transfer of Richardson number from 1D model to Obukhov length caused
69! devision by zero in neutral stratification
70!
71! 1691 2015-10-26 16:17:44Z maronga
72! Call to init_surface_layer added. rif is replaced by ol and zeta.
73!
74! 1682 2015-10-07 23:56:08Z knoop
75! Code annotations made doxygen readable
76!
77! 1615 2015-07-08 18:49:19Z suehring
78! Enable turbulent inflow for passive_scalar and humidity
79!
80! 1585 2015-04-30 07:05:52Z maronga
81! Initialization of radiation code is now done after LSM initializtion
82!
83! 1575 2015-03-27 09:56:27Z raasch
84! adjustments for psolver-queries
85!
86! 1551 2015-03-03 14:18:16Z maronga
87! Allocation of land surface arrays is now done in the subroutine lsm_init_arrays,
88! which is part of land_surface_model.
89!
90! 1507 2014-12-10 12:14:18Z suehring
91! Bugfix: set horizontal velocity components to zero inside topography
92!
93! 1496 2014-12-02 17:25:50Z maronga
94! Added initialization of the land surface and radiation schemes
95!
96! 1484 2014-10-21 10:53:05Z kanani
97! Changes due to new module structure of the plant canopy model:
98! canopy-related initialization (e.g. lad and canopy_heat_flux) moved to new
99! subroutine init_plant_canopy within the module plant_canopy_model_mod,
100! call of subroutine init_plant_canopy added.
101!
102! 1431 2014-07-15 14:47:17Z suehring
103! var_d added, in order to normalize spectra.
104!
105! 1429 2014-07-15 12:53:45Z knoop
106! Ensemble run capability added to parallel random number generator
107!
108! 1411 2014-05-16 18:01:51Z suehring
109! Initial horizontal velocity profiles were not set to zero at the first vertical
110! grid level in case of non-cyclic lateral boundary conditions.
111!
112! 1406 2014-05-16 13:47:01Z raasch
113! bugfix: setting of initial velocities at k=1 to zero not in case of a
114! no-slip boundary condition for uv
115!
116! 1402 2014-05-09 14:25:13Z raasch
117! location messages modified
118!
119! 1400 2014-05-09 14:03:54Z knoop
120! Parallel random number generator added
121!
122! 1384 2014-05-02 14:31:06Z raasch
123! location messages added
124!
125! 1361 2014-04-16 15:17:48Z hoffmann
126! tend_* removed
127! Bugfix: w_subs is not allocated anymore if it is already allocated
128!
129! 1359 2014-04-11 17:15:14Z hoffmann
130! module lpm_init_mod added to use statements, because lpm_init has become a
131! module
132!
133! 1353 2014-04-08 15:21:23Z heinze
134! REAL constants provided with KIND-attribute
135!
136! 1340 2014-03-25 19:45:13Z kanani
137! REAL constants defined as wp-kind
138!
139! 1322 2014-03-20 16:38:49Z raasch
140! REAL constants defined as wp-kind
141! module interfaces removed
142!
143! 1320 2014-03-20 08:40:49Z raasch
144! ONLY-attribute added to USE-statements,
145! kind-parameters added to all INTEGER and REAL declaration statements,
146! kinds are defined in new module kinds,
147! revision history before 2012 removed,
148! comment fields (!:) to be used for variable explanations added to
149! all variable declaration statements
150!
151! 1316 2014-03-17 07:44:59Z heinze
152! Bugfix: allocation of w_subs
153!
154! 1299 2014-03-06 13:15:21Z heinze
155! Allocate w_subs due to extension of large scale subsidence in combination
156! with large scale forcing data (LSF_DATA)
157!
158! 1241 2013-10-30 11:36:58Z heinze
159! Overwrite initial profiles in case of nudging
160! Inititialize shf and qsws in case of large_scale_forcing
161!
162! 1221 2013-09-10 08:59:13Z raasch
163! +rflags_s_inner in copyin statement, use copyin for most arrays instead of
164! copy
165!
166! 1212 2013-08-15 08:46:27Z raasch
167! array tri is allocated and included in data copy statement
168!
169! 1195 2013-07-01 12:27:57Z heinze
170! Bugfix: move allocation of ref_state to parin.f90 and read_var_list.f90
171!
172! 1179 2013-06-14 05:57:58Z raasch
173! allocate and set ref_state to be used in buoyancy terms
174!
175! 1171 2013-05-30 11:27:45Z raasch
176! diss array is allocated with full size if accelerator boards are used
177!
178! 1159 2013-05-21 11:58:22Z fricke
179! -bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir
180!
181! 1153 2013-05-10 14:33:08Z raasch
182! diss array is allocated with dummy elements even if it is not needed
183! (required by PGI 13.4 / CUDA 5.0)
184!
185! 1115 2013-03-26 18:16:16Z hoffmann
186! unused variables removed
187!
188! 1113 2013-03-10 02:48:14Z raasch
189! openACC directive modified
190!
191! 1111 2013-03-08 23:54:10Z raasch
192! openACC directives added for pres
193! array diss allocated only if required
194!
195! 1092 2013-02-02 11:24:22Z raasch
196! unused variables removed
197!
198! 1065 2012-11-22 17:42:36Z hoffmann
199! allocation of diss (dissipation rate) in case of turbulence = .TRUE. added
200!
201! 1053 2012-11-13 17:11:03Z hoffmann
202! allocation and initialisation of necessary data arrays for the two-moment
203! cloud physics scheme the two new prognostic equations (nr, qr):
204! +dr, lambda_r, mu_r, sed_*, xr, *s, *sws, *swst, *, *_p, t*_m, *_1, *_2, *_3,
205! +tend_*, prr
206!
207! 1036 2012-10-22 13:43:42Z raasch
208! code put under GPL (PALM 3.9)
209!
210! 1032 2012-10-21 13:03:21Z letzel
211! save memory by not allocating pt_2 in case of neutral = .T.
212!
213! 1025 2012-10-07 16:04:41Z letzel
214! bugfix: swap indices of mask for ghost boundaries
215!
216! 1015 2012-09-27 09:23:24Z raasch
217! mask is set to zero for ghost boundaries
218!
219! 1010 2012-09-20 07:59:54Z raasch
220! cpp switch __nopointer added for pointer free version
221!
222! 1003 2012-09-14 14:35:53Z raasch
223! nxra,nyna, nzta replaced ny nxr, nyn, nzt
224!
225! 1001 2012-09-13 14:08:46Z raasch
226! all actions concerning leapfrog scheme removed
227!
228! 996 2012-09-07 10:41:47Z raasch
229! little reformatting
230!
231! 978 2012-08-09 08:28:32Z fricke
232! outflow damping layer removed
233! roughness length for scalar quantites z0h added
234! damping zone for the potential temperatur in case of non-cyclic lateral
235! boundaries added
236! initialization of ptdf_x, ptdf_y
237! initialization of c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l
238!
239! 849 2012-03-15 10:35:09Z raasch
240! init_particles renamed lpm_init
241!
242! 825 2012-02-19 03:03:44Z raasch
243! wang_collision_kernel renamed wang_kernel
244!
245! Revision 1.1  1998/03/09 16:22:22  raasch
246! Initial revision
247!
248!
249! Description:
250! ------------
251!> Allocation of arrays and initialization of the 3D model via
252!> a) pre-run the 1D model
253!> or
254!> b) pre-set constant linear profiles
255!> or
256!> c) read values of a previous run
257!------------------------------------------------------------------------------!
258 SUBROUTINE init_3d_model
259 
260
261    USE advec_ws
262
263    USE arrays_3d
264   
265    USE cloud_parameters,                                                      &
266        ONLY:  collision_turbulence, nc_const, precipitation_amount,           &
267               precipitation_rate, prr
268   
269    USE constants,                                                             &
270        ONLY:  pi
271   
272    USE control_parameters
273   
274    USE grid_variables,                                                        &
275        ONLY:  dx, dy
276   
277    USE indices
278
279    USE lpm_init_mod,                                                          &
280        ONLY:  lpm_init
281   
282    USE kinds
283
284    USE land_surface_model_mod,                                                &
285        ONLY:  lsm_init, lsm_init_arrays, land_surface
286 
287    USE ls_forcing_mod
288   
289    USE model_1d,                                                              &
290        ONLY:  e1d, kh1d, km1d, l1d, rif1d, u1d, us1d, usws1d, v1d, vsws1d 
291   
292    USE netcdf_interface,                                                      &
293        ONLY:  dots_max, dots_num
294   
295    USE particle_attributes,                                                   &
296        ONLY:  particle_advection, use_sgs_for_particles, wang_kernel
297   
298    USE pegrid
299   
300    USE plant_canopy_model_mod,                                                &
301        ONLY:  pcm_init, plant_canopy
302
303    USE radiation_model_mod,                                                   &
304        ONLY:  radiation_init, radiation
305   
306    USE random_function_mod 
307   
308    USE random_generator_parallel,                                             &
309        ONLY:  random_number_parallel, random_seed_parallel, random_dummy,     &
310               id_random_array, seq_random_array
311   
312    USE statistics,                                                            &
313        ONLY:  hom, hom_sum, mean_surface_level_height, pr_palm, rmask,        &
314               statistic_regions, sums, sums_divnew_l, sums_divold_l, sums_l,  &
315               sums_l_l, sums_up_fraction_l, sums_wsts_bc_l, ts_value,         &
316               weight_pres, weight_substep
317 
318    USE surface_layer_fluxes_mod,                                              &
319        ONLY:  init_surface_layer_fluxes
320   
321    USE transpose_indices 
322
323    IMPLICIT NONE
324
325    INTEGER(iwp) ::  i             !<
326    INTEGER(iwp) ::  ind_array(1)  !<
327    INTEGER(iwp) ::  j             !<
328    INTEGER(iwp) ::  k             !<
329    INTEGER(iwp) ::  sr            !<
330
331    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  ngp_2dh_l  !<
332
333    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l    !<
334    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner_l  !<
335
336    REAL(wp), DIMENSION(1:3) ::  volume_flow_area_l     !<
337    REAL(wp), DIMENSION(1:3) ::  volume_flow_initial_l  !<
338
339    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height_l    !<
340    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l    !<
341    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_tmp  !<
342
343
344    CALL location_message( 'allocating arrays', .FALSE. )
345!
346!-- Allocate arrays
347    ALLOCATE( mean_surface_level_height(0:statistic_regions),                  &
348              mean_surface_level_height_l(0:statistic_regions),                &
349              ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions),    &
350              ngp_3d(0:statistic_regions),                                     &
351              ngp_3d_inner(0:statistic_regions),                               &
352              ngp_3d_inner_l(0:statistic_regions),                             &
353              ngp_3d_inner_tmp(0:statistic_regions),                           &
354              sums_divnew_l(0:statistic_regions),                              &
355              sums_divold_l(0:statistic_regions) )
356    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
357    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                    &
358              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),                  &
359              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),                  &
360              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),                &
361              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),                  &
362              sums(nzb:nzt+1,pr_palm+max_pr_user),                             &
363              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),      &
364              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1),    &
365              sums_up_fraction_l(10,3,0:statistic_regions),                    &
366              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                   &
367              ts_value(dots_max,0:statistic_regions) )
368    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
369
370    ALLOCATE( ol(nysg:nyng,nxlg:nxrg), shf(nysg:nyng,nxlg:nxrg),               &
371              ts(nysg:nyng,nxlg:nxrg), tswst(nysg:nyng,nxlg:nxrg),             &
372              us(nysg:nyng,nxlg:nxrg), usws(nysg:nyng,nxlg:nxrg),              &
373              uswst(nysg:nyng,nxlg:nxrg), vsws(nysg:nyng,nxlg:nxrg),           &
374              vswst(nysg:nyng,nxlg:nxrg), z0(nysg:nyng,nxlg:nxrg),             &
375              z0h(nysg:nyng,nxlg:nxrg), z0q(nysg:nyng,nxlg:nxrg) )
376
377    ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr),                                    &
378              kh(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
379              km(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
380              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
381              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
382
383#if defined( __nopointer )
384    ALLOCATE( e(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
385              e_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
386              pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
387              pt_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
388              u(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
389              u_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
390              v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
391              v_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
392              w(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
393              w_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
394              te_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
395              tpt_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
396              tu_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
397              tv_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
398              tw_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
399#else
400    ALLOCATE( e_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
401              e_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
402              e_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
403              pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
404              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
405              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
406              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
407              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
408              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
409              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
410              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
411              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
412              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
413              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
414    IF (  .NOT.  neutral )  THEN
415       ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
416    ENDIF
417#endif
418
419!
420!-- Following array is required for perturbation pressure within the iterative
421!-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds
422!-- the weighted average of the substeps and cannot be used in the Poisson
423!-- solver.
424    IF ( psolver == 'sor' )  THEN
425       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
426    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
427!
428!--    For performance reasons, multigrid is using one ghost layer only
429       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
430    ENDIF
431
432!
433!-- Array for storing constant coeffficients of the tridiagonal solver
434    IF ( psolver == 'poisfft' )  THEN
435       ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) )
436       ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) )
437    ENDIF
438
439    IF ( humidity  .OR.  passive_scalar )  THEN
440!
441!--    2D-humidity/scalar arrays
442       ALLOCATE ( qs(nysg:nyng,nxlg:nxrg),                                     &
443                  qsws(nysg:nyng,nxlg:nxrg),                                   &
444                  qswst(nysg:nyng,nxlg:nxrg) )
445
446!
447!--    3D-humidity/scalar arrays
448#if defined( __nopointer )
449       ALLOCATE( q(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
450                 q_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
451                 tq_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
452#else
453       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
454                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
455                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
456#endif
457
458!
459!--    3D-arrays needed for humidity only
460       IF ( humidity )  THEN
461#if defined( __nopointer )
462          ALLOCATE( vpt(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
463#else
464          ALLOCATE( vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
465#endif
466
467          IF ( cloud_physics )  THEN
468
469!
470!--          Liquid water content
471#if defined( __nopointer )
472             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
473#else
474             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
475#endif
476!
477!--          Precipitation amount and rate (only needed if output is switched)
478             ALLOCATE( precipitation_amount(nysg:nyng,nxlg:nxrg),              &
479                       precipitation_rate(nysg:nyng,nxlg:nxrg) )
480
481!
482!--          1D-arrays
483             ALLOCATE ( nc_1d(nzb:nzt+1), pt_1d(nzb:nzt+1),                 &
484                        q_1d(nzb:nzt+1), qc_1d(nzb:nzt+1) )
485!
486!--          3D-cloud water content
487#if defined( __nopointer )
488             ALLOCATE( qc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
489#else
490             ALLOCATE( qc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
491#endif
492!
493!--          3d-precipitation rate
494             ALLOCATE( prr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
495
496             IF ( microphysics_seifert )  THEN
497!
498!--             1D-arrays
499                ALLOCATE ( nr_1d(nzb:nzt+1), qr_1d(nzb:nzt+1) )
500
501!
502!--             2D-rain water content and rain drop concentration arrays
503                ALLOCATE ( qrs(nysg:nyng,nxlg:nxrg),                        &
504                           qrsws(nysg:nyng,nxlg:nxrg),                      &
505                           qrswst(nysg:nyng,nxlg:nxrg),                     &
506                           nrs(nysg:nyng,nxlg:nxrg),                        &
507                           nrsws(nysg:nyng,nxlg:nxrg),                      &
508                           nrswst(nysg:nyng,nxlg:nxrg) )
509!
510!--             3D-rain water content, rain drop concentration arrays
511#if defined( __nopointer )
512                ALLOCATE( nr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
513                          nr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
514                          qr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
515                          qr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
516                          tnr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),             &
517                          tqr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
518#else
519                ALLOCATE( nr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
520                          nr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
521                          nr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
522                          qr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
523                          qr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
524                          qr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
525#endif
526             ENDIF
527
528          ENDIF
529
530          IF ( cloud_droplets )  THEN
531!
532!--          Liquid water content, change in liquid water content
533#if defined( __nopointer )
534             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
535                        ql_c(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
536#else
537             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
538                        ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
539#endif
540!
541!--          Real volume of particles (with weighting), volume of particles
542             ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
543                        ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
544          ENDIF
545
546       ENDIF
547
548    ENDIF
549
550    IF ( ocean )  THEN
551       ALLOCATE( saswsb(nysg:nyng,nxlg:nxrg),                                  &
552                 saswst(nysg:nyng,nxlg:nxrg) )
553#if defined( __nopointer )
554       ALLOCATE( prho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
555                 rho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
556                 sa(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
557                 sa_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
558                 tsa_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
559#else
560       ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
561                 rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                         &
562                 sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
563                 sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
564                 sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
565       prho => prho_1
566       rho  => rho_1  ! routines calc_mean_profile and diffusion_e require
567                      ! density to be apointer
568#endif
569       IF ( humidity_remote )  THEN
570          ALLOCATE( qswst_remote(nysg:nyng,nxlg:nxrg))
571          qswst_remote = 0.0_wp
572       ENDIF
573    ENDIF
574
575!
576!-- 3D-array for storing the dissipation, needed for calculating the sgs
577!-- particle velocities
578    IF ( use_sgs_for_particles  .OR.  wang_kernel  .OR.  collision_turbulence  &
579         .OR.  num_acc_per_node > 0 )  THEN
580       ALLOCATE( diss(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
581    ENDIF
582
583!
584!-- 1D-array for large scale subsidence velocity
585    IF ( .NOT. ALLOCATED( w_subs ) )  THEN
586       ALLOCATE ( w_subs(nzb:nzt+1) )
587       w_subs = 0.0_wp
588    ENDIF
589
590!
591!-- ID-array and state-space-array for the parallel random number generator
592    IF ( random_generator == 'random-parallel' )  THEN
593       ALLOCATE ( seq_random_array(5,nysg:nyng,nxlg:nxrg) )
594       ALLOCATE ( id_random_array(0:ny,0:nx) )
595       seq_random_array = 0
596       id_random_array  = 0
597    ENDIF
598   
599!
600!-- 4D-array for storing the Rif-values at vertical walls
601    IF ( topography /= 'flat' )  THEN
602       ALLOCATE( rif_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg,1:4) )
603       rif_wall = 0.0_wp
604    ENDIF
605
606!
607!-- Arrays to store velocity data from t-dt and the phase speeds which
608!-- are needed for radiation boundary conditions
609    IF ( outflow_l )  THEN
610       ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2),                               &
611                 v_m_l(nzb:nzt+1,nysg:nyng,0:1),                               &
612                 w_m_l(nzb:nzt+1,nysg:nyng,0:1) )
613    ENDIF
614    IF ( outflow_r )  THEN
615       ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
616                 v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
617                 w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) )
618    ENDIF
619    IF ( outflow_l  .OR.  outflow_r )  THEN
620       ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng),           &
621                 c_w(nzb:nzt+1,nysg:nyng) )
622    ENDIF
623    IF ( outflow_s )  THEN
624       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg),                               &
625                 v_m_s(nzb:nzt+1,1:2,nxlg:nxrg),                               &
626                 w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) )
627    ENDIF
628    IF ( outflow_n )  THEN
629       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
630                 v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
631                 w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) )
632    ENDIF
633    IF ( outflow_s  .OR.  outflow_n )  THEN
634       ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg),           &
635                 c_w(nzb:nzt+1,nxlg:nxrg) )
636    ENDIF
637    IF ( outflow_l  .OR.  outflow_r  .OR.  outflow_s  .OR.  outflow_n )  THEN
638       ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) )                   
639       ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) )
640    ENDIF
641
642
643#if ! defined( __nopointer )
644!
645!-- Initial assignment of the pointers
646    e  => e_1;   e_p  => e_2;   te_m  => e_3
647    IF ( .NOT. neutral )  THEN
648       pt => pt_1;  pt_p => pt_2;  tpt_m => pt_3
649    ELSE
650       pt => pt_1;  pt_p => pt_1;  tpt_m => pt_3
651    ENDIF
652    u  => u_1;   u_p  => u_2;   tu_m  => u_3
653    v  => v_1;   v_p  => v_2;   tv_m  => v_3
654    w  => w_1;   w_p  => w_2;   tw_m  => w_3
655
656    IF ( humidity  .OR.  passive_scalar )  THEN
657       q => q_1;  q_p => q_2;  tq_m => q_3
658       IF ( humidity )  THEN
659          vpt  => vpt_1   
660          IF ( cloud_physics )  THEN
661             ql => ql_1
662             qc => qc_1
663             IF ( microphysics_seifert )  THEN
664                qr => qr_1;  qr_p  => qr_2;  tqr_m  => qr_3
665                nr => nr_1;  nr_p  => nr_2;  tnr_m  => nr_3
666             ENDIF
667          ENDIF
668       ENDIF
669       IF ( cloud_droplets )  THEN
670          ql   => ql_1
671          ql_c => ql_2
672       ENDIF
673    ENDIF
674
675    IF ( ocean )  THEN
676       sa => sa_1;  sa_p => sa_2;  tsa_m => sa_3
677    ENDIF
678#endif
679
680!
681!-- Allocate land surface model arrays
682    IF ( land_surface )  THEN
683       CALL lsm_init_arrays
684    ENDIF
685
686!
687!-- Allocate arrays containing the RK coefficient for calculation of
688!-- perturbation pressure and turbulent fluxes. At this point values are
689!-- set for pressure calculation during initialization (where no timestep
690!-- is done). Further below the values needed within the timestep scheme
691!-- will be set.
692!-- TO_DO: zeroth element added to weight_pres because in case of nesting
693!--        pres may need to be called  outside the RK-substeps. Antti will
694!--        check if this is really required.
695    ALLOCATE( weight_substep(1:intermediate_timestep_count_max),               &
696              weight_pres(0:intermediate_timestep_count_max) )
697    weight_substep = 1.0_wp
698    weight_pres    = 1.0_wp
699    intermediate_timestep_count = 1  ! needed when simulated_time = 0.0
700       
701    CALL location_message( 'finished', .TRUE. )
702!
703!-- Initialize model variables
704    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
705         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
706!
707!--    First model run of a possible job queue.
708!--    Initial profiles of the variables must be computes.
709       IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
710
711          CALL location_message( 'initializing with 1D model profiles', .FALSE. )
712!
713!--       Use solutions of the 1D model as initial profiles,
714!--       start 1D model
715          CALL init_1d_model
716!
717!--       Transfer initial profiles to the arrays of the 3D model
718          DO  i = nxlg, nxrg
719             DO  j = nysg, nyng
720                e(:,j,i)  = e1d
721                kh(:,j,i) = kh1d
722                km(:,j,i) = km1d
723                pt(:,j,i) = pt_init
724                u(:,j,i)  = u1d
725                v(:,j,i)  = v1d
726             ENDDO
727          ENDDO
728
729          IF ( humidity  .OR.  passive_scalar )  THEN
730             DO  i = nxlg, nxrg
731                DO  j = nysg, nyng
732                   q(:,j,i) = q_init
733                ENDDO
734             ENDDO
735             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
736                DO  i = nxlg, nxrg
737                   DO  j = nysg, nyng
738                      qr(:,j,i) = 0.0_wp
739                      nr(:,j,i) = 0.0_wp
740                   ENDDO
741                ENDDO
742!
743!--             Initialze nc_1d with default value
744                nc_1d(:) = nc_const
745
746             ENDIF
747          ENDIF
748
749          IF ( .NOT. constant_diffusion )  THEN
750             DO  i = nxlg, nxrg
751                DO  j = nysg, nyng
752                   e(:,j,i)  = e1d
753                ENDDO
754             ENDDO
755!
756!--          Store initial profiles for output purposes etc.
757             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
758
759             IF ( constant_flux_layer )  THEN
760                ol   = ( zu(nzb+1) - zw(nzb) ) / ( rif1d(nzb+1) + 1.0E-20_wp )
761                ts   = 0.0_wp  ! could actually be computed more accurately in the
762                               ! 1D model. Update when opportunity arises.
763                us   = us1d
764                usws = usws1d
765                vsws = vsws1d
766             ELSE
767                ts   = 0.0_wp  ! must be set, because used in
768                ol   = ( zu(nzb+1) - zw(nzb) ) / zeta_min  ! flowste
769                us   = 0.0_wp
770                usws = 0.0_wp
771                vsws = 0.0_wp
772             ENDIF
773
774          ELSE
775             e    = 0.0_wp  ! must be set, because used in
776             ol   = ( zu(nzb+1) - zw(nzb) ) / zeta_min  ! flowste
777             ts   = 0.0_wp
778             us   = 0.0_wp
779             usws = 0.0_wp
780             vsws = 0.0_wp
781          ENDIF
782          uswst = top_momentumflux_u
783          vswst = top_momentumflux_v
784
785!
786!--       In every case qs = 0.0 (see also pt)
787!--       This could actually be computed more accurately in the 1D model.
788!--       Update when opportunity arises!
789          IF ( humidity  .OR.  passive_scalar )  THEN
790             qs = 0.0_wp
791             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
792                qrs = 0.0_wp
793                nrs = 0.0_wp
794             ENDIF
795          ENDIF
796
797!
798!--       Inside buildings set velocities back to zero
799          IF ( topography /= 'flat' )  THEN
800             DO  i = nxlg, nxrg
801                DO  j = nysg, nyng
802                   u(nzb:nzb_u_inner(j,i),j,i) = 0.0_wp
803                   v(nzb:nzb_v_inner(j,i),j,i) = 0.0_wp
804                ENDDO
805             ENDDO
806             
807!
808!--          WARNING: The extra boundary conditions set after running the
809!--          -------  1D model impose an error on the divergence one layer
810!--                   below the topography; need to correct later
811!--          ATTENTION: Provisional correction for Piacsek & Williams
812!--          ---------  advection scheme: keep u and v zero one layer below
813!--                     the topography.
814             IF ( ibc_uv_b == 1 )  THEN
815!
816!--             Neumann condition
817                DO  i = nxl-1, nxr+1
818                   DO  j = nys-1, nyn+1
819                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i)
820                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i)
821                   ENDDO
822                ENDDO
823
824             ENDIF
825
826          ENDIF
827
828          CALL location_message( 'finished', .TRUE. )
829
830       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )    &
831       THEN
832
833          CALL location_message( 'initializing with constant profiles', .FALSE. )
834!
835!--       Overwrite initial profiles in case of nudging
836          IF ( nudging )  THEN
837             pt_init = ptnudge(:,1)
838             u_init  = unudge(:,1)
839             v_init  = vnudge(:,1)
840             IF ( humidity  .OR.  passive_scalar )  THEN
841                q_init = qnudge(:,1)
842             ENDIF
843
844             WRITE( message_string, * ) 'Initial profiles of u, v and ',       &
845                 'scalars from NUDGING_DATA are used.'
846             CALL message( 'init_3d_model', 'PA0370', 0, 0, 0, 6, 0 )
847          ENDIF
848
849!
850!--       Use constructed initial profiles (velocity constant with height,
851!--       temperature profile with constant gradient)
852          DO  i = nxlg, nxrg
853             DO  j = nysg, nyng
854                pt(:,j,i) = pt_init
855                u(:,j,i)  = u_init
856                v(:,j,i)  = v_init
857             ENDDO
858          ENDDO
859
860!
861!--       Set initial horizontal velocities at the lowest computational grid
862!--       levels to zero in order to avoid too small time steps caused by the
863!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
864!--       in the limiting formula!).
865          IF ( ibc_uv_b /= 1 )  THEN
866             DO  i = nxlg, nxrg
867                DO  j = nysg, nyng
868                   u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0_wp
869                   v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0_wp
870                ENDDO
871             ENDDO
872          ENDIF
873
874          IF ( humidity  .OR.  passive_scalar )  THEN
875             DO  i = nxlg, nxrg
876                DO  j = nysg, nyng
877                   q(:,j,i) = q_init
878                ENDDO
879             ENDDO
880             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
881!
882!--             Initialze nc_1d with default value
883                nc_1d(:) = nc_const
884
885                DO  i = nxlg, nxrg
886                   DO  j = nysg, nyng
887                      qr(:,j,i) = 0.0_wp
888                      nr(:,j,i) = 0.0_wp
889                   ENDDO
890                ENDDO
891
892             ENDIF
893          ENDIF
894
895          IF ( ocean )  THEN
896             DO  i = nxlg, nxrg
897                DO  j = nysg, nyng
898                   sa(:,j,i) = sa_init
899                ENDDO
900             ENDDO
901          ENDIF
902         
903          IF ( constant_diffusion )  THEN
904             km   = km_constant
905             kh   = km / prandtl_number
906             e    = 0.0_wp
907          ELSEIF ( e_init > 0.0_wp )  THEN
908             DO  k = nzb+1, nzt
909                km(k,:,:) = 0.1_wp * l_grid(k) * SQRT( e_init )
910             ENDDO
911             km(nzb,:,:)   = km(nzb+1,:,:)
912             km(nzt+1,:,:) = km(nzt,:,:)
913             kh   = km / prandtl_number
914             e    = e_init
915          ELSE
916             IF ( .NOT. ocean )  THEN
917                kh   = 0.01_wp   ! there must exist an initial diffusion, because
918                km   = 0.01_wp   ! otherwise no TKE would be produced by the
919                              ! production terms, as long as not yet
920                              ! e = (u*/cm)**2 at k=nzb+1
921             ELSE
922                kh   = 0.00001_wp
923                km   = 0.00001_wp
924             ENDIF
925             e    = 0.0_wp
926          ENDIF
927          ol    = ( zu(nzb+1) - zw(nzb) ) / zeta_min
928          ts    = 0.0_wp
929          us    = 0.0_wp
930          usws  = 0.0_wp
931          uswst = top_momentumflux_u
932          vsws  = 0.0_wp
933          vswst = top_momentumflux_v
934          IF ( humidity  .OR.  passive_scalar ) qs = 0.0_wp
935
936!
937!--       Compute initial temperature field and other constants used in case
938!--       of a sloping surface
939          IF ( sloping_surface )  CALL init_slope
940
941          CALL location_message( 'finished', .TRUE. )
942
943       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )                  &
944       THEN
945
946          CALL location_message( 'initializing by user', .FALSE. )
947!
948!--       Initialization will completely be done by the user
949          CALL user_init_3d_model
950
951          CALL location_message( 'finished', .TRUE. )
952
953       ENDIF
954
955       CALL location_message( 'initializing statistics, boundary conditions, etc.', &
956                              .FALSE. )
957
958!
959!--    Bottom boundary
960       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
961          u(nzb,:,:) = 0.0_wp
962          v(nzb,:,:) = 0.0_wp
963       ENDIF
964
965!
966!--    Apply channel flow boundary condition
967       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
968          u(nzt+1,:,:) = 0.0_wp
969          v(nzt+1,:,:) = 0.0_wp
970       ENDIF
971
972!
973!--    Calculate virtual potential temperature
974       IF ( humidity ) vpt = pt * ( 1.0_wp + 0.61_wp * q )
975
976!
977!--    Store initial profiles for output purposes etc.
978       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
979       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
980       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
981          hom(nzb,1,5,:) = 0.0_wp
982          hom(nzb,1,6,:) = 0.0_wp
983       ENDIF
984       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
985       hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 )
986       hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 )
987
988       IF ( ocean )  THEN
989!
990!--       Store initial salinity profile
991          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
992       ENDIF
993
994       IF ( humidity )  THEN
995!
996!--       Store initial profile of total water content, virtual potential
997!--       temperature
998          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
999          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
1000          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
1001!
1002!--          Store initial profile of specific humidity and potential
1003!--          temperature
1004             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1005             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1006          ENDIF
1007       ENDIF
1008
1009       IF ( passive_scalar )  THEN
1010!
1011!--       Store initial scalar profile
1012          hom(:,1,26,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1013       ENDIF
1014
1015!
1016!--    Initialize the random number generators (from numerical recipes)
1017       CALL random_function_ini
1018       
1019       IF ( random_generator == 'random-parallel' )  THEN
1020!--       Asigning an ID to every vertical gridpoint column
1021!--       dependig on the ensemble run number.
1022          random_dummy=1
1023          DO j=0,ny
1024             DO i=0,nx
1025                id_random_array(j,i) = random_dummy + 1E6                      &
1026                                       * ( ensemble_member_nr - 1000 )
1027                random_dummy = random_dummy + 1
1028             END DO
1029          ENDDO
1030!--       Initializing with random_seed_parallel for every vertical
1031!--       gridpoint column.
1032          random_dummy=0
1033          DO j = nysg, nyng
1034             DO i = nxlg, nxrg
1035                CALL random_seed_parallel (random_sequence=id_random_array(j, i))
1036                CALL random_number_parallel (random_dummy)
1037                CALL random_seed_parallel (get=seq_random_array(:, j, i))
1038             END DO
1039          ENDDO
1040       ENDIF
1041
1042!
1043!--    Initialize fluxes at bottom surface
1044       IF ( use_surface_fluxes )  THEN
1045
1046          IF ( constant_heatflux )  THEN
1047!
1048!--          Heat flux is prescribed
1049             IF ( random_heatflux )  THEN
1050                CALL disturb_heatflux
1051             ELSE
1052                shf = surface_heatflux
1053!
1054!--             Initialize shf with data from external file LSF_DATA
1055                IF ( large_scale_forcing .AND. lsf_surf )  THEN
1056                   CALL ls_forcing_surf ( simulated_time )
1057                ENDIF
1058
1059!
1060!--             Over topography surface_heatflux is replaced by wall_heatflux(0)
1061                IF ( TRIM( topography ) /= 'flat' )  THEN
1062                   DO  i = nxlg, nxrg
1063                      DO  j = nysg, nyng
1064                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
1065                            shf(j,i) = wall_heatflux(0)
1066                         ENDIF
1067                      ENDDO
1068                   ENDDO
1069                ENDIF
1070             ENDIF
1071          ENDIF
1072
1073!
1074!--       Determine the near-surface water flux
1075          IF ( humidity  .OR.  passive_scalar )  THEN
1076             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1077                qrsws = 0.0_wp
1078                nrsws = 0.0_wp
1079             ENDIF
1080             IF ( constant_waterflux )  THEN
1081                qsws   = surface_waterflux
1082!
1083!--             Over topography surface_waterflux is replaced by
1084!--             wall_humidityflux(0)
1085                IF ( TRIM( topography ) /= 'flat' )  THEN
1086                   wall_qflux = wall_humidityflux
1087                   DO  i = nxlg, nxrg
1088                      DO  j = nysg, nyng
1089                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
1090                            qsws(j,i) = wall_qflux(0)
1091                         ENDIF
1092                      ENDDO
1093                   ENDDO
1094                ENDIF
1095             ENDIF
1096          ENDIF
1097
1098       ENDIF
1099
1100!
1101!--    Initialize fluxes at top surface
1102!--    Currently, only the heatflux and salinity flux can be prescribed.
1103!--    The latent flux is zero in this case!
1104       IF ( use_top_fluxes )  THEN
1105
1106          IF ( constant_top_heatflux )  THEN
1107!
1108!--          Heat flux is prescribed
1109             tswst = top_heatflux
1110
1111             IF ( humidity  .OR.  passive_scalar )  THEN
1112                qswst = 0.0_wp
1113                IF ( cloud_physics  .AND.  microphysics_seifert ) 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.  microphysics_seifert )  THEN
1154                qrsws = 0.0_wp
1155                nrsws = 0.0_wp
1156             ENDIF
1157          ENDIF
1158
1159       ENDIF
1160
1161!
1162!--    Set the reference state to be used in the buoyancy terms (for ocean runs
1163!--    the reference state will be set (overwritten) in init_ocean)
1164       IF ( use_single_reference_value )  THEN
1165          IF (  .NOT.  humidity )  THEN
1166             ref_state(:) = pt_reference
1167          ELSE
1168             ref_state(:) = vpt_reference
1169          ENDIF
1170       ELSE
1171          IF (  .NOT.  humidity )  THEN
1172             ref_state(:) = pt_init(:)
1173          ELSE
1174             ref_state(:) = vpt(:,nys,nxl)
1175          ENDIF
1176       ENDIF
1177
1178!
1179!--    For the moment, vertical velocity is zero
1180       w = 0.0_wp
1181
1182!
1183!--    Initialize array sums (must be defined in first call of pres)
1184       sums = 0.0_wp
1185
1186!
1187!--    In case of iterative solvers, p must get an initial value
1188       IF ( psolver(1:9) == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0_wp
1189
1190!
1191!--    Treating cloud physics, liquid water content and precipitation amount
1192!--    are zero at beginning of the simulation
1193       IF ( cloud_physics )  THEN
1194          ql = 0.0_wp
1195          qc = 0.0_wp
1196
1197          precipitation_amount = 0.0_wp
1198
1199          IF ( microphysics_seifert )  THEN
1200             nc_1d = nc_const
1201          ENDIF
1202       ENDIF
1203!
1204!--    Impose vortex with vertical axis on the initial velocity profile
1205       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
1206          CALL init_rankine
1207       ENDIF
1208
1209!
1210!--    Impose temperature anomaly (advection test only)
1211       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
1212          CALL init_pt_anomaly
1213       ENDIF
1214
1215!
1216!--    If required, change the surface temperature at the start of the 3D run
1217       IF ( pt_surface_initial_change /= 0.0_wp )  THEN
1218          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
1219       ENDIF
1220
1221!
1222!--    If required, change the surface humidity/scalar at the start of the 3D
1223!--    run
1224       IF ( ( humidity .OR. passive_scalar ) .AND.                             &
1225            q_surface_initial_change /= 0.0_wp )  THEN
1226          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
1227       ENDIF
1228
1229!
1230!--    Initialize old and new time levels.
1231       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
1232       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
1233
1234       IF ( humidity  .OR.  passive_scalar )  THEN
1235          tq_m = 0.0_wp
1236          q_p = q
1237          IF ( cloud_physics  .AND.  microphysics_seifert )  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.  microphysics_seifert )  THEN
1423             qr_p = qr
1424             nr_p = nr
1425          ENDIF
1426       ENDIF
1427       IF ( ocean )  sa_p = sa
1428
1429!
1430!--    Allthough tendency arrays are set in prognostic_equations, they have
1431!--    have to be predefined here because they are used (but multiplied with 0)
1432!--    there before they are set.
1433       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
1434       IF ( humidity  .OR.  passive_scalar )  THEN
1435          tq_m = 0.0_wp
1436          IF ( cloud_physics  .AND.  microphysics_seifert )  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    CALL location_message( 'initializing plant canopy model', .FALSE. )
1635    IF ( plant_canopy )  CALL pcm_init
1636    CALL location_message( 'finished', .TRUE. )
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 lsm_init
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 radiation_init
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.