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

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