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

Last change on this file since 2182 was 2173, checked in by knoop, 8 years ago

last commit documented

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