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