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

Last change on this file since 2320 was 2320, checked in by suehring, 7 years ago

large-scale forcing and nudging modularized

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