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

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

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