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

Last change on this file since 1822 was 1822, checked in by hoffmann, 8 years ago

changes in LPM and bulk cloud microphysics

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