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

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

Merge branch salsa with trunk

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