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

Last change on this file since 1952 was 1921, checked in by suehring, 8 years ago

last commit documented

  • Property svn:keywords set to Id
  • Property svn:mergeinfo set to False
    /palm/branches/forwind/SOURCE/init_3d_model.f901564-1913
File size: 70.7 KB
Line 
1!> @file init_3d_model.f90
2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2016 Leibniz Universitaet Hannover
17!--------------------------------------------------------------------------------!
18!
19! Current revisions:
20! ------------------
21!
22!
23! Former revisions:
24! -----------------
25! $Id: init_3d_model.f90 1921 2016-05-30 10:52:07Z suehring $
26!
27! 1920 2016-05-30 10:50:15Z suehring
28! Initialize us with very small number to avoid segmentation fault during
29! calculation of Obukhov length
30!
31! 1918 2016-05-27 14:35:57Z raasch
32! intermediate_timestep_count is set 0 instead 1 for first call of pres,
33! bugfix: initialization of local sum arrays are moved to the beginning of the
34!         routine because otherwise results from pres are overwritten
35!
36! 1914 2016-05-26 14:44:07Z witha
37! Added initialization of the wind turbine model
38!
39! 1878 2016-04-19 12:30:36Z hellstea
40! The zeroth element of weight_pres removed as unnecessary
41!
42! 1849 2016-04-08 11:33:18Z hoffmann
43! Adapted for modularization of microphysics.
44! precipitation_amount, precipitation_rate, prr moved to arrays_3d.
45! Initialization of nc_1d, nr_1d, pt_1d, qc_1d, qr_1d, q_1d moved to
46! microphysics_init.
47!
48! 1845 2016-04-08 08:29:13Z raasch
49! nzb_2d replaced by nzb_u|v_inner
50!
51! 1833 2016-04-07 14:23:03Z raasch
52! initialization of spectra quantities moved to spectra_mod
53!
54! 1831 2016-04-07 13:15:51Z hoffmann
55! turbulence renamed collision_turbulence
56!
57! 1826 2016-04-07 12:01:39Z maronga
58! Renamed radiation calls.
59! Renamed canopy model calls.
60!
61! 1822 2016-04-07 07:49:42Z hoffmann
62! icloud_scheme replaced by microphysics_*
63!
64! 1817 2016-04-06 15:44:20Z maronga
65! Renamed lsm calls.
66!
67! 1815 2016-04-06 13:49:59Z raasch
68! zero-settings for velocities inside topography re-activated (was deactivated
69! in r1762)
70!
71! 1788 2016-03-10 11:01:04Z maronga
72! Added z0q.
73! Syntax layout improved.
74!
75! 1783 2016-03-06 18:36:17Z raasch
76! netcdf module name changed + related changes
77!
78! 1764 2016-02-28 12:45:19Z raasch
79! bugfix: increase size of volume_flow_area_l and volume_flow_initial_l by 1
80!
81! 1762 2016-02-25 12:31:13Z hellstea
82! Introduction of nested domain feature
83!
84! 1738 2015-12-18 13:56:05Z raasch
85! calculate mean surface level height for each statistic region
86!
87! 1734 2015-12-02 12:17:12Z raasch
88! no initial disturbances in case that the disturbance energy limit has been
89! set zero
90!
91! 1707 2015-11-02 15:24:52Z maronga
92! Bugfix: transfer of Richardson number from 1D model to Obukhov length caused
93! devision by zero in neutral stratification
94!
95! 1691 2015-10-26 16:17:44Z maronga
96! Call to init_surface_layer added. rif is replaced by ol and zeta.
97!
98! 1682 2015-10-07 23:56:08Z knoop
99! Code annotations made doxygen readable
100!
101! 1615 2015-07-08 18:49:19Z suehring
102! Enable turbulent inflow for passive_scalar and humidity
103!
104! 1585 2015-04-30 07:05:52Z maronga
105! Initialization of radiation code is now done after LSM initializtion
106!
107! 1575 2015-03-27 09:56:27Z raasch
108! adjustments for psolver-queries
109!
110! 1551 2015-03-03 14:18:16Z maronga
111! Allocation of land surface arrays is now done in the subroutine lsm_init_arrays,
112! which is part of land_surface_model.
113!
114! 1507 2014-12-10 12:14:18Z suehring
115! Bugfix: set horizontal velocity components to zero inside topography
116!
117! 1496 2014-12-02 17:25:50Z maronga
118! Added initialization of the land surface and radiation schemes
119!
120! 1484 2014-10-21 10:53:05Z kanani
121! Changes due to new module structure of the plant canopy model:
122! canopy-related initialization (e.g. lad and canopy_heat_flux) moved to new
123! subroutine init_plant_canopy within the module plant_canopy_model_mod,
124! call of subroutine init_plant_canopy added.
125!
126! 1431 2014-07-15 14:47:17Z suehring
127! var_d added, in order to normalize spectra.
128!
129! 1429 2014-07-15 12:53:45Z knoop
130! Ensemble run capability added to parallel random number generator
131!
132! 1411 2014-05-16 18:01:51Z suehring
133! Initial horizontal velocity profiles were not set to zero at the first vertical
134! grid level in case of non-cyclic lateral boundary conditions.
135!
136! 1406 2014-05-16 13:47:01Z raasch
137! bugfix: setting of initial velocities at k=1 to zero not in case of a
138! no-slip boundary condition for uv
139!
140! 1402 2014-05-09 14:25:13Z raasch
141! location messages modified
142!
143! 1400 2014-05-09 14:03:54Z knoop
144! Parallel random number generator added
145!
146! 1384 2014-05-02 14:31:06Z raasch
147! location messages added
148!
149! 1361 2014-04-16 15:17:48Z hoffmann
150! tend_* removed
151! Bugfix: w_subs is not allocated anymore if it is already allocated
152!
153! 1359 2014-04-11 17:15:14Z hoffmann
154! module lpm_init_mod added to use statements, because lpm_init has become a
155! module
156!
157! 1353 2014-04-08 15:21:23Z heinze
158! REAL constants provided with KIND-attribute
159!
160! 1340 2014-03-25 19:45:13Z kanani
161! REAL constants defined as wp-kind
162!
163! 1322 2014-03-20 16:38:49Z raasch
164! REAL constants defined as wp-kind
165! module interfaces removed
166!
167! 1320 2014-03-20 08:40:49Z raasch
168! ONLY-attribute added to USE-statements,
169! kind-parameters added to all INTEGER and REAL declaration statements,
170! kinds are defined in new module kinds,
171! revision history before 2012 removed,
172! comment fields (!:) to be used for variable explanations added to
173! all variable declaration statements
174!
175! 1316 2014-03-17 07:44:59Z heinze
176! Bugfix: allocation of w_subs
177!
178! 1299 2014-03-06 13:15:21Z heinze
179! Allocate w_subs due to extension of large scale subsidence in combination
180! with large scale forcing data (LSF_DATA)
181!
182! 1241 2013-10-30 11:36:58Z heinze
183! Overwrite initial profiles in case of nudging
184! Inititialize shf and qsws in case of large_scale_forcing
185!
186! 1221 2013-09-10 08:59:13Z raasch
187! +rflags_s_inner in copyin statement, use copyin for most arrays instead of
188! copy
189!
190! 1212 2013-08-15 08:46:27Z raasch
191! array tri is allocated and included in data copy statement
192!
193! 1195 2013-07-01 12:27:57Z heinze
194! Bugfix: move allocation of ref_state to parin.f90 and read_var_list.f90
195!
196! 1179 2013-06-14 05:57:58Z raasch
197! allocate and set ref_state to be used in buoyancy terms
198!
199! 1171 2013-05-30 11:27:45Z raasch
200! diss array is allocated with full size if accelerator boards are used
201!
202! 1159 2013-05-21 11:58:22Z fricke
203! -bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir
204!
205! 1153 2013-05-10 14:33:08Z raasch
206! diss array is allocated with dummy elements even if it is not needed
207! (required by PGI 13.4 / CUDA 5.0)
208!
209! 1115 2013-03-26 18:16:16Z hoffmann
210! unused variables removed
211!
212! 1113 2013-03-10 02:48:14Z raasch
213! openACC directive modified
214!
215! 1111 2013-03-08 23:54:10Z raasch
216! openACC directives added for pres
217! array diss allocated only if required
218!
219! 1092 2013-02-02 11:24:22Z raasch
220! unused variables removed
221!
222! 1065 2012-11-22 17:42:36Z hoffmann
223! allocation of diss (dissipation rate) in case of turbulence = .TRUE. added
224!
225! 1053 2012-11-13 17:11:03Z hoffmann
226! allocation and initialisation of necessary data arrays for the two-moment
227! cloud physics scheme the two new prognostic equations (nr, qr):
228! +dr, lambda_r, mu_r, sed_*, xr, *s, *sws, *swst, *, *_p, t*_m, *_1, *_2, *_3,
229! +tend_*, prr
230!
231! 1036 2012-10-22 13:43:42Z raasch
232! code put under GPL (PALM 3.9)
233!
234! 1032 2012-10-21 13:03:21Z letzel
235! save memory by not allocating pt_2 in case of neutral = .T.
236!
237! 1025 2012-10-07 16:04:41Z letzel
238! bugfix: swap indices of mask for ghost boundaries
239!
240! 1015 2012-09-27 09:23:24Z raasch
241! mask is set to zero for ghost boundaries
242!
243! 1010 2012-09-20 07:59:54Z raasch
244! cpp switch __nopointer added for pointer free version
245!
246! 1003 2012-09-14 14:35:53Z raasch
247! nxra,nyna, nzta replaced ny nxr, nyn, nzt
248!
249! 1001 2012-09-13 14:08:46Z raasch
250! all actions concerning leapfrog scheme removed
251!
252! 996 2012-09-07 10:41:47Z raasch
253! little reformatting
254!
255! 978 2012-08-09 08:28:32Z fricke
256! outflow damping layer removed
257! roughness length for scalar quantites z0h added
258! damping zone for the potential temperatur in case of non-cyclic lateral
259! boundaries added
260! initialization of ptdf_x, ptdf_y
261! initialization of c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l
262!
263! 849 2012-03-15 10:35:09Z raasch
264! init_particles renamed lpm_init
265!
266! 825 2012-02-19 03:03:44Z raasch
267! wang_collision_kernel renamed wang_kernel
268!
269! Revision 1.1  1998/03/09 16:22:22  raasch
270! Initial revision
271!
272!
273! Description:
274! ------------
275!> Allocation of arrays and initialization of the 3D model via
276!> a) pre-run the 1D model
277!> or
278!> b) pre-set constant linear profiles
279!> or
280!> c) read values of a previous run
281!------------------------------------------------------------------------------!
282 SUBROUTINE init_3d_model
283 
284
285    USE advec_ws
286
287    USE arrays_3d
288
289    USE constants,                                                             &
290        ONLY:  pi
291   
292    USE control_parameters
293   
294    USE grid_variables,                                                        &
295        ONLY:  dx, dy
296   
297    USE indices
298
299    USE lpm_init_mod,                                                          &
300        ONLY:  lpm_init
301   
302    USE kinds
303
304    USE land_surface_model_mod,                                                &
305        ONLY:  lsm_init, lsm_init_arrays, land_surface
306 
307    USE ls_forcing_mod
308
309    USE microphysics_mod,                                                      &
310        ONLY:  collision_turbulence, microphysics_init
311
312    USE model_1d,                                                              &
313        ONLY:  e1d, kh1d, km1d, l1d, rif1d, u1d, us1d, usws1d, v1d, vsws1d 
314   
315    USE netcdf_interface,                                                      &
316        ONLY:  dots_max, dots_num
317   
318    USE particle_attributes,                                                   &
319        ONLY:  particle_advection, use_sgs_for_particles, wang_kernel
320   
321    USE pegrid
322   
323    USE plant_canopy_model_mod,                                                &
324        ONLY:  pcm_init, plant_canopy
325
326    USE radiation_model_mod,                                                   &
327        ONLY:  radiation_init, radiation
328   
329    USE random_function_mod 
330   
331    USE random_generator_parallel,                                             &
332        ONLY:  random_number_parallel, random_seed_parallel, random_dummy,     &
333               id_random_array, seq_random_array
334   
335    USE statistics,                                                            &
336        ONLY:  hom, hom_sum, mean_surface_level_height, pr_palm, rmask,        &
337               statistic_regions, sums, sums_divnew_l, sums_divold_l, sums_l,  &
338               sums_l_l, sums_up_fraction_l, sums_wsts_bc_l, ts_value,         &
339               weight_pres, weight_substep
340 
341    USE surface_layer_fluxes_mod,                                              &
342        ONLY:  init_surface_layer_fluxes
343   
344    USE transpose_indices 
345
346    USE wind_turbine_model_mod,                                                &
347        ONLY:  wtm_init, wtm_init_arrays, wind_turbine
348
349    IMPLICIT NONE
350
351    INTEGER(iwp) ::  i             !<
352    INTEGER(iwp) ::  ind_array(1)  !<
353    INTEGER(iwp) ::  j             !<
354    INTEGER(iwp) ::  k             !<
355    INTEGER(iwp) ::  sr            !<
356
357    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  ngp_2dh_l  !<
358
359    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer_l    !<
360    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner_l  !<
361
362    REAL(wp), DIMENSION(1:3) ::  volume_flow_area_l     !<
363    REAL(wp), DIMENSION(1:3) ::  volume_flow_initial_l  !<
364
365    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height_l    !<
366    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_l    !<
367    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner_tmp  !<
368
369
370    CALL location_message( 'allocating arrays', .FALSE. )
371!
372!-- Allocate arrays
373    ALLOCATE( mean_surface_level_height(0:statistic_regions),                  &
374              mean_surface_level_height_l(0:statistic_regions),                &
375              ngp_2dh(0:statistic_regions), ngp_2dh_l(0:statistic_regions),    &
376              ngp_3d(0:statistic_regions),                                     &
377              ngp_3d_inner(0:statistic_regions),                               &
378              ngp_3d_inner_l(0:statistic_regions),                             &
379              ngp_3d_inner_tmp(0:statistic_regions),                           &
380              sums_divnew_l(0:statistic_regions),                              &
381              sums_divold_l(0:statistic_regions) )
382    ALLOCATE( dp_smooth_factor(nzb:nzt), rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt) )
383    ALLOCATE( ngp_2dh_outer(nzb:nzt+1,0:statistic_regions),                    &
384              ngp_2dh_outer_l(nzb:nzt+1,0:statistic_regions),                  &
385              ngp_2dh_s_inner(nzb:nzt+1,0:statistic_regions),                  &
386              ngp_2dh_s_inner_l(nzb:nzt+1,0:statistic_regions),                &
387              rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions),                  &
388              sums(nzb:nzt+1,pr_palm+max_pr_user),                             &
389              sums_l(nzb:nzt+1,pr_palm+max_pr_user,0:threads_per_task-1),      &
390              sums_l_l(nzb:nzt+1,0:statistic_regions,0:threads_per_task-1),    &
391              sums_up_fraction_l(10,3,0:statistic_regions),                    &
392              sums_wsts_bc_l(nzb:nzt+1,0:statistic_regions),                   &
393              ts_value(dots_max,0:statistic_regions) )
394    ALLOCATE( ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng) )
395
396    ALLOCATE( ol(nysg:nyng,nxlg:nxrg), shf(nysg:nyng,nxlg:nxrg),               &
397              ts(nysg:nyng,nxlg:nxrg), tswst(nysg:nyng,nxlg:nxrg),             &
398              us(nysg:nyng,nxlg:nxrg), usws(nysg:nyng,nxlg:nxrg),              &
399              uswst(nysg:nyng,nxlg:nxrg), vsws(nysg:nyng,nxlg:nxrg),           &
400              vswst(nysg:nyng,nxlg:nxrg), z0(nysg:nyng,nxlg:nxrg),             &
401              z0h(nysg:nyng,nxlg:nxrg), z0q(nysg:nyng,nxlg:nxrg) )
402
403    ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr),                                    &
404              kh(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
405              km(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
406              p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
407              tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
408
409#if defined( __nopointer )
410    ALLOCATE( e(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
411              e_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
412              pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                               &
413              pt_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
414              u(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
415              u_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
416              v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
417              v_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
418              w(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                                &
419              w_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
420              te_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
421              tpt_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
422              tu_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
423              tv_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
424              tw_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
425#else
426    ALLOCATE( e_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
427              e_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
428              e_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
429              pt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
430              pt_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
431              u_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
432              u_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
433              u_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
434              v_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
435              v_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
436              v_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
437              w_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
438              w_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                              &
439              w_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
440    IF (  .NOT.  neutral )  THEN
441       ALLOCATE( pt_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
442    ENDIF
443#endif
444
445!
446!-- Following array is required for perturbation pressure within the iterative
447!-- pressure solvers. For the multistep schemes (Runge-Kutta), array p holds
448!-- the weighted average of the substeps and cannot be used in the Poisson
449!-- solver.
450    IF ( psolver == 'sor' )  THEN
451       ALLOCATE( p_loc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
452    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
453!
454!--    For performance reasons, multigrid is using one ghost layer only
455       ALLOCATE( p_loc(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) )
456    ENDIF
457
458!
459!-- Array for storing constant coeffficients of the tridiagonal solver
460    IF ( psolver == 'poisfft' )  THEN
461       ALLOCATE( tri(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1,2) )
462       ALLOCATE( tric(nxl_z:nxr_z,nys_z:nyn_z,0:nz-1) )
463    ENDIF
464
465    IF ( humidity  .OR.  passive_scalar )  THEN
466!
467!--    2D-humidity/scalar arrays
468       ALLOCATE ( qs(nysg:nyng,nxlg:nxrg),                                     &
469                  qsws(nysg:nyng,nxlg:nxrg),                                   &
470                  qswst(nysg:nyng,nxlg:nxrg) )
471
472!
473!--    3D-humidity/scalar arrays
474#if defined( __nopointer )
475       ALLOCATE( q(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                             &
476                 q_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
477                 tq_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
478#else
479       ALLOCATE( q_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
480                 q_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
481                 q_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
482#endif
483
484!
485!--    3D-arrays needed for humidity only
486       IF ( humidity )  THEN
487#if defined( __nopointer )
488          ALLOCATE( vpt(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
489#else
490          ALLOCATE( vpt_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
491#endif
492
493          IF ( cloud_physics )  THEN
494
495!
496!--          Liquid water content
497#if defined( __nopointer )
498             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
499#else
500             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
501#endif
502!
503!--          Precipitation amount and rate (only needed if output is switched)
504             ALLOCATE( precipitation_amount(nysg:nyng,nxlg:nxrg),              &
505                       precipitation_rate(nysg:nyng,nxlg:nxrg) )
506
507!
508!--          3D-cloud water content
509#if defined( __nopointer )
510             ALLOCATE( qc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
511#else
512             ALLOCATE( qc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
513#endif
514!
515!--          3d-precipitation rate
516             ALLOCATE( prr(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
517
518             IF ( microphysics_seifert )  THEN
519!
520!--             2D-rain water content and rain drop concentration arrays
521                ALLOCATE ( qrs(nysg:nyng,nxlg:nxrg),                        &
522                           qrsws(nysg:nyng,nxlg:nxrg),                      &
523                           qrswst(nysg:nyng,nxlg:nxrg),                     &
524                           nrs(nysg:nyng,nxlg:nxrg),                        &
525                           nrsws(nysg:nyng,nxlg:nxrg),                      &
526                           nrswst(nysg:nyng,nxlg:nxrg) )
527!
528!--             3D-rain water content, rain drop concentration arrays
529#if defined( __nopointer )
530                ALLOCATE( nr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
531                          nr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
532                          qr(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                &
533                          qr_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
534                          tnr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg),             &
535                          tqr_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
536#else
537                ALLOCATE( nr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
538                          nr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
539                          nr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
540                          qr_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
541                          qr_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),              &
542                          qr_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
543#endif
544             ENDIF
545
546          ENDIF
547
548          IF ( cloud_droplets )  THEN
549!
550!--          Liquid water content, change in liquid water content
551#if defined( __nopointer )
552             ALLOCATE ( ql(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
553                        ql_c(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
554#else
555             ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
556                        ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
557#endif
558!
559!--          Real volume of particles (with weighting), volume of particles
560             ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                   &
561                        ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
562          ENDIF
563
564       ENDIF
565
566    ENDIF
567
568    IF ( ocean )  THEN
569       ALLOCATE( saswsb(nysg:nyng,nxlg:nxrg),                                  &
570                 saswst(nysg:nyng,nxlg:nxrg) )
571#if defined( __nopointer )
572       ALLOCATE( prho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
573                 rho(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                           &
574                 sa(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                            &
575                 sa_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
576                 tsa_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
577#else
578       ALLOCATE( prho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
579                 rho_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                         &
580                 sa_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
581                 sa_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                          &
582                 sa_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
583       prho => prho_1
584       rho  => rho_1  ! routines calc_mean_profile and diffusion_e require
585                      ! density to be apointer
586#endif
587       IF ( humidity_remote )  THEN
588          ALLOCATE( qswst_remote(nysg:nyng,nxlg:nxrg))
589          qswst_remote = 0.0_wp
590       ENDIF
591    ENDIF
592
593!
594!-- 3D-array for storing the dissipation, needed for calculating the sgs
595!-- particle velocities
596    IF ( use_sgs_for_particles  .OR.  wang_kernel  .OR.  collision_turbulence  &
597         .OR.  num_acc_per_node > 0 )  THEN
598       ALLOCATE( diss(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
599    ENDIF
600
601!
602!-- 1D-array for large scale subsidence velocity
603    IF ( .NOT. ALLOCATED( w_subs ) )  THEN
604       ALLOCATE ( w_subs(nzb:nzt+1) )
605       w_subs = 0.0_wp
606    ENDIF
607
608!
609!-- ID-array and state-space-array for the parallel random number generator
610    IF ( random_generator == 'random-parallel' )  THEN
611       ALLOCATE ( seq_random_array(5,nysg:nyng,nxlg:nxrg) )
612       ALLOCATE ( id_random_array(0:ny,0:nx) )
613       seq_random_array = 0
614       id_random_array  = 0
615    ENDIF
616   
617!
618!-- 4D-array for storing the Rif-values at vertical walls
619    IF ( topography /= 'flat' )  THEN
620       ALLOCATE( rif_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg,1:4) )
621       rif_wall = 0.0_wp
622    ENDIF
623
624!
625!-- Arrays to store velocity data from t-dt and the phase speeds which
626!-- are needed for radiation boundary conditions
627    IF ( outflow_l )  THEN
628       ALLOCATE( u_m_l(nzb:nzt+1,nysg:nyng,1:2),                               &
629                 v_m_l(nzb:nzt+1,nysg:nyng,0:1),                               &
630                 w_m_l(nzb:nzt+1,nysg:nyng,0:1) )
631    ENDIF
632    IF ( outflow_r )  THEN
633       ALLOCATE( u_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
634                 v_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx),                           &
635                 w_m_r(nzb:nzt+1,nysg:nyng,nx-1:nx) )
636    ENDIF
637    IF ( outflow_l  .OR.  outflow_r )  THEN
638       ALLOCATE( c_u(nzb:nzt+1,nysg:nyng), c_v(nzb:nzt+1,nysg:nyng),           &
639                 c_w(nzb:nzt+1,nysg:nyng) )
640    ENDIF
641    IF ( outflow_s )  THEN
642       ALLOCATE( u_m_s(nzb:nzt+1,0:1,nxlg:nxrg),                               &
643                 v_m_s(nzb:nzt+1,1:2,nxlg:nxrg),                               &
644                 w_m_s(nzb:nzt+1,0:1,nxlg:nxrg) )
645    ENDIF
646    IF ( outflow_n )  THEN
647       ALLOCATE( u_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
648                 v_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg),                           &
649                 w_m_n(nzb:nzt+1,ny-1:ny,nxlg:nxrg) )
650    ENDIF
651    IF ( outflow_s  .OR.  outflow_n )  THEN
652       ALLOCATE( c_u(nzb:nzt+1,nxlg:nxrg), c_v(nzb:nzt+1,nxlg:nxrg),           &
653                 c_w(nzb:nzt+1,nxlg:nxrg) )
654    ENDIF
655    IF ( outflow_l  .OR.  outflow_r  .OR.  outflow_s  .OR.  outflow_n )  THEN
656       ALLOCATE( c_u_m_l(nzb:nzt+1), c_v_m_l(nzb:nzt+1), c_w_m_l(nzb:nzt+1) )                   
657       ALLOCATE( c_u_m(nzb:nzt+1), c_v_m(nzb:nzt+1), c_w_m(nzb:nzt+1) )
658    ENDIF
659
660
661#if ! defined( __nopointer )
662!
663!-- Initial assignment of the pointers
664    e  => e_1;   e_p  => e_2;   te_m  => e_3
665    IF ( .NOT. neutral )  THEN
666       pt => pt_1;  pt_p => pt_2;  tpt_m => pt_3
667    ELSE
668       pt => pt_1;  pt_p => pt_1;  tpt_m => pt_3
669    ENDIF
670    u  => u_1;   u_p  => u_2;   tu_m  => u_3
671    v  => v_1;   v_p  => v_2;   tv_m  => v_3
672    w  => w_1;   w_p  => w_2;   tw_m  => w_3
673
674    IF ( humidity  .OR.  passive_scalar )  THEN
675       q => q_1;  q_p => q_2;  tq_m => q_3
676       IF ( humidity )  THEN
677          vpt  => vpt_1   
678          IF ( cloud_physics )  THEN
679             ql => ql_1
680             qc => qc_1
681             IF ( microphysics_seifert )  THEN
682                qr => qr_1;  qr_p  => qr_2;  tqr_m  => qr_3
683                nr => nr_1;  nr_p  => nr_2;  tnr_m  => nr_3
684             ENDIF
685          ENDIF
686       ENDIF
687       IF ( cloud_droplets )  THEN
688          ql   => ql_1
689          ql_c => ql_2
690       ENDIF
691    ENDIF
692
693    IF ( ocean )  THEN
694       sa => sa_1;  sa_p => sa_2;  tsa_m => sa_3
695    ENDIF
696#endif
697
698!
699!-- Allocate land surface model arrays
700    IF ( land_surface )  THEN
701       CALL lsm_init_arrays
702    ENDIF
703
704!
705!-- Allocate wind turbine model arrays
706    IF ( wind_turbine )  THEN
707       CALL wtm_init_arrays
708    ENDIF
709
710!
711!-- Allocate arrays containing the RK coefficient for calculation of
712!-- perturbation pressure and turbulent fluxes. At this point values are
713!-- set for pressure calculation during initialization (where no timestep
714!-- is done). Further below the values needed within the timestep scheme
715!-- will be set.
716    ALLOCATE( weight_substep(1:intermediate_timestep_count_max),               &
717              weight_pres(1:intermediate_timestep_count_max) )
718    weight_substep = 1.0_wp
719    weight_pres    = 1.0_wp
720    intermediate_timestep_count = 0  ! needed when simulated_time = 0.0
721       
722    CALL location_message( 'finished', .TRUE. )
723
724!
725!-- Initialize local summation arrays for routine flow_statistics.
726!-- This is necessary because they may not yet have been initialized when they
727!-- are called from flow_statistics (or - depending on the chosen model run -
728!-- are never initialized)
729    sums_divnew_l      = 0.0_wp
730    sums_divold_l      = 0.0_wp
731    sums_l_l           = 0.0_wp
732    sums_up_fraction_l = 0.0_wp
733    sums_wsts_bc_l     = 0.0_wp
734
735
736!
737!-- Initialize model variables
738    IF ( TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
739         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
740!
741!--    First model run of a possible job queue.
742!--    Initial profiles of the variables must be computes.
743       IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 )  THEN
744
745          CALL location_message( 'initializing with 1D model profiles', .FALSE. )
746!
747!--       Use solutions of the 1D model as initial profiles,
748!--       start 1D model
749          CALL init_1d_model
750!
751!--       Transfer initial profiles to the arrays of the 3D model
752          DO  i = nxlg, nxrg
753             DO  j = nysg, nyng
754                e(:,j,i)  = e1d
755                kh(:,j,i) = kh1d
756                km(:,j,i) = km1d
757                pt(:,j,i) = pt_init
758                u(:,j,i)  = u1d
759                v(:,j,i)  = v1d
760             ENDDO
761          ENDDO
762
763          IF ( humidity  .OR.  passive_scalar )  THEN
764             DO  i = nxlg, nxrg
765                DO  j = nysg, nyng
766                   q(:,j,i) = q_init
767                ENDDO
768             ENDDO
769             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
770                DO  i = nxlg, nxrg
771                   DO  j = nysg, nyng
772                      qr(:,j,i) = 0.0_wp
773                      nr(:,j,i) = 0.0_wp
774                   ENDDO
775                ENDDO
776
777             ENDIF
778          ENDIF
779
780          IF ( .NOT. constant_diffusion )  THEN
781             DO  i = nxlg, nxrg
782                DO  j = nysg, nyng
783                   e(:,j,i)  = e1d
784                ENDDO
785             ENDDO
786!
787!--          Store initial profiles for output purposes etc.
788             hom(:,1,25,:) = SPREAD( l1d, 2, statistic_regions+1 )
789
790             IF ( constant_flux_layer )  THEN
791                ol   = ( zu(nzb+1) - zw(nzb) ) / ( rif1d(nzb+1) + 1.0E-20_wp )
792                ts   = 0.0_wp  ! could actually be computed more accurately in the
793                               ! 1D model. Update when opportunity arises.
794                us   = us1d
795                usws = usws1d
796                vsws = vsws1d
797             ELSE
798                ts   = 0.0_wp  ! must be set, because used in
799                ol   = ( zu(nzb+1) - zw(nzb) ) / zeta_min  ! flowste
800                us   = 0.0_wp
801                usws = 0.0_wp
802                vsws = 0.0_wp
803             ENDIF
804
805          ELSE
806             e    = 0.0_wp  ! must be set, because used in
807             ol   = ( zu(nzb+1) - zw(nzb) ) / zeta_min  ! flowste
808             ts   = 0.0_wp
809             us   = 0.0_wp
810             usws = 0.0_wp
811             vsws = 0.0_wp
812          ENDIF
813          uswst = top_momentumflux_u
814          vswst = top_momentumflux_v
815
816!
817!--       In every case qs = 0.0 (see also pt)
818!--       This could actually be computed more accurately in the 1D model.
819!--       Update when opportunity arises!
820          IF ( humidity  .OR.  passive_scalar )  THEN
821             qs = 0.0_wp
822             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
823                qrs = 0.0_wp
824                nrs = 0.0_wp
825             ENDIF
826          ENDIF
827
828!
829!--       Inside buildings set velocities back to zero
830          IF ( topography /= 'flat' )  THEN
831             DO  i = nxlg, nxrg
832                DO  j = nysg, nyng
833                   u(nzb:nzb_u_inner(j,i),j,i) = 0.0_wp
834                   v(nzb:nzb_v_inner(j,i),j,i) = 0.0_wp
835                ENDDO
836             ENDDO
837             
838!
839!--          WARNING: The extra boundary conditions set after running the
840!--          -------  1D model impose an error on the divergence one layer
841!--                   below the topography; need to correct later
842!--          ATTENTION: Provisional correction for Piacsek & Williams
843!--          ---------  advection scheme: keep u and v zero one layer below
844!--                     the topography.
845             IF ( ibc_uv_b == 1 )  THEN
846!
847!--             Neumann condition
848                DO  i = nxl-1, nxr+1
849                   DO  j = nys-1, nyn+1
850                      IF ( nzb_u_inner(j,i) == 0 ) u(0,j,i) = u(1,j,i)
851                      IF ( nzb_v_inner(j,i) == 0 ) v(0,j,i) = v(1,j,i)
852                   ENDDO
853                ENDDO
854
855             ENDIF
856
857          ENDIF
858
859          CALL location_message( 'finished', .TRUE. )
860
861       ELSEIF ( INDEX(initializing_actions, 'set_constant_profiles') /= 0 )    &
862       THEN
863
864          CALL location_message( 'initializing with constant profiles', .FALSE. )
865!
866!--       Overwrite initial profiles in case of nudging
867          IF ( nudging )  THEN
868             pt_init = ptnudge(:,1)
869             u_init  = unudge(:,1)
870             v_init  = vnudge(:,1)
871             IF ( humidity  .OR.  passive_scalar )  THEN
872                q_init = qnudge(:,1)
873             ENDIF
874
875             WRITE( message_string, * ) 'Initial profiles of u, v and ',       &
876                 'scalars from NUDGING_DATA are used.'
877             CALL message( 'init_3d_model', 'PA0370', 0, 0, 0, 6, 0 )
878          ENDIF
879
880!
881!--       Use constructed initial profiles (velocity constant with height,
882!--       temperature profile with constant gradient)
883          DO  i = nxlg, nxrg
884             DO  j = nysg, nyng
885                pt(:,j,i) = pt_init
886                u(:,j,i)  = u_init
887                v(:,j,i)  = v_init
888             ENDDO
889          ENDDO
890
891!
892!--       Set initial horizontal velocities at the lowest computational grid
893!--       levels to zero in order to avoid too small time steps caused by the
894!--       diffusion limit in the initial phase of a run (at k=1, dz/2 occurs
895!--       in the limiting formula!).
896          IF ( ibc_uv_b /= 1 )  THEN
897             DO  i = nxlg, nxrg
898                DO  j = nysg, nyng
899                   u(nzb:nzb_u_inner(j,i)+1,j,i) = 0.0_wp
900                   v(nzb:nzb_v_inner(j,i)+1,j,i) = 0.0_wp
901                ENDDO
902             ENDDO
903          ENDIF
904
905          IF ( humidity  .OR.  passive_scalar )  THEN
906             DO  i = nxlg, nxrg
907                DO  j = nysg, nyng
908                   q(:,j,i) = q_init
909                ENDDO
910             ENDDO
911             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
912
913                DO  i = nxlg, nxrg
914                   DO  j = nysg, nyng
915                      qr(:,j,i) = 0.0_wp
916                      nr(:,j,i) = 0.0_wp
917                   ENDDO
918                ENDDO
919
920             ENDIF
921          ENDIF
922
923          IF ( ocean )  THEN
924             DO  i = nxlg, nxrg
925                DO  j = nysg, nyng
926                   sa(:,j,i) = sa_init
927                ENDDO
928             ENDDO
929          ENDIF
930         
931          IF ( constant_diffusion )  THEN
932             km   = km_constant
933             kh   = km / prandtl_number
934             e    = 0.0_wp
935          ELSEIF ( e_init > 0.0_wp )  THEN
936             DO  k = nzb+1, nzt
937                km(k,:,:) = 0.1_wp * l_grid(k) * SQRT( e_init )
938             ENDDO
939             km(nzb,:,:)   = km(nzb+1,:,:)
940             km(nzt+1,:,:) = km(nzt,:,:)
941             kh   = km / prandtl_number
942             e    = e_init
943          ELSE
944             IF ( .NOT. ocean )  THEN
945                kh   = 0.01_wp   ! there must exist an initial diffusion, because
946                km   = 0.01_wp   ! otherwise no TKE would be produced by the
947                              ! production terms, as long as not yet
948                              ! e = (u*/cm)**2 at k=nzb+1
949             ELSE
950                kh   = 0.00001_wp
951                km   = 0.00001_wp
952             ENDIF
953             e    = 0.0_wp
954          ENDIF
955          ol    = ( zu(nzb+1) - zw(nzb) ) / zeta_min
956          ts    = 0.0_wp
957!
958!--       Very small number is required for calculation of Obukhov length
959!--       at first timestep     
960          us    = 1E-30_wp 
961          usws  = 0.0_wp
962          uswst = top_momentumflux_u
963          vsws  = 0.0_wp
964          vswst = top_momentumflux_v
965          IF ( humidity  .OR.  passive_scalar ) qs = 0.0_wp
966
967!
968!--       Compute initial temperature field and other constants used in case
969!--       of a sloping surface
970          IF ( sloping_surface )  CALL init_slope
971
972          CALL location_message( 'finished', .TRUE. )
973
974       ELSEIF ( INDEX(initializing_actions, 'by_user') /= 0 )                  &
975       THEN
976
977          CALL location_message( 'initializing by user', .FALSE. )
978!
979!--       Initialization will completely be done by the user
980          CALL user_init_3d_model
981
982          CALL location_message( 'finished', .TRUE. )
983
984       ENDIF
985
986       CALL location_message( 'initializing statistics, boundary conditions, etc.', &
987                              .FALSE. )
988
989!
990!--    Bottom boundary
991       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2  )  THEN
992          u(nzb,:,:) = 0.0_wp
993          v(nzb,:,:) = 0.0_wp
994       ENDIF
995
996!
997!--    Apply channel flow boundary condition
998       IF ( TRIM( bc_uv_t ) == 'dirichlet_0' )  THEN
999          u(nzt+1,:,:) = 0.0_wp
1000          v(nzt+1,:,:) = 0.0_wp
1001       ENDIF
1002
1003!
1004!--    Calculate virtual potential temperature
1005       IF ( humidity ) vpt = pt * ( 1.0_wp + 0.61_wp * q )
1006
1007!
1008!--    Store initial profiles for output purposes etc.
1009       hom(:,1,5,:) = SPREAD( u(:,nys,nxl), 2, statistic_regions+1 )
1010       hom(:,1,6,:) = SPREAD( v(:,nys,nxl), 2, statistic_regions+1 )
1011       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2)  THEN
1012          hom(nzb,1,5,:) = 0.0_wp
1013          hom(nzb,1,6,:) = 0.0_wp
1014       ENDIF
1015       hom(:,1,7,:)  = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1016       hom(:,1,23,:) = SPREAD( km(:,nys,nxl), 2, statistic_regions+1 )
1017       hom(:,1,24,:) = SPREAD( kh(:,nys,nxl), 2, statistic_regions+1 )
1018
1019       IF ( ocean )  THEN
1020!
1021!--       Store initial salinity profile
1022          hom(:,1,26,:)  = SPREAD( sa(:,nys,nxl), 2, statistic_regions+1 )
1023       ENDIF
1024
1025       IF ( humidity )  THEN
1026!
1027!--       Store initial profile of total water content, virtual potential
1028!--       temperature
1029          hom(:,1,26,:) = SPREAD(   q(:,nys,nxl), 2, statistic_regions+1 )
1030          hom(:,1,29,:) = SPREAD( vpt(:,nys,nxl), 2, statistic_regions+1 )
1031          IF ( cloud_physics  .OR.  cloud_droplets ) THEN
1032!
1033!--          Store initial profile of specific humidity and potential
1034!--          temperature
1035             hom(:,1,27,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1036             hom(:,1,28,:) = SPREAD( pt(:,nys,nxl), 2, statistic_regions+1 )
1037          ENDIF
1038       ENDIF
1039
1040       IF ( passive_scalar )  THEN
1041!
1042!--       Store initial scalar profile
1043          hom(:,1,26,:) = SPREAD(  q(:,nys,nxl), 2, statistic_regions+1 )
1044       ENDIF
1045
1046!
1047!--    Initialize the random number generators (from numerical recipes)
1048       CALL random_function_ini
1049       
1050       IF ( random_generator == 'random-parallel' )  THEN
1051!--       Asigning an ID to every vertical gridpoint column
1052!--       dependig on the ensemble run number.
1053          random_dummy=1
1054          DO j=0,ny
1055             DO i=0,nx
1056                id_random_array(j,i) = random_dummy + 1E6                      &
1057                                       * ( ensemble_member_nr - 1000 )
1058                random_dummy = random_dummy + 1
1059             END DO
1060          ENDDO
1061!--       Initializing with random_seed_parallel for every vertical
1062!--       gridpoint column.
1063          random_dummy=0
1064          DO j = nysg, nyng
1065             DO i = nxlg, nxrg
1066                CALL random_seed_parallel (random_sequence=id_random_array(j, i))
1067                CALL random_number_parallel (random_dummy)
1068                CALL random_seed_parallel (get=seq_random_array(:, j, i))
1069             END DO
1070          ENDDO
1071       ENDIF
1072
1073!
1074!--    Initialize fluxes at bottom surface
1075       IF ( use_surface_fluxes )  THEN
1076
1077          IF ( constant_heatflux )  THEN
1078!
1079!--          Heat flux is prescribed
1080             IF ( random_heatflux )  THEN
1081                CALL disturb_heatflux
1082             ELSE
1083                shf = surface_heatflux
1084!
1085!--             Initialize shf with data from external file LSF_DATA
1086                IF ( large_scale_forcing .AND. lsf_surf )  THEN
1087                   CALL ls_forcing_surf ( simulated_time )
1088                ENDIF
1089
1090!
1091!--             Over topography surface_heatflux is replaced by wall_heatflux(0)
1092                IF ( TRIM( topography ) /= 'flat' )  THEN
1093                   DO  i = nxlg, nxrg
1094                      DO  j = nysg, nyng
1095                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
1096                            shf(j,i) = wall_heatflux(0)
1097                         ENDIF
1098                      ENDDO
1099                   ENDDO
1100                ENDIF
1101             ENDIF
1102          ENDIF
1103
1104!
1105!--       Determine the near-surface water flux
1106          IF ( humidity  .OR.  passive_scalar )  THEN
1107             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1108                qrsws = 0.0_wp
1109                nrsws = 0.0_wp
1110             ENDIF
1111             IF ( constant_waterflux )  THEN
1112                qsws   = surface_waterflux
1113!
1114!--             Over topography surface_waterflux is replaced by
1115!--             wall_humidityflux(0)
1116                IF ( TRIM( topography ) /= 'flat' )  THEN
1117                   wall_qflux = wall_humidityflux
1118                   DO  i = nxlg, nxrg
1119                      DO  j = nysg, nyng
1120                         IF ( nzb_s_inner(j,i) /= 0 )  THEN
1121                            qsws(j,i) = wall_qflux(0)
1122                         ENDIF
1123                      ENDDO
1124                   ENDDO
1125                ENDIF
1126             ENDIF
1127          ENDIF
1128
1129       ENDIF
1130
1131!
1132!--    Initialize fluxes at top surface
1133!--    Currently, only the heatflux and salinity flux can be prescribed.
1134!--    The latent flux is zero in this case!
1135       IF ( use_top_fluxes )  THEN
1136
1137          IF ( constant_top_heatflux )  THEN
1138!
1139!--          Heat flux is prescribed
1140             tswst = top_heatflux
1141
1142             IF ( humidity  .OR.  passive_scalar )  THEN
1143                qswst = 0.0_wp
1144                IF ( cloud_physics  .AND.  microphysics_seifert ) THEN
1145                   nrswst = 0.0_wp
1146                   qrswst = 0.0_wp
1147                ENDIF
1148             ENDIF
1149
1150             IF ( ocean )  THEN
1151                saswsb = bottom_salinityflux
1152                saswst = top_salinityflux
1153             ENDIF
1154          ENDIF
1155
1156!
1157!--       Initialization in case of a coupled model run
1158          IF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
1159             tswst = 0.0_wp
1160          ENDIF
1161
1162       ENDIF
1163
1164!
1165!--    Initialize Prandtl layer quantities
1166       IF ( constant_flux_layer )  THEN
1167
1168          z0 = roughness_length
1169          z0h = z0h_factor * z0
1170          z0q = z0h_factor * z0
1171
1172          IF ( .NOT. constant_heatflux )  THEN 
1173!
1174!--          Surface temperature is prescribed. Here the heat flux cannot be
1175!--          simply estimated, because therefore ol, u* and theta* would have
1176!--          to be computed by iteration. This is why the heat flux is assumed
1177!--          to be zero before the first time step. It approaches its correct
1178!--          value in the course of the first few time steps.
1179             shf   = 0.0_wp
1180          ENDIF
1181
1182          IF ( humidity  .OR.  passive_scalar )  THEN
1183             IF (  .NOT.  constant_waterflux )  qsws   = 0.0_wp
1184             IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1185                qrsws = 0.0_wp
1186                nrsws = 0.0_wp
1187             ENDIF
1188          ENDIF
1189
1190       ENDIF
1191
1192!
1193!--    Set the reference state to be used in the buoyancy terms (for ocean runs
1194!--    the reference state will be set (overwritten) in init_ocean)
1195       IF ( use_single_reference_value )  THEN
1196          IF (  .NOT.  humidity )  THEN
1197             ref_state(:) = pt_reference
1198          ELSE
1199             ref_state(:) = vpt_reference
1200          ENDIF
1201       ELSE
1202          IF (  .NOT.  humidity )  THEN
1203             ref_state(:) = pt_init(:)
1204          ELSE
1205             ref_state(:) = vpt(:,nys,nxl)
1206          ENDIF
1207       ENDIF
1208
1209!
1210!--    For the moment, vertical velocity is zero
1211       w = 0.0_wp
1212
1213!
1214!--    Initialize array sums (must be defined in first call of pres)
1215       sums = 0.0_wp
1216
1217!
1218!--    In case of iterative solvers, p must get an initial value
1219       IF ( psolver(1:9) == 'multigrid'  .OR.  psolver == 'sor' )  p = 0.0_wp
1220
1221!
1222!--    Treating cloud physics, liquid water content and precipitation amount
1223!--    are zero at beginning of the simulation
1224       IF ( cloud_physics )  THEN
1225          ql = 0.0_wp
1226          qc = 0.0_wp
1227
1228          precipitation_amount = 0.0_wp
1229       ENDIF
1230!
1231!--    Impose vortex with vertical axis on the initial velocity profile
1232       IF ( INDEX( initializing_actions, 'initialize_vortex' ) /= 0 )  THEN
1233          CALL init_rankine
1234       ENDIF
1235
1236!
1237!--    Impose temperature anomaly (advection test only)
1238       IF ( INDEX( initializing_actions, 'initialize_ptanom' ) /= 0 )  THEN
1239          CALL init_pt_anomaly
1240       ENDIF
1241
1242!
1243!--    If required, change the surface temperature at the start of the 3D run
1244       IF ( pt_surface_initial_change /= 0.0_wp )  THEN
1245          pt(nzb,:,:) = pt(nzb,:,:) + pt_surface_initial_change
1246       ENDIF
1247
1248!
1249!--    If required, change the surface humidity/scalar at the start of the 3D
1250!--    run
1251       IF ( ( humidity .OR. passive_scalar ) .AND.                             &
1252            q_surface_initial_change /= 0.0_wp )  THEN
1253          q(nzb,:,:) = q(nzb,:,:) + q_surface_initial_change
1254       ENDIF
1255
1256!
1257!--    Initialize old and new time levels.
1258       te_m = 0.0_wp; tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1259       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
1260
1261       IF ( humidity  .OR.  passive_scalar )  THEN
1262          tq_m = 0.0_wp
1263          q_p = q
1264          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1265             tqr_m = 0.0_wp
1266             qr_p  = qr
1267             tnr_m = 0.0_wp
1268             nr_p  = nr
1269          ENDIF
1270       ENDIF
1271
1272       IF ( ocean )  THEN
1273          tsa_m = 0.0_wp
1274          sa_p  = sa
1275       ENDIF
1276       
1277       CALL location_message( 'finished', .TRUE. )
1278
1279    ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .OR.         &
1280         TRIM( initializing_actions ) == 'cyclic_fill' )                       &
1281    THEN
1282
1283       CALL location_message( 'initializing in case of restart / cyclic_fill', &
1284                              .FALSE. )
1285!
1286!--    When reading data for cyclic fill of 3D prerun data files, read
1287!--    some of the global variables from the restart file which are required
1288!--    for initializing the inflow
1289       IF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1290
1291          DO  i = 0, io_blocks-1
1292             IF ( i == io_group )  THEN
1293                CALL read_parts_of_var_list
1294                CALL close_file( 13 )
1295             ENDIF
1296#if defined( __parallel )
1297             CALL MPI_BARRIER( comm2d, ierr )
1298#endif
1299          ENDDO
1300
1301       ENDIF
1302
1303!
1304!--    Read binary data from restart file
1305       DO  i = 0, io_blocks-1
1306          IF ( i == io_group )  THEN
1307             CALL read_3d_binary
1308          ENDIF
1309#if defined( __parallel )
1310          CALL MPI_BARRIER( comm2d, ierr )
1311#endif
1312       ENDDO
1313
1314!
1315!--    Initialization of the turbulence recycling method
1316       IF ( TRIM( initializing_actions ) == 'cyclic_fill'  .AND.               &
1317            turbulent_inflow )  THEN
1318!
1319!--       First store the profiles to be used at the inflow.
1320!--       These profiles are the (temporally) and horizontally averaged vertical
1321!--       profiles from the prerun. Alternatively, prescribed profiles
1322!--       for u,v-components can be used.
1323          ALLOCATE( mean_inflow_profiles(nzb:nzt+1,6) )
1324
1325          IF ( use_prescribed_profile_data )  THEN
1326             mean_inflow_profiles(:,1) = u_init            ! u
1327             mean_inflow_profiles(:,2) = v_init            ! v
1328          ELSE
1329             mean_inflow_profiles(:,1) = hom_sum(:,1,0)    ! u
1330             mean_inflow_profiles(:,2) = hom_sum(:,2,0)    ! v
1331          ENDIF
1332          mean_inflow_profiles(:,4) = hom_sum(:,4,0)       ! pt
1333          mean_inflow_profiles(:,5) = hom_sum(:,8,0)       ! e
1334          mean_inflow_profiles(:,6) = hom_sum(:,41,0)      ! q
1335
1336!
1337!--       If necessary, adjust the horizontal flow field to the prescribed
1338!--       profiles
1339          IF ( use_prescribed_profile_data )  THEN
1340             DO  i = nxlg, nxrg
1341                DO  j = nysg, nyng
1342                   DO  k = nzb, nzt+1
1343                      u(k,j,i) = u(k,j,i) - hom_sum(k,1,0) + u_init(k)
1344                      v(k,j,i) = v(k,j,i) - hom_sum(k,2,0) + v_init(k)
1345                   ENDDO
1346                ENDDO
1347             ENDDO
1348          ENDIF
1349
1350!
1351!--       Use these mean profiles at the inflow (provided that Dirichlet
1352!--       conditions are used)
1353          IF ( inflow_l )  THEN
1354             DO  j = nysg, nyng
1355                DO  k = nzb, nzt+1
1356                   u(k,j,nxlg:-1)  = mean_inflow_profiles(k,1)
1357                   v(k,j,nxlg:-1)  = mean_inflow_profiles(k,2)
1358                   w(k,j,nxlg:-1)  = 0.0_wp
1359                   pt(k,j,nxlg:-1) = mean_inflow_profiles(k,4)
1360                   e(k,j,nxlg:-1)  = mean_inflow_profiles(k,5)
1361                   IF ( humidity  .OR.  passive_scalar )                       &
1362                      q(k,j,nxlg:-1)  = mean_inflow_profiles(k,6)
1363                ENDDO
1364             ENDDO
1365          ENDIF
1366
1367!
1368!--       Calculate the damping factors to be used at the inflow. For a
1369!--       turbulent inflow the turbulent fluctuations have to be limited
1370!--       vertically because otherwise the turbulent inflow layer will grow
1371!--       in time.
1372          IF ( inflow_damping_height == 9999999.9_wp )  THEN
1373!
1374!--          Default: use the inversion height calculated by the prerun; if
1375!--          this is zero, inflow_damping_height must be explicitly
1376!--          specified.
1377             IF ( hom_sum(nzb+6,pr_palm,0) /= 0.0_wp )  THEN
1378                inflow_damping_height = hom_sum(nzb+6,pr_palm,0)
1379             ELSE
1380                WRITE( message_string, * ) 'inflow_damping_height must be ',   &
1381                     'explicitly specified because&the inversion height ',     &
1382                     'calculated by the prerun is zero.'
1383                CALL message( 'init_3d_model', 'PA0318', 1, 2, 0, 6, 0 )
1384             ENDIF
1385
1386          ENDIF
1387
1388          IF ( inflow_damping_width == 9999999.9_wp )  THEN
1389!
1390!--          Default for the transition range: one tenth of the undamped
1391!--          layer
1392             inflow_damping_width = 0.1_wp * inflow_damping_height
1393
1394          ENDIF
1395
1396          ALLOCATE( inflow_damping_factor(nzb:nzt+1) )
1397
1398          DO  k = nzb, nzt+1
1399
1400             IF ( zu(k) <= inflow_damping_height )  THEN
1401                inflow_damping_factor(k) = 1.0_wp
1402             ELSEIF ( zu(k) <= ( inflow_damping_height + inflow_damping_width ) )  THEN
1403                inflow_damping_factor(k) = 1.0_wp -                            &
1404                                           ( zu(k) - inflow_damping_height ) / &
1405                                           inflow_damping_width
1406             ELSE
1407                inflow_damping_factor(k) = 0.0_wp
1408             ENDIF
1409
1410          ENDDO
1411
1412       ENDIF
1413
1414!
1415!--    Inside buildings set velocities and TKE back to zero
1416       IF ( TRIM( initializing_actions ) == 'cyclic_fill' .AND.                &
1417            topography /= 'flat' )  THEN
1418!
1419!--       Inside buildings set velocities and TKE back to zero.
1420!--       Other scalars (pt, q, s, km, kh, p, sa, ...) are ignored at present,
1421!--       maybe revise later.
1422          DO  i = nxlg, nxrg
1423             DO  j = nysg, nyng
1424                u  (nzb:nzb_u_inner(j,i),j,i)   = 0.0_wp
1425                v  (nzb:nzb_v_inner(j,i),j,i)   = 0.0_wp
1426                w  (nzb:nzb_w_inner(j,i),j,i)   = 0.0_wp
1427                e  (nzb:nzb_w_inner(j,i),j,i)   = 0.0_wp
1428                tu_m(nzb:nzb_u_inner(j,i),j,i)  = 0.0_wp
1429                tv_m(nzb:nzb_v_inner(j,i),j,i)  = 0.0_wp
1430                tw_m(nzb:nzb_w_inner(j,i),j,i)  = 0.0_wp
1431                te_m(nzb:nzb_w_inner(j,i),j,i)  = 0.0_wp
1432                tpt_m(nzb:nzb_w_inner(j,i),j,i) = 0.0_wp
1433             ENDDO
1434          ENDDO
1435
1436       ENDIF
1437
1438!
1439!--    Calculate initial temperature field and other constants used in case
1440!--    of a sloping surface
1441       IF ( sloping_surface )  CALL init_slope
1442
1443!
1444!--    Initialize new time levels (only done in order to set boundary values
1445!--    including ghost points)
1446       e_p = e; pt_p = pt; u_p = u; v_p = v; w_p = w
1447       IF ( humidity  .OR.  passive_scalar )  THEN
1448          q_p = q
1449          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1450             qr_p = qr
1451             nr_p = nr
1452          ENDIF
1453       ENDIF
1454       IF ( ocean )  sa_p = sa
1455
1456!
1457!--    Allthough tendency arrays are set in prognostic_equations, they have
1458!--    have to be predefined here because they are used (but multiplied with 0)
1459!--    there before they are set.
1460       te_m = 0.0_wp; tpt_m = 0.0_wp; tu_m = 0.0_wp; tv_m = 0.0_wp; tw_m = 0.0_wp
1461       IF ( humidity  .OR.  passive_scalar )  THEN
1462          tq_m = 0.0_wp
1463          IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
1464             tqr_m = 0.0_wp
1465             tnr_m = 0.0_wp
1466          ENDIF
1467       ENDIF
1468       IF ( ocean )  tsa_m = 0.0_wp
1469
1470       CALL location_message( 'finished', .TRUE. )
1471
1472    ELSE
1473!
1474!--    Actually this part of the programm should not be reached
1475       message_string = 'unknown initializing problem'
1476       CALL message( 'init_3d_model', 'PA0193', 1, 2, 0, 6, 0 )
1477    ENDIF
1478
1479
1480    IF (  TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1481!
1482!--    Initialize old timelevels needed for radiation boundary conditions
1483       IF ( outflow_l )  THEN
1484          u_m_l(:,:,:) = u(:,:,1:2)
1485          v_m_l(:,:,:) = v(:,:,0:1)
1486          w_m_l(:,:,:) = w(:,:,0:1)
1487       ENDIF
1488       IF ( outflow_r )  THEN
1489          u_m_r(:,:,:) = u(:,:,nx-1:nx)
1490          v_m_r(:,:,:) = v(:,:,nx-1:nx)
1491          w_m_r(:,:,:) = w(:,:,nx-1:nx)
1492       ENDIF
1493       IF ( outflow_s )  THEN
1494          u_m_s(:,:,:) = u(:,0:1,:)
1495          v_m_s(:,:,:) = v(:,1:2,:)
1496          w_m_s(:,:,:) = w(:,0:1,:)
1497       ENDIF
1498       IF ( outflow_n )  THEN
1499          u_m_n(:,:,:) = u(:,ny-1:ny,:)
1500          v_m_n(:,:,:) = v(:,ny-1:ny,:)
1501          w_m_n(:,:,:) = w(:,ny-1:ny,:)
1502       ENDIF
1503       
1504    ENDIF
1505
1506!
1507!-- Calculate the initial volume flow at the right and north boundary
1508    IF ( conserve_volume_flow )  THEN
1509
1510       IF ( use_prescribed_profile_data )  THEN
1511
1512          volume_flow_initial_l = 0.0_wp
1513          volume_flow_area_l    = 0.0_wp
1514
1515          IF ( nxr == nx )  THEN
1516             DO  j = nys, nyn
1517                DO  k = nzb_u_inner(j,nx)+1, nzt
1518                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1519                                              u_init(k) * dzw(k)
1520                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1521                ENDDO
1522             ENDDO
1523          ENDIF
1524         
1525          IF ( nyn == ny )  THEN
1526             DO  i = nxl, nxr
1527                DO  k = nzb_v_inner(ny,i)+1, nzt
1528                   volume_flow_initial_l(2) = volume_flow_initial_l(2) + &
1529                                              v_init(k) * dzw(k)
1530                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1531                ENDDO
1532             ENDDO
1533          ENDIF
1534
1535#if defined( __parallel )
1536          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1537                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1538          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1539                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1540
1541#else
1542          volume_flow_initial = volume_flow_initial_l
1543          volume_flow_area    = volume_flow_area_l
1544#endif 
1545
1546       ELSEIF ( TRIM( initializing_actions ) == 'cyclic_fill' )  THEN
1547
1548          volume_flow_initial_l = 0.0_wp
1549          volume_flow_area_l    = 0.0_wp
1550
1551          IF ( nxr == nx )  THEN
1552             DO  j = nys, nyn
1553                DO  k = nzb_u_inner(j,nx)+1, nzt
1554                   volume_flow_initial_l(1) = volume_flow_initial_l(1) +       &
1555                                              hom_sum(k,1,0) * dzw(k)
1556                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1557                ENDDO
1558             ENDDO
1559          ENDIF
1560         
1561          IF ( nyn == ny )  THEN
1562             DO  i = nxl, nxr
1563                DO  k = nzb_v_inner(ny,i)+1, nzt
1564                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1565                                              hom_sum(k,2,0) * dzw(k)
1566                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1567                ENDDO
1568             ENDDO
1569          ENDIF
1570
1571#if defined( __parallel )
1572          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1573                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1574          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1575                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1576
1577#else
1578          volume_flow_initial = volume_flow_initial_l
1579          volume_flow_area    = volume_flow_area_l
1580#endif 
1581
1582       ELSEIF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
1583
1584          volume_flow_initial_l = 0.0_wp
1585          volume_flow_area_l    = 0.0_wp
1586
1587          IF ( nxr == nx )  THEN
1588             DO  j = nys, nyn
1589                DO  k = nzb_u_inner(j,nx)+1, nzt
1590                   volume_flow_initial_l(1) = volume_flow_initial_l(1) + &
1591                                              u(k,j,nx) * dzw(k)
1592                   volume_flow_area_l(1)    = volume_flow_area_l(1) + dzw(k)
1593                ENDDO
1594             ENDDO
1595          ENDIF
1596         
1597          IF ( nyn == ny )  THEN
1598             DO  i = nxl, nxr
1599                DO  k = nzb_v_inner(ny,i)+1, nzt
1600                   volume_flow_initial_l(2) = volume_flow_initial_l(2) +       &
1601                                              v(k,ny,i) * dzw(k)
1602                   volume_flow_area_l(2)    = volume_flow_area_l(2) + dzw(k)
1603                ENDDO
1604             ENDDO
1605          ENDIF
1606
1607#if defined( __parallel )
1608          CALL MPI_ALLREDUCE( volume_flow_initial_l(1), volume_flow_initial(1),&
1609                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1610          CALL MPI_ALLREDUCE( volume_flow_area_l(1), volume_flow_area(1),      &
1611                              2, MPI_REAL, MPI_SUM, comm2d, ierr )
1612
1613#else
1614          volume_flow_initial = volume_flow_initial_l
1615          volume_flow_area    = volume_flow_area_l
1616#endif 
1617
1618       ENDIF
1619
1620!
1621!--    In case of 'bulk_velocity' mode, volume_flow_initial is calculated
1622!--    from u|v_bulk instead
1623       IF ( TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' )  THEN
1624          volume_flow_initial(1) = u_bulk * volume_flow_area(1)
1625          volume_flow_initial(2) = v_bulk * volume_flow_area(2)
1626       ENDIF
1627
1628    ENDIF
1629
1630!
1631!-- Initialize quantities for special advections schemes
1632    CALL init_advec
1633
1634!
1635!-- Impose random perturbation on the horizontal velocity field and then
1636!-- remove the divergences from the velocity field at the initial stage
1637    IF ( create_disturbances  .AND.  disturbance_energy_limit /= 0.0_wp  .AND. &
1638         TRIM( initializing_actions ) /= 'read_restart_data'  .AND.            &
1639         TRIM( initializing_actions ) /= 'cyclic_fill' )  THEN
1640
1641       CALL location_message( 'creating initial disturbances', .FALSE. )
1642       CALL disturb_field( nzb_u_inner, tend, u )
1643       CALL disturb_field( nzb_v_inner, tend, v )
1644       CALL location_message( 'finished', .TRUE. )
1645
1646       CALL location_message( 'calling pressure solver', .FALSE. )
1647       n_sor = nsor_ini
1648       !$acc data copyin( d, ddzu, ddzw, nzb_s_inner, nzb_u_inner )            &
1649       !$acc      copyin( nzb_v_inner, nzb_w_inner, p, rflags_s_inner, tend )  &
1650       !$acc      copyin( weight_pres, weight_substep )                        &
1651       !$acc      copy( tri, tric, u, v, w )
1652       CALL pres
1653       !$acc end data
1654       n_sor = nsor
1655       CALL location_message( 'finished', .TRUE. )
1656
1657    ENDIF
1658
1659!
1660!-- If required, initialize quantities needed for the plant canopy model
1661    CALL location_message( 'initializing plant canopy model', .FALSE. )
1662    IF ( plant_canopy )  CALL pcm_init
1663    CALL location_message( 'finished', .TRUE. )
1664
1665!
1666!-- If required, initialize dvrp-software
1667    IF ( dt_dvrp /= 9999999.9_wp )  CALL init_dvrp
1668
1669    IF ( ocean )  THEN
1670!
1671!--    Initialize quantities needed for the ocean model
1672       CALL init_ocean
1673
1674    ELSE
1675!
1676!--    Initialize quantities for handling cloud physics
1677!--    This routine must be called before lpm_init, because
1678!--    otherwise, array pt_d_t, needed in data_output_dvrp (called by
1679!--    lpm_init) is not defined.
1680       CALL init_cloud_physics
1681!
1682!--    Initialize bulk cloud microphysics
1683       CALL microphysics_init
1684    ENDIF
1685
1686!
1687!-- If required, initialize particles
1688    IF ( particle_advection )  CALL lpm_init
1689
1690!
1691!-- If required, initialize quantities needed for the LSM
1692    IF ( land_surface )  THEN
1693       CALL location_message( 'initializing land surface model', .FALSE. )
1694       CALL lsm_init
1695       CALL location_message( 'finished', .TRUE. )
1696    ENDIF
1697
1698!
1699!-- Initialize surface layer (done after LSM as roughness length are required
1700!-- for initialization
1701    IF ( constant_flux_layer )  THEN
1702       CALL location_message( 'initializing surface layer', .FALSE. )
1703       CALL init_surface_layer_fluxes
1704       CALL location_message( 'finished', .TRUE. )
1705    ENDIF
1706
1707!
1708!-- If required, initialize radiation model
1709    IF ( radiation )  THEN
1710       CALL location_message( 'initializing radiation model', .FALSE. )
1711       CALL radiation_init
1712       CALL location_message( 'finished', .TRUE. )
1713    ENDIF
1714   
1715!
1716!-- If required, initialize quantities needed for the wind turbine model
1717    IF ( wind_turbine )  THEN
1718       CALL location_message( 'initializing wind turbine model', .FALSE. )
1719       CALL wtm_init
1720       CALL location_message( 'finished', .TRUE. )
1721    ENDIF
1722
1723
1724!
1725!-- Initialize the ws-scheme.   
1726    IF ( ws_scheme_sca .OR. ws_scheme_mom )  CALL ws_init       
1727
1728!
1729!-- Setting weighting factors for calculation of perturbation pressure
1730!-- and turbulent quantities from the RK substeps
1731    IF ( TRIM(timestep_scheme) == 'runge-kutta-3' )  THEN      ! for RK3-method
1732
1733       weight_substep(1) = 1._wp/6._wp
1734       weight_substep(2) = 3._wp/10._wp
1735       weight_substep(3) = 8._wp/15._wp
1736
1737       weight_pres(1)    = 1._wp/3._wp
1738       weight_pres(2)    = 5._wp/12._wp
1739       weight_pres(3)    = 1._wp/4._wp
1740
1741    ELSEIF ( TRIM(timestep_scheme) == 'runge-kutta-2' )  THEN  ! for RK2-method
1742
1743       weight_substep(1) = 1._wp/2._wp
1744       weight_substep(2) = 1._wp/2._wp
1745         
1746       weight_pres(1)    = 1._wp/2._wp
1747       weight_pres(2)    = 1._wp/2._wp       
1748
1749    ELSE                                     ! for Euler-method
1750
1751       weight_substep(1) = 1.0_wp     
1752       weight_pres(1)    = 1.0_wp                   
1753
1754    ENDIF
1755
1756!
1757!-- Initialize Rayleigh damping factors
1758    rdf    = 0.0_wp
1759    rdf_sc = 0.0_wp
1760    IF ( rayleigh_damping_factor /= 0.0_wp )  THEN
1761       IF (  .NOT.  ocean )  THEN
1762          DO  k = nzb+1, nzt
1763             IF ( zu(k) >= rayleigh_damping_height )  THEN
1764                rdf(k) = rayleigh_damping_factor *                             &
1765                      ( SIN( pi * 0.5_wp * ( zu(k) - rayleigh_damping_height ) &
1766                             / ( zu(nzt) - rayleigh_damping_height ) )         &
1767                      )**2
1768             ENDIF
1769          ENDDO
1770       ELSE
1771          DO  k = nzt, nzb+1, -1
1772             IF ( zu(k) <= rayleigh_damping_height )  THEN
1773                rdf(k) = rayleigh_damping_factor *                             &
1774                      ( SIN( pi * 0.5_wp * ( rayleigh_damping_height - zu(k) ) &
1775                             / ( rayleigh_damping_height - zu(nzb+1) ) )       &
1776                      )**2
1777             ENDIF
1778          ENDDO
1779       ENDIF
1780    ENDIF
1781    IF ( scalar_rayleigh_damping )  rdf_sc = rdf
1782
1783!
1784!-- Initialize the starting level and the vertical smoothing factor used for
1785!-- the external pressure gradient
1786    dp_smooth_factor = 1.0_wp
1787    IF ( dp_external )  THEN
1788!
1789!--    Set the starting level dp_level_ind_b only if it has not been set before
1790!--    (e.g. in init_grid).
1791       IF ( dp_level_ind_b == 0 )  THEN
1792          ind_array = MINLOC( ABS( dp_level_b - zu ) )
1793          dp_level_ind_b = ind_array(1) - 1 + nzb 
1794                                        ! MINLOC uses lower array bound 1
1795       ENDIF
1796       IF ( dp_smooth )  THEN
1797          dp_smooth_factor(:dp_level_ind_b) = 0.0_wp
1798          DO  k = dp_level_ind_b+1, nzt
1799             dp_smooth_factor(k) = 0.5_wp * ( 1.0_wp + SIN( pi *               &
1800                        ( REAL( k - dp_level_ind_b, KIND=wp ) /                &
1801                          REAL( nzt - dp_level_ind_b, KIND=wp ) - 0.5_wp ) ) )
1802          ENDDO
1803       ENDIF
1804    ENDIF
1805
1806!
1807!-- Initialize damping zone for the potential temperature in case of
1808!-- non-cyclic lateral boundaries. The damping zone has the maximum value
1809!-- at the inflow boundary and decreases to zero at pt_damping_width.
1810    ptdf_x = 0.0_wp
1811    ptdf_y = 0.0_wp
1812    IF ( bc_lr_dirrad )  THEN
1813       DO  i = nxl, nxr
1814          IF ( ( i * dx ) < pt_damping_width )  THEN
1815             ptdf_x(i) = pt_damping_factor * ( SIN( pi * 0.5_wp *              &
1816                            REAL( pt_damping_width - i * dx, KIND=wp ) / (     &
1817                            REAL( pt_damping_width, KIND=wp ) ) ) )**2 
1818          ENDIF
1819       ENDDO
1820    ELSEIF ( bc_lr_raddir )  THEN
1821       DO  i = nxl, nxr
1822          IF ( ( i * dx ) > ( nx * dx - pt_damping_width ) )  THEN
1823             ptdf_x(i) = pt_damping_factor *                                   &
1824                         SIN( pi * 0.5_wp *                                    &
1825                                 ( ( i - nx ) * dx + pt_damping_width ) /      &
1826                                 REAL( pt_damping_width, KIND=wp ) )**2
1827          ENDIF
1828       ENDDO 
1829    ELSEIF ( bc_ns_dirrad )  THEN
1830       DO  j = nys, nyn
1831          IF ( ( j * dy ) > ( ny * dy - pt_damping_width ) )  THEN
1832             ptdf_y(j) = pt_damping_factor *                                   &
1833                         SIN( pi * 0.5_wp *                                    &
1834                                 ( ( j - ny ) * dy + pt_damping_width ) /      &
1835                                 REAL( pt_damping_width, KIND=wp ) )**2
1836          ENDIF
1837       ENDDO 
1838    ELSEIF ( bc_ns_raddir )  THEN
1839       DO  j = nys, nyn
1840          IF ( ( j * dy ) < pt_damping_width )  THEN
1841             ptdf_y(j) = pt_damping_factor *                                   &
1842                         SIN( pi * 0.5_wp *                                    &
1843                                ( pt_damping_width - j * dy ) /                &
1844                                REAL( pt_damping_width, KIND=wp ) )**2
1845          ENDIF
1846       ENDDO
1847    ENDIF
1848
1849!
1850!-- Pre-set masks for regional statistics. Default is the total model domain.
1851!-- Ghost points are excluded because counting values at the ghost boundaries
1852!-- would bias the statistics
1853    rmask = 1.0_wp
1854    rmask(:,nxlg:nxl-1,:) = 0.0_wp;  rmask(:,nxr+1:nxrg,:) = 0.0_wp
1855    rmask(nysg:nys-1,:,:) = 0.0_wp;  rmask(nyn+1:nyng,:,:) = 0.0_wp
1856
1857!
1858!-- User-defined initializing actions. Check afterwards, if maximum number
1859!-- of allowed timeseries is exceeded
1860    CALL user_init
1861
1862    IF ( dots_num > dots_max )  THEN
1863       WRITE( message_string, * ) 'number of time series quantities exceeds',  &
1864                                  ' its maximum of dots_max = ', dots_max,     &
1865                                  ' &Please increase dots_max in modules.f90.'
1866       CALL message( 'init_3d_model', 'PA0194', 1, 2, 0, 6, 0 )   
1867    ENDIF
1868
1869!
1870!-- Input binary data file is not needed anymore. This line must be placed
1871!-- after call of user_init!
1872    CALL close_file( 13 )
1873
1874!
1875!-- Compute total sum of active mask grid points
1876!-- and the mean surface level height for each statistic region
1877!-- ngp_2dh: number of grid points of a horizontal cross section through the
1878!--          total domain
1879!-- ngp_3d:  number of grid points of the total domain
1880    ngp_2dh_outer_l   = 0
1881    ngp_2dh_outer     = 0
1882    ngp_2dh_s_inner_l = 0
1883    ngp_2dh_s_inner   = 0
1884    ngp_2dh_l         = 0
1885    ngp_2dh           = 0
1886    ngp_3d_inner_l    = 0.0_wp
1887    ngp_3d_inner      = 0
1888    ngp_3d            = 0
1889    ngp_sums          = ( nz + 2 ) * ( pr_palm + max_pr_user )
1890
1891    mean_surface_level_height   = 0.0_wp
1892    mean_surface_level_height_l = 0.0_wp
1893
1894    DO  sr = 0, statistic_regions
1895       DO  i = nxl, nxr
1896          DO  j = nys, nyn
1897             IF ( rmask(j,i,sr) == 1.0_wp )  THEN
1898!
1899!--             All xy-grid points
1900                ngp_2dh_l(sr) = ngp_2dh_l(sr) + 1
1901                mean_surface_level_height_l(sr) = mean_surface_level_height_l(sr) &
1902                                                  + zw(nzb_s_inner(j,i))
1903!
1904!--             xy-grid points above topography
1905                DO  k = nzb_s_outer(j,i), nz + 1
1906                   ngp_2dh_outer_l(k,sr) = ngp_2dh_outer_l(k,sr) + 1
1907                ENDDO
1908                DO  k = nzb_s_inner(j,i), nz + 1
1909                   ngp_2dh_s_inner_l(k,sr) = ngp_2dh_s_inner_l(k,sr) + 1
1910                ENDDO
1911!
1912!--             All grid points of the total domain above topography
1913                ngp_3d_inner_l(sr) = ngp_3d_inner_l(sr)                        &
1914                                     + ( nz - nzb_s_inner(j,i) + 2 )
1915             ENDIF
1916          ENDDO
1917       ENDDO
1918    ENDDO
1919
1920    sr = statistic_regions + 1
1921#if defined( __parallel )
1922    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1923    CALL MPI_ALLREDUCE( ngp_2dh_l(0), ngp_2dh(0), sr, MPI_INTEGER, MPI_SUM,    &
1924                        comm2d, ierr )
1925    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1926    CALL MPI_ALLREDUCE( ngp_2dh_outer_l(0,0), ngp_2dh_outer(0,0), (nz+2)*sr,   &
1927                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
1928    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1929    CALL MPI_ALLREDUCE( ngp_2dh_s_inner_l(0,0), ngp_2dh_s_inner(0,0),          &
1930                        (nz+2)*sr, MPI_INTEGER, MPI_SUM, comm2d, ierr )
1931    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1932    CALL MPI_ALLREDUCE( ngp_3d_inner_l(0), ngp_3d_inner_tmp(0), sr, MPI_REAL,  &
1933                        MPI_SUM, comm2d, ierr )
1934    ngp_3d_inner = INT( ngp_3d_inner_tmp, KIND = SELECTED_INT_KIND( 18 ) )
1935    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1936    CALL MPI_ALLREDUCE( mean_surface_level_height_l(0),                        &
1937                        mean_surface_level_height(0), sr, MPI_REAL,            &
1938                        MPI_SUM, comm2d, ierr )
1939    mean_surface_level_height = mean_surface_level_height / REAL( ngp_2dh )
1940#else
1941    ngp_2dh         = ngp_2dh_l
1942    ngp_2dh_outer   = ngp_2dh_outer_l
1943    ngp_2dh_s_inner = ngp_2dh_s_inner_l
1944    ngp_3d_inner    = INT( ngp_3d_inner_l, KIND = SELECTED_INT_KIND( 18 ) )
1945    mean_surface_level_height = mean_surface_level_height_l / REAL( ngp_2dh_l )
1946#endif
1947
1948    ngp_3d = INT ( ngp_2dh, KIND = SELECTED_INT_KIND( 18 ) ) * &
1949             INT ( (nz + 2 ), KIND = SELECTED_INT_KIND( 18 ) )
1950
1951!
1952!-- Set a lower limit of 1 in order to avoid zero divisions in flow_statistics,
1953!-- buoyancy, etc. A zero value will occur for cases where all grid points of
1954!-- the respective subdomain lie below the surface topography
1955    ngp_2dh_outer   = MAX( 1, ngp_2dh_outer(:,:)   ) 
1956    ngp_3d_inner    = MAX( INT(1, KIND = SELECTED_INT_KIND( 18 )),             &
1957                           ngp_3d_inner(:) )
1958    ngp_2dh_s_inner = MAX( 1, ngp_2dh_s_inner(:,:) ) 
1959
1960    DEALLOCATE( mean_surface_level_height_l, ngp_2dh_l, ngp_2dh_outer_l,       &
1961                ngp_3d_inner_l, ngp_3d_inner_tmp )
1962
1963    CALL location_message( 'leaving init_3d_model', .TRUE. )
1964
1965 END SUBROUTINE init_3d_model
Note: See TracBrowser for help on using the repository browser.