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

Last change on this file since 3946 was 3939, checked in by suehring, 6 years ago

Remove unused variable

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