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