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

Last change on this file since 1833 was 1833, checked in by raasch, 8 years ago

spectrum renamed spactra_par and further modularized, POINTER-attributes added in coupler-routines to avoid gfortran error messages

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