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

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

Framework to input single surface variables independent on land/urban-surface model via static input file provided. This way, input and initialization of heterogeneous roughness was already added.

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