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

Last change on this file since 2001 was 2001, checked in by knoop, 8 years ago

last commit documented

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