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

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

Adjustments according new topography and surface-modelling concept implemented

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