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

Last change on this file since 1682 was 1682, checked in by knoop, 9 years ago

Code annotations made doxygen readable

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