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

Last change on this file since 3454 was 3448, checked in by kanani, 6 years ago

Implementation of human thermal indices (from branch biomet_p2 at r3444)

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