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

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

last commit documented

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