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

Last change on this file since 1818 was 1818, checked in by maronga, 8 years ago

last commit documented / copyright update

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