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

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

last commit documented

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