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