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

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

flight module added

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