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