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

Last change on this file since 1858 was 1852, checked in by hoffmann, 9 years ago

last commit documented

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