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

Last change on this file since 2604 was 2564, checked in by Giersch, 7 years ago

Bugfix: Variable wind_turbine changed to module control_parameters

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