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