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

Last change on this file since 3608 was 3608, checked in by suehring, 5 years ago

Bugfix in initialization of surfaces in cyclic-fill case + bugfix in radiation output

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