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

Last change on this file since 1764 was 1764, checked in by raasch, 9 years ago

update of the nested domain system + some bugfixes

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