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

Last change on this file since 1361 was 1361, checked in by hoffmann, 10 years ago

improved version of two-moment cloud physics

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