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

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