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

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

Separate balance equations for humidity and passive_scalar

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