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

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

Avoid segmentation fault (see change in 1915) by different initialization of us instead of adding a very small number in the denominator

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