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

Last change on this file since 2331 was 2329, checked in by knoop, 7 years ago

Bugfix for topography usage with anelastic approximation and boussinesq approximation with air density != 1

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