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

Last change on this file since 2713 was 2705, checked in by maronga, 7 years ago

minor bugfix for restarts

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