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

Last change on this file since 1831 was 1831, checked in by hoffmann, 8 years ago

cloud physics variables renamed

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