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

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

bugfixes for calculating run control quantities, bugfix for calculating pressure with fft-method in case of Neumann conditions both at bottom and top, steering of pres modified, ocean mode now uses initial density profile as reference in the buoyancy term

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