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

Last change on this file since 1903 was 1879, checked in by hellstea, 9 years ago

last commit documented

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