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

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