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

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

NetCDF routines modularized; new parameter netcdf_deflate; further changes in the pmc

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