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

Last change on this file since 3294 was 3294, checked in by raasch, 6 years ago

modularization of the ocean code

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