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

Last change on this file since 3570 was 3569, checked in by kanani, 6 years ago

Fix for biomet output (ticket:757), merge of uv_exposure into biometeorology_mod

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