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

Last change on this file since 2864 was 2864, checked in by suehring, 6 years ago

Bugfix, move call of user_init in front of initialization of grid-point arrays

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