[1682] | 1 | !> @file time_integration.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 | ! |
---|
[484] | 20 | ! Current revisions: |
---|
[1092] | 21 | ! ------------------ |
---|
[1928] | 22 | ! |
---|
[2119] | 23 | ! |
---|
[1366] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: time_integration.f90 2119 2017-01-17 16:51:50Z scharf $ |
---|
| 27 | ! |
---|
[2119] | 28 | ! 2118 2017-01-17 16:38:49Z raasch |
---|
| 29 | ! OpenACC directives and related code removed |
---|
| 30 | ! |
---|
[2051] | 31 | ! 2050 2016-11-08 15:00:55Z gronemeier |
---|
| 32 | ! Implement turbulent outflow condition |
---|
| 33 | ! |
---|
[2032] | 34 | ! 2031 2016-10-21 15:11:58Z knoop |
---|
| 35 | ! renamed variable rho to rho_ocean |
---|
| 36 | ! |
---|
[2012] | 37 | ! 2011 2016-09-19 17:29:57Z kanani |
---|
| 38 | ! Flag urban_surface is now defined in module control_parameters, |
---|
| 39 | ! removed commented CALLs of global_min_max. |
---|
| 40 | ! |
---|
[2008] | 41 | ! 2007 2016-08-24 15:47:17Z kanani |
---|
| 42 | ! Added CALLs for new urban surface model |
---|
| 43 | ! |
---|
[2001] | 44 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
| 45 | ! Forced header and separation lines into 80 columns |
---|
| 46 | ! |
---|
[1977] | 47 | ! 1976 2016-07-27 13:28:04Z maronga |
---|
| 48 | ! Simplified calls to radiation model |
---|
| 49 | ! |
---|
[1961] | 50 | ! 1960 2016-07-12 16:34:24Z suehring |
---|
| 51 | ! Separate humidity and passive scalar |
---|
| 52 | ! |
---|
[1958] | 53 | ! 1957 2016-07-07 10:43:48Z suehring |
---|
| 54 | ! flight module added |
---|
| 55 | ! |
---|
[1933] | 56 | ! 1919 2016-05-27 14:51:23Z raasch |
---|
| 57 | ! Initial version of purely vertical nesting introduced. |
---|
[1928] | 58 | ! |
---|
[1919] | 59 | ! 1918 2016-05-27 14:35:57Z raasch |
---|
| 60 | ! determination of time step moved to the end of the time step loop, |
---|
| 61 | ! the first time step is now always calculated before the time step loop (i.e. |
---|
| 62 | ! also in case of restart runs) |
---|
| 63 | ! |
---|
[1917] | 64 | ! 1914 2016-05-26 14:44:07Z witha |
---|
| 65 | ! Added call for wind turbine model |
---|
| 66 | ! |
---|
[1879] | 67 | ! 1878 2016-04-19 12:30:36Z hellstea |
---|
| 68 | ! Synchronization for nested runs rewritten |
---|
| 69 | ! |
---|
[1854] | 70 | ! 1853 2016-04-11 09:00:35Z maronga |
---|
| 71 | ! Adjusted for use with radiation_scheme = constant |
---|
| 72 | ! |
---|
[1851] | 73 | ! 1849 2016-04-08 11:33:18Z hoffmann |
---|
| 74 | ! Adapted for modularization of microphysics |
---|
| 75 | ! |
---|
[1834] | 76 | ! 1833 2016-04-07 14:23:03Z raasch |
---|
| 77 | ! spectrum renamed spectra_mod, spectra related variables moved to spectra_mod |
---|
| 78 | ! |
---|
[1832] | 79 | ! 1831 2016-04-07 13:15:51Z hoffmann |
---|
| 80 | ! turbulence renamed collision_turbulence |
---|
| 81 | ! |
---|
[1823] | 82 | ! 1822 2016-04-07 07:49:42Z hoffmann |
---|
| 83 | ! icloud_scheme replaced by microphysics_* |
---|
| 84 | ! |
---|
[1809] | 85 | ! 1808 2016-04-05 19:44:00Z raasch |
---|
| 86 | ! output message in case unscheduled radiation calls removed |
---|
| 87 | ! |
---|
[1798] | 88 | ! 1797 2016-03-21 16:50:28Z raasch |
---|
| 89 | ! introduction of different datatransfer modes |
---|
| 90 | ! |
---|
[1792] | 91 | ! 1791 2016-03-11 10:41:25Z raasch |
---|
| 92 | ! call of pmci_update_new removed |
---|
| 93 | ! |
---|
[1787] | 94 | ! 1786 2016-03-08 05:49:27Z raasch |
---|
| 95 | ! +module spectrum |
---|
| 96 | ! |
---|
[1784] | 97 | ! 1783 2016-03-06 18:36:17Z raasch |
---|
| 98 | ! switch back of netcdf data format for mask output moved to the mask output |
---|
| 99 | ! routine |
---|
| 100 | ! |
---|
[1782] | 101 | ! 1781 2016-03-03 15:12:23Z raasch |
---|
| 102 | ! some pmc calls removed at the beginning (before timeloop), |
---|
| 103 | ! pmc initialization moved to the main program |
---|
| 104 | ! |
---|
[1765] | 105 | ! 1764 2016-02-28 12:45:19Z raasch |
---|
| 106 | ! PMC_ACTIVE flags removed, |
---|
| 107 | ! bugfix: nest synchronization after first call of timestep |
---|
| 108 | ! |
---|
[1763] | 109 | ! 1762 2016-02-25 12:31:13Z hellstea |
---|
| 110 | ! Introduction of nested domain feature |
---|
| 111 | ! |
---|
[1737] | 112 | ! 1736 2015-12-04 08:56:33Z raasch |
---|
| 113 | ! no perturbations added to total domain if energy limit has been set zero |
---|
| 114 | ! |
---|
[1692] | 115 | ! 1691 2015-10-26 16:17:44Z maronga |
---|
| 116 | ! Added option for spin-ups without land surface and radiation models. Moved calls |
---|
| 117 | ! for radiation and lan surface schemes. |
---|
| 118 | ! |
---|
[1683] | 119 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
| 120 | ! Code annotations made doxygen readable |
---|
| 121 | ! |
---|
[1672] | 122 | ! 1671 2015-09-25 03:29:37Z raasch |
---|
| 123 | ! bugfix: ghostpoint exchange for array diss in case that sgs velocities are used |
---|
| 124 | ! for particles |
---|
| 125 | ! |
---|
[1586] | 126 | ! 1585 2015-04-30 07:05:52Z maronga |
---|
| 127 | ! Moved call of radiation scheme. Added support for RRTM |
---|
| 128 | ! |
---|
[1552] | 129 | ! 1551 2015-03-03 14:18:16Z maronga |
---|
| 130 | ! Added interface for different radiation schemes. |
---|
| 131 | ! |
---|
[1497] | 132 | ! 1496 2014-12-02 17:25:50Z maronga |
---|
| 133 | ! Added calls for the land surface model and radiation scheme |
---|
| 134 | ! |
---|
[1403] | 135 | ! 1402 2014-05-09 14:25:13Z raasch |
---|
| 136 | ! location messages modified |
---|
| 137 | ! |
---|
[1385] | 138 | ! 1384 2014-05-02 14:31:06Z raasch |
---|
| 139 | ! location messages added |
---|
| 140 | ! |
---|
[1381] | 141 | ! 1380 2014-04-28 12:40:45Z heinze |
---|
| 142 | ! CALL of nudge_ref added |
---|
| 143 | ! bc_pt_t_val and bc_q_t_val are updated in case nudging is used |
---|
| 144 | ! |
---|
[1366] | 145 | ! 1365 2014-04-22 15:03:56Z boeske |
---|
[1365] | 146 | ! Reset sums_ls_l to zero at each timestep |
---|
| 147 | ! +sums_ls_l |
---|
| 148 | ! Calculation of reference state (previously in subroutine calc_mean_profile) |
---|
| 149 | |
---|
[1343] | 150 | ! 1342 2014-03-26 17:04:47Z kanani |
---|
| 151 | ! REAL constants defined as wp-kind |
---|
| 152 | ! |
---|
[1321] | 153 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 154 | ! ONLY-attribute added to USE-statements, |
---|
| 155 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 156 | ! kinds are defined in new module kinds, |
---|
| 157 | ! old module precision_kind is removed, |
---|
| 158 | ! revision history before 2012 removed, |
---|
| 159 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 160 | ! all variable declaration statements |
---|
[1319] | 161 | ! 1318 2014-03-17 13:35:16Z raasch |
---|
| 162 | ! module interfaces removed |
---|
| 163 | ! |
---|
[1309] | 164 | ! 1308 2014-03-13 14:58:42Z fricke |
---|
| 165 | ! +netcdf_data_format_save |
---|
| 166 | ! For masked data, parallel netcdf output is not tested so far, hence |
---|
| 167 | ! netcdf_data_format is switched back to non-paralell output. |
---|
| 168 | ! |
---|
[1277] | 169 | ! 1276 2014-01-15 13:40:41Z heinze |
---|
| 170 | ! Use LSF_DATA also in case of Dirichlet bottom boundary condition for scalars |
---|
| 171 | ! |
---|
[1258] | 172 | ! 1257 2013-11-08 15:18:40Z raasch |
---|
| 173 | ! acc-update-host directive for timestep removed |
---|
| 174 | ! |
---|
[1242] | 175 | ! 1241 2013-10-30 11:36:58Z heinze |
---|
| 176 | ! Generalize calc_mean_profile for wider use |
---|
| 177 | ! Determine shf and qsws in dependence on data from LSF_DATA |
---|
| 178 | ! Determine ug and vg in dependence on data from LSF_DATA |
---|
[1222] | 179 | ! 1221 2013-09-10 08:59:13Z raasch |
---|
| 180 | ! host update of arrays before timestep is called |
---|
| 181 | ! |
---|
[1182] | 182 | ! 1179 2013-06-14 05:57:58Z raasch |
---|
| 183 | ! mean profiles for reference state are only calculated if required, |
---|
| 184 | ! small bugfix for background communication |
---|
| 185 | ! |
---|
[1172] | 186 | ! 1171 2013-05-30 11:27:45Z raasch |
---|
| 187 | ! split of prognostic_equations deactivated (comment lines), for the time being |
---|
| 188 | ! |
---|
[1132] | 189 | ! 1128 2013-04-12 06:19:32Z raasch |
---|
[1128] | 190 | ! asynchronous transfer of ghost point data realized for acc-optimized version: |
---|
| 191 | ! prognostic_equations are first called two times for those points required for |
---|
| 192 | ! the left-right and north-south exchange, respectively, and then for the |
---|
| 193 | ! remaining points, |
---|
| 194 | ! those parts requiring global communication moved from prognostic_equations to |
---|
| 195 | ! here |
---|
[392] | 196 | ! |
---|
[1116] | 197 | ! 1115 2013-03-26 18:16:16Z hoffmann |
---|
| 198 | ! calculation of qr and nr is restricted to precipitation |
---|
| 199 | ! |
---|
[1114] | 200 | ! 1113 2013-03-10 02:48:14Z raasch |
---|
| 201 | ! GPU-porting of boundary conditions, |
---|
| 202 | ! openACC directives updated |
---|
| 203 | ! formal parameter removed from routine boundary_conds |
---|
| 204 | ! |
---|
[1112] | 205 | ! 1111 2013-03-08 23:54:10Z raasch |
---|
| 206 | ! +internal timestep counter for cpu statistics added, |
---|
| 207 | ! openACC directives updated |
---|
| 208 | ! |
---|
[1093] | 209 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
| 210 | ! unused variables removed |
---|
| 211 | ! |
---|
[1066] | 212 | ! 1065 2012-11-22 17:42:36Z hoffmann |
---|
| 213 | ! exchange of diss (dissipation rate) in case of turbulence = .TRUE. added |
---|
| 214 | ! |
---|
[1054] | 215 | ! 1053 2012-11-13 17:11:03Z hoffmann |
---|
| 216 | ! exchange of ghost points for nr, qr added |
---|
| 217 | ! |
---|
[1037] | 218 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 219 | ! code put under GPL (PALM 3.9) |
---|
| 220 | ! |
---|
[1020] | 221 | ! 1019 2012-09-28 06:46:45Z raasch |
---|
| 222 | ! non-optimized version of prognostic_equations removed |
---|
| 223 | ! |
---|
[1017] | 224 | ! 1015 2012-09-27 09:23:24Z raasch |
---|
| 225 | ! +call of prognostic_equations_acc |
---|
| 226 | ! |
---|
[1002] | 227 | ! 1001 2012-09-13 14:08:46Z raasch |
---|
| 228 | ! all actions concerning leapfrog- and upstream-spline-scheme removed |
---|
| 229 | ! |
---|
[850] | 230 | ! 849 2012-03-15 10:35:09Z raasch |
---|
| 231 | ! advec_particles renamed lpm, first_call_advec_particles renamed first_call_lpm |
---|
| 232 | ! |
---|
[826] | 233 | ! 825 2012-02-19 03:03:44Z raasch |
---|
| 234 | ! wang_collision_kernel renamed wang_kernel |
---|
| 235 | ! |
---|
[1] | 236 | ! Revision 1.1 1997/08/11 06:19:04 raasch |
---|
| 237 | ! Initial revision |
---|
| 238 | ! |
---|
| 239 | ! |
---|
| 240 | ! Description: |
---|
| 241 | ! ------------ |
---|
[1682] | 242 | !> Integration in time of the model equations, statistical analysis and graphic |
---|
| 243 | !> output |
---|
[1] | 244 | !------------------------------------------------------------------------------! |
---|
[1682] | 245 | SUBROUTINE time_integration |
---|
| 246 | |
---|
[1] | 247 | |
---|
[1320] | 248 | USE advec_ws, & |
---|
| 249 | ONLY: ws_statistics |
---|
| 250 | |
---|
| 251 | USE arrays_3d, & |
---|
[1762] | 252 | ONLY: diss, dzu, e, e_p, nr_p, prho, pt, pt_p, pt_init, q_init, q, & |
---|
[2031] | 253 | ql, ql_c, ql_v, ql_vp, qr_p, q_p, ref_state, rho_ocean, s, s_p, sa_p, & |
---|
[1960] | 254 | tend, u, u_p, v, vpt, v_p, w, w_p |
---|
[1320] | 255 | |
---|
[1365] | 256 | USE calc_mean_profile_mod, & |
---|
[1320] | 257 | ONLY: calc_mean_profile |
---|
| 258 | |
---|
| 259 | USE control_parameters, & |
---|
| 260 | ONLY: advected_distance_x, advected_distance_y, average_count_3d, & |
---|
[1833] | 261 | averaging_interval, averaging_interval_pr, & |
---|
| 262 | bc_lr_cyc, bc_ns_cyc, bc_pt_t_val, & |
---|
[1380] | 263 | bc_q_t_val, call_psolver_at_all_substeps, cloud_droplets, & |
---|
[1691] | 264 | cloud_physics, constant_flux_layer, constant_heatflux, & |
---|
| 265 | create_disturbances, dopr_n, constant_diffusion, coupling_mode, & |
---|
| 266 | coupling_start_time, current_timestep_number, & |
---|
| 267 | disturbance_created, disturbance_energy_limit, dist_range, & |
---|
| 268 | do_sum, dt_3d, dt_averaging_input, dt_averaging_input_pr, & |
---|
| 269 | dt_coupling, dt_data_output_av, dt_disturb, dt_do2d_xy, & |
---|
| 270 | dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_domask,dt_dopts, dt_dopr, & |
---|
[1833] | 271 | dt_dopr_listing, dt_dots, dt_dvrp, dt_run_control, & |
---|
[1320] | 272 | end_time, first_call_lpm, galilei_transformation, humidity, & |
---|
[1822] | 273 | intermediate_timestep_count, & |
---|
[1320] | 274 | intermediate_timestep_count_max, large_scale_forcing, & |
---|
[1822] | 275 | loop_optimization, lsf_surf, lsf_vert, masks, & |
---|
| 276 | microphysics_seifert, mid, nest_domain, & |
---|
[1783] | 277 | neutral, nr_timesteps_this_run, nudging, & |
---|
[2118] | 278 | ocean, passive_scalar, & |
---|
[1320] | 279 | prho_reference, pt_reference, pt_slope_offset, random_heatflux, & |
---|
| 280 | run_coupled, simulated_time, simulated_time_chr, & |
---|
| 281 | skip_time_do2d_xy, skip_time_do2d_xz, skip_time_do2d_yz, & |
---|
| 282 | skip_time_do3d, skip_time_domask, skip_time_dopr, & |
---|
[1833] | 283 | skip_time_data_output_av, sloping_surface, & |
---|
[1320] | 284 | stop_dt, terminate_coupled, terminate_run, timestep_scheme, & |
---|
| 285 | time_coupling, time_do2d_xy, time_do2d_xz, time_do2d_yz, & |
---|
| 286 | time_do3d, time_domask, time_dopr, time_dopr_av, & |
---|
| 287 | time_dopr_listing, time_dopts, time_dosp, time_dosp_av, & |
---|
| 288 | time_dots, time_do_av, time_do_sla, time_disturb, time_dvrp, & |
---|
[1496] | 289 | time_run_control, time_since_reference_point, & |
---|
[2050] | 290 | turbulent_inflow, turbulent_outflow, urban_surface, & |
---|
[2011] | 291 | use_initial_profile_as_reference, & |
---|
[1957] | 292 | use_single_reference_value, u_gtrans, v_gtrans, virtual_flight, & |
---|
| 293 | ws_scheme_mom, ws_scheme_sca |
---|
[1320] | 294 | |
---|
| 295 | USE cpulog, & |
---|
| 296 | ONLY: cpu_log, log_point, log_point_s |
---|
| 297 | |
---|
[1957] | 298 | USE flight_mod, & |
---|
| 299 | ONLY: flight_measurement |
---|
| 300 | |
---|
| 301 | |
---|
[1320] | 302 | USE indices, & |
---|
[2118] | 303 | ONLY: nbgp, nx, nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb, nzt, & |
---|
| 304 | nzb_u_inner, nzb_v_inner |
---|
[1320] | 305 | |
---|
| 306 | USE interaction_droplets_ptq_mod, & |
---|
| 307 | ONLY: interaction_droplets_ptq |
---|
| 308 | |
---|
[1918] | 309 | USE interfaces |
---|
| 310 | |
---|
[1320] | 311 | USE kinds |
---|
| 312 | |
---|
[1496] | 313 | USE land_surface_model_mod, & |
---|
[1691] | 314 | ONLY: land_surface, lsm_energy_balance, lsm_soil_model, & |
---|
| 315 | skip_time_do_lsm |
---|
[1496] | 316 | |
---|
[1320] | 317 | USE ls_forcing_mod, & |
---|
| 318 | ONLY: ls_forcing_surf, ls_forcing_vert |
---|
| 319 | |
---|
[1849] | 320 | USE microphysics_mod, & |
---|
| 321 | ONLY: collision_turbulence |
---|
| 322 | |
---|
[1365] | 323 | USE nudge_mod, & |
---|
[1380] | 324 | ONLY: calc_tnudge, nudge_ref |
---|
[1365] | 325 | |
---|
[1320] | 326 | USE particle_attributes, & |
---|
[1671] | 327 | ONLY: particle_advection, particle_advection_start, & |
---|
| 328 | use_sgs_for_particles, wang_kernel |
---|
[1320] | 329 | |
---|
[1] | 330 | USE pegrid |
---|
| 331 | |
---|
[1762] | 332 | USE pmc_interface, & |
---|
[1933] | 333 | ONLY: nested_run, nesting_mode, pmci_datatrans, & |
---|
| 334 | pmci_ensure_nest_mass_conservation, pmci_synchronize |
---|
[1762] | 335 | |
---|
[1320] | 336 | USE production_e_mod, & |
---|
| 337 | ONLY: production_e_init |
---|
| 338 | |
---|
[1402] | 339 | USE progress_bar, & |
---|
| 340 | ONLY: finish_progress_bar, output_progress_bar |
---|
| 341 | |
---|
[1320] | 342 | USE prognostic_equations_mod, & |
---|
[2118] | 343 | ONLY: prognostic_equations_cache, prognostic_equations_vector |
---|
[1320] | 344 | |
---|
[1496] | 345 | USE radiation_model_mod, & |
---|
[1976] | 346 | ONLY: dt_radiation, force_radiation_call, radiation, radiation_control,& |
---|
| 347 | skip_time_do_radiation, time_radiation |
---|
[1496] | 348 | |
---|
[1833] | 349 | USE spectra_mod, & |
---|
| 350 | ONLY: average_count_sp, averaging_interval_sp, calc_spectra, dt_dosp, & |
---|
| 351 | skip_time_dosp |
---|
[1786] | 352 | |
---|
[1320] | 353 | USE statistics, & |
---|
[1918] | 354 | ONLY: flow_statistics_called, hom, pr_palm, sums_ls_l, u_max, & |
---|
| 355 | u_max_ijk, v_max, v_max_ijk, w_max, w_max_ijk |
---|
[1320] | 356 | |
---|
[1691] | 357 | USE surface_layer_fluxes_mod, & |
---|
| 358 | ONLY: surface_layer_fluxes |
---|
| 359 | |
---|
[2007] | 360 | USE urban_surface_mod, & |
---|
[2011] | 361 | ONLY: usm_material_heat_model, usm_material_model, & |
---|
[2007] | 362 | usm_radiation, usm_surface_energy_balance |
---|
| 363 | |
---|
[1320] | 364 | USE user_actions_mod, & |
---|
| 365 | ONLY: user_actions |
---|
| 366 | |
---|
[1914] | 367 | USE wind_turbine_model_mod, & |
---|
| 368 | ONLY: wind_turbine, wtm_forces |
---|
| 369 | |
---|
[1] | 370 | IMPLICIT NONE |
---|
| 371 | |
---|
[1682] | 372 | CHARACTER (LEN=9) :: time_to_string !< |
---|
[1] | 373 | |
---|
[1918] | 374 | REAL(wp) :: dt_3d_old !< temporary storage of timestep to be used for |
---|
| 375 | !< steering of run control output interval |
---|
| 376 | |
---|
[1] | 377 | ! |
---|
[1918] | 378 | !-- At beginning determine the first time step |
---|
| 379 | CALL timestep |
---|
[667] | 380 | |
---|
[1764] | 381 | ! |
---|
| 382 | !-- Synchronize the timestep in case of nested run. |
---|
| 383 | IF ( nested_run ) THEN |
---|
[1878] | 384 | ! |
---|
| 385 | !-- Synchronization by unifying the time step. |
---|
| 386 | !-- Global minimum of all time-steps is used for all. |
---|
| 387 | CALL pmci_synchronize |
---|
[1764] | 388 | ENDIF |
---|
| 389 | |
---|
[1918] | 390 | ! |
---|
| 391 | !-- Determine and print out the run control quantities before the first time |
---|
| 392 | !-- step of this run. For the initial run, some statistics (e.g. divergence) |
---|
| 393 | !-- need to be determined first. |
---|
| 394 | IF ( simulated_time == 0.0_wp ) CALL flow_statistics |
---|
[1] | 395 | CALL run_control |
---|
| 396 | |
---|
[108] | 397 | ! |
---|
| 398 | !-- Data exchange between coupled models in case that a call has been omitted |
---|
| 399 | !-- at the end of the previous run of a job chain. |
---|
[291] | 400 | IF ( coupling_mode /= 'uncoupled' .AND. run_coupled ) THEN |
---|
[108] | 401 | ! |
---|
| 402 | !-- In case of model termination initiated by the local model the coupler |
---|
| 403 | !-- must not be called because this would again cause an MPI hang. |
---|
[1918] | 404 | DO WHILE ( time_coupling >= dt_coupling .AND. terminate_coupled == 0 ) |
---|
[108] | 405 | CALL surface_coupler |
---|
| 406 | time_coupling = time_coupling - dt_coupling |
---|
| 407 | ENDDO |
---|
[1918] | 408 | IF (time_coupling == 0.0_wp .AND. & |
---|
| 409 | time_since_reference_point < dt_coupling ) & |
---|
[348] | 410 | THEN |
---|
| 411 | time_coupling = time_since_reference_point |
---|
| 412 | ENDIF |
---|
[108] | 413 | ENDIF |
---|
| 414 | |
---|
[1] | 415 | #if defined( __dvrp_graphics ) |
---|
| 416 | ! |
---|
| 417 | !-- Time measurement with dvrp software |
---|
| 418 | CALL DVRP_LOG_EVENT( 2, current_timestep_number ) |
---|
| 419 | #endif |
---|
| 420 | |
---|
[1402] | 421 | CALL location_message( 'start with time-stepping', .TRUE. ) |
---|
[1] | 422 | ! |
---|
| 423 | !-- Start of the time loop |
---|
| 424 | DO WHILE ( simulated_time < end_time .AND. .NOT. stop_dt .AND. & |
---|
| 425 | .NOT. terminate_run ) |
---|
| 426 | |
---|
| 427 | CALL cpu_log( log_point_s(10), 'timesteps', 'start' ) |
---|
[1221] | 428 | |
---|
[1] | 429 | ! |
---|
[1241] | 430 | !-- Determine ug, vg and w_subs in dependence on data from external file |
---|
| 431 | !-- LSF_DATA |
---|
[1365] | 432 | IF ( large_scale_forcing .AND. lsf_vert ) THEN |
---|
[1241] | 433 | CALL ls_forcing_vert ( simulated_time ) |
---|
[1365] | 434 | sums_ls_l = 0.0_wp |
---|
[1241] | 435 | ENDIF |
---|
| 436 | |
---|
| 437 | ! |
---|
[1380] | 438 | !-- Set pt_init and q_init to the current profiles taken from |
---|
| 439 | !-- NUDGING_DATA |
---|
| 440 | IF ( nudging ) THEN |
---|
| 441 | CALL nudge_ref ( simulated_time ) |
---|
| 442 | ! |
---|
| 443 | !-- Store temperature gradient at the top boundary for possible Neumann |
---|
| 444 | !-- boundary condition |
---|
| 445 | bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1) |
---|
| 446 | bc_q_t_val = ( q_init(nzt+1) - q_init(nzt) ) / dzu(nzt+1) |
---|
| 447 | ENDIF |
---|
| 448 | |
---|
| 449 | ! |
---|
[1] | 450 | !-- Execute the user-defined actions |
---|
| 451 | CALL user_actions( 'before_timestep' ) |
---|
| 452 | |
---|
| 453 | ! |
---|
[1914] | 454 | !-- Calculate forces by wind turbines |
---|
| 455 | IF ( wind_turbine ) THEN |
---|
| 456 | |
---|
| 457 | CALL cpu_log( log_point(55), 'wind_turbine', 'start' ) |
---|
| 458 | |
---|
| 459 | CALL wtm_forces |
---|
| 460 | |
---|
| 461 | CALL cpu_log( log_point(55), 'wind_turbine', 'stop' ) |
---|
| 462 | |
---|
| 463 | ENDIF |
---|
| 464 | |
---|
| 465 | ! |
---|
[1] | 466 | !-- Start of intermediate step loop |
---|
| 467 | intermediate_timestep_count = 0 |
---|
| 468 | DO WHILE ( intermediate_timestep_count < & |
---|
| 469 | intermediate_timestep_count_max ) |
---|
| 470 | |
---|
| 471 | intermediate_timestep_count = intermediate_timestep_count + 1 |
---|
| 472 | |
---|
| 473 | ! |
---|
| 474 | !-- Set the steering factors for the prognostic equations which depend |
---|
| 475 | !-- on the timestep scheme |
---|
| 476 | CALL timestep_scheme_steering |
---|
| 477 | |
---|
| 478 | ! |
---|
[1128] | 479 | !-- Calculate those variables needed in the tendency terms which need |
---|
| 480 | !-- global communication |
---|
[1179] | 481 | IF ( .NOT. use_single_reference_value .AND. & |
---|
| 482 | .NOT. use_initial_profile_as_reference ) THEN |
---|
| 483 | ! |
---|
| 484 | !-- Horizontally averaged profiles to be used as reference state in |
---|
| 485 | !-- buoyancy terms (WARNING: only the respective last call of |
---|
| 486 | !-- calc_mean_profile defines the reference state!) |
---|
[1365] | 487 | IF ( .NOT. neutral ) THEN |
---|
| 488 | CALL calc_mean_profile( pt, 4 ) |
---|
| 489 | ref_state(:) = hom(:,1,4,0) ! this is used in the buoyancy term |
---|
| 490 | ENDIF |
---|
| 491 | IF ( ocean ) THEN |
---|
[2031] | 492 | CALL calc_mean_profile( rho_ocean, 64 ) |
---|
[1365] | 493 | ref_state(:) = hom(:,1,64,0) |
---|
| 494 | ENDIF |
---|
| 495 | IF ( humidity ) THEN |
---|
| 496 | CALL calc_mean_profile( vpt, 44 ) |
---|
| 497 | ref_state(:) = hom(:,1,44,0) |
---|
| 498 | ENDIF |
---|
| 499 | |
---|
[1179] | 500 | ENDIF |
---|
| 501 | |
---|
[1128] | 502 | IF ( .NOT. constant_diffusion ) CALL production_e_init |
---|
| 503 | IF ( ( ws_scheme_mom .OR. ws_scheme_sca ) .AND. & |
---|
| 504 | intermediate_timestep_count == 1 ) CALL ws_statistics |
---|
[1365] | 505 | ! |
---|
| 506 | !-- In case of nudging calculate current nudging time scale and horizontal |
---|
[1380] | 507 | !-- means of u, v, pt and q |
---|
[1365] | 508 | IF ( nudging ) THEN |
---|
| 509 | CALL calc_tnudge( simulated_time ) |
---|
| 510 | CALL calc_mean_profile( u, 1 ) |
---|
| 511 | CALL calc_mean_profile( v, 2 ) |
---|
| 512 | CALL calc_mean_profile( pt, 4 ) |
---|
| 513 | CALL calc_mean_profile( q, 41 ) |
---|
| 514 | ENDIF |
---|
[1128] | 515 | |
---|
| 516 | ! |
---|
[1] | 517 | !-- Solve the prognostic equations. A fast cache optimized version with |
---|
| 518 | !-- only one single loop is used in case of Piascek-Williams advection |
---|
| 519 | !-- scheme. NEC vector machines use a different version, because |
---|
| 520 | !-- in the other versions a good vectorization is prohibited due to |
---|
| 521 | !-- inlining problems. |
---|
[1019] | 522 | IF ( loop_optimization == 'cache' ) THEN |
---|
| 523 | CALL prognostic_equations_cache |
---|
| 524 | ELSEIF ( loop_optimization == 'vector' ) THEN |
---|
[63] | 525 | CALL prognostic_equations_vector |
---|
[1] | 526 | ENDIF |
---|
| 527 | |
---|
| 528 | ! |
---|
[849] | 529 | !-- Particle transport/physics with the Lagrangian particle model |
---|
| 530 | !-- (only once during intermediate steps, because it uses an Euler-step) |
---|
[1128] | 531 | !-- ### particle model should be moved before prognostic_equations, in order |
---|
| 532 | !-- to regard droplet interactions directly |
---|
[63] | 533 | IF ( particle_advection .AND. & |
---|
| 534 | simulated_time >= particle_advection_start .AND. & |
---|
[1] | 535 | intermediate_timestep_count == 1 ) THEN |
---|
[849] | 536 | CALL lpm |
---|
| 537 | first_call_lpm = .FALSE. |
---|
[1] | 538 | ENDIF |
---|
| 539 | |
---|
| 540 | ! |
---|
| 541 | !-- Interaction of droplets with temperature and specific humidity. |
---|
| 542 | !-- Droplet condensation and evaporation is calculated within |
---|
| 543 | !-- advec_particles. |
---|
| 544 | IF ( cloud_droplets .AND. & |
---|
| 545 | intermediate_timestep_count == intermediate_timestep_count_max )& |
---|
| 546 | THEN |
---|
| 547 | CALL interaction_droplets_ptq |
---|
| 548 | ENDIF |
---|
| 549 | |
---|
| 550 | ! |
---|
| 551 | !-- Exchange of ghost points (lateral boundary conditions) |
---|
[2118] | 552 | CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'start' ) |
---|
[1113] | 553 | |
---|
[2118] | 554 | CALL exchange_horiz( u_p, nbgp ) |
---|
| 555 | CALL exchange_horiz( v_p, nbgp ) |
---|
| 556 | CALL exchange_horiz( w_p, nbgp ) |
---|
| 557 | CALL exchange_horiz( pt_p, nbgp ) |
---|
| 558 | IF ( .NOT. constant_diffusion ) CALL exchange_horiz( e_p, nbgp ) |
---|
| 559 | IF ( ocean ) THEN |
---|
| 560 | CALL exchange_horiz( sa_p, nbgp ) |
---|
| 561 | CALL exchange_horiz( rho_ocean, nbgp ) |
---|
| 562 | CALL exchange_horiz( prho, nbgp ) |
---|
| 563 | ENDIF |
---|
| 564 | IF ( humidity ) THEN |
---|
| 565 | CALL exchange_horiz( q_p, nbgp ) |
---|
| 566 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
| 567 | CALL exchange_horiz( qr_p, nbgp ) |
---|
| 568 | CALL exchange_horiz( nr_p, nbgp ) |
---|
[1053] | 569 | ENDIF |
---|
[2118] | 570 | ENDIF |
---|
| 571 | IF ( cloud_droplets ) THEN |
---|
| 572 | CALL exchange_horiz( ql, nbgp ) |
---|
| 573 | CALL exchange_horiz( ql_c, nbgp ) |
---|
| 574 | CALL exchange_horiz( ql_v, nbgp ) |
---|
| 575 | CALL exchange_horiz( ql_vp, nbgp ) |
---|
| 576 | ENDIF |
---|
| 577 | IF ( wang_kernel .OR. collision_turbulence .OR. & |
---|
| 578 | use_sgs_for_particles ) THEN |
---|
| 579 | CALL exchange_horiz( diss, nbgp ) |
---|
| 580 | ENDIF |
---|
| 581 | IF ( passive_scalar ) CALL exchange_horiz( s_p, nbgp ) |
---|
[1] | 582 | |
---|
[2118] | 583 | CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'stop' ) |
---|
[1128] | 584 | |
---|
[1] | 585 | ! |
---|
| 586 | !-- Boundary conditions for the prognostic quantities (except of the |
---|
| 587 | !-- velocities at the outflow in case of a non-cyclic lateral wall) |
---|
[1113] | 588 | CALL boundary_conds |
---|
[1] | 589 | |
---|
| 590 | ! |
---|
[73] | 591 | !-- Swap the time levels in preparation for the next time step. |
---|
| 592 | CALL swap_timelevel |
---|
| 593 | |
---|
[1764] | 594 | IF ( nested_run ) THEN |
---|
[1797] | 595 | |
---|
[1764] | 596 | CALL cpu_log( log_point(60), 'nesting', 'start' ) |
---|
[1762] | 597 | ! |
---|
[1933] | 598 | !-- Domain nesting. The data transfer subroutines pmci_parent_datatrans |
---|
| 599 | !-- and pmci_child_datatrans are called inside the wrapper |
---|
[1797] | 600 | !-- subroutine pmci_datatrans according to the control parameters |
---|
| 601 | !-- nesting_mode and nesting_datatransfer_mode. |
---|
| 602 | !-- TO_DO: why is nesting_mode given as a parameter here? |
---|
| 603 | CALL pmci_datatrans( nesting_mode ) |
---|
[1762] | 604 | |
---|
[1933] | 605 | IF ( TRIM( nesting_mode ) == 'two-way' .OR. & |
---|
| 606 | nesting_mode == 'vertical' ) THEN |
---|
[1762] | 607 | ! |
---|
[1933] | 608 | !-- Exchange_horiz is needed for all parent-domains after the |
---|
[1764] | 609 | !-- anterpolation |
---|
| 610 | CALL exchange_horiz( u, nbgp ) |
---|
| 611 | CALL exchange_horiz( v, nbgp ) |
---|
| 612 | CALL exchange_horiz( w, nbgp ) |
---|
[1960] | 613 | IF ( .NOT. neutral ) CALL exchange_horiz( pt, nbgp ) |
---|
| 614 | IF ( humidity ) CALL exchange_horiz( q, nbgp ) |
---|
| 615 | IF ( passive_scalar ) CALL exchange_horiz( s, nbgp ) |
---|
| 616 | IF ( .NOT. constant_diffusion ) CALL exchange_horiz( e, nbgp ) |
---|
[1762] | 617 | ENDIF |
---|
| 618 | ! |
---|
[1764] | 619 | !-- Correct the w top-BC in nest domains to ensure mass conservation. |
---|
[1933] | 620 | !-- This action must never be done for the root domain. Vertical |
---|
| 621 | !-- nesting implies mass conservation. |
---|
[1764] | 622 | IF ( nest_domain ) THEN |
---|
| 623 | CALL pmci_ensure_nest_mass_conservation |
---|
| 624 | ENDIF |
---|
| 625 | |
---|
| 626 | CALL cpu_log( log_point(60), 'nesting', 'stop' ) |
---|
| 627 | |
---|
[1762] | 628 | ENDIF |
---|
| 629 | |
---|
| 630 | ! |
---|
[1] | 631 | !-- Temperature offset must be imposed at cyclic boundaries in x-direction |
---|
| 632 | !-- when a sloping surface is used |
---|
| 633 | IF ( sloping_surface ) THEN |
---|
[707] | 634 | IF ( nxl == 0 ) pt(:,:,nxlg:nxl-1) = pt(:,:,nxlg:nxl-1) - & |
---|
| 635 | pt_slope_offset |
---|
| 636 | IF ( nxr == nx ) pt(:,:,nxr+1:nxrg) = pt(:,:,nxr+1:nxrg) + & |
---|
| 637 | pt_slope_offset |
---|
[1] | 638 | ENDIF |
---|
| 639 | |
---|
| 640 | ! |
---|
[151] | 641 | !-- Impose a turbulent inflow using the recycling method |
---|
| 642 | IF ( turbulent_inflow ) CALL inflow_turbulence |
---|
| 643 | |
---|
| 644 | ! |
---|
[2050] | 645 | !-- Set values at outflow boundary using the special outflow condition |
---|
| 646 | IF ( turbulent_outflow ) CALL outflow_turbulence |
---|
| 647 | |
---|
| 648 | ! |
---|
[1] | 649 | !-- Impose a random perturbation on the horizontal velocity field |
---|
[106] | 650 | IF ( create_disturbances .AND. & |
---|
| 651 | ( call_psolver_at_all_substeps .AND. & |
---|
[1] | 652 | intermediate_timestep_count == intermediate_timestep_count_max )& |
---|
[106] | 653 | .OR. ( .NOT. call_psolver_at_all_substeps .AND. & |
---|
| 654 | intermediate_timestep_count == 1 ) ) & |
---|
[1] | 655 | THEN |
---|
| 656 | time_disturb = time_disturb + dt_3d |
---|
| 657 | IF ( time_disturb >= dt_disturb ) THEN |
---|
[1736] | 658 | IF ( disturbance_energy_limit /= 0.0_wp .AND. & |
---|
| 659 | hom(nzb+5,1,pr_palm,0) < disturbance_energy_limit ) THEN |
---|
[75] | 660 | CALL disturb_field( nzb_u_inner, tend, u ) |
---|
| 661 | CALL disturb_field( nzb_v_inner, tend, v ) |
---|
[707] | 662 | ELSEIF ( .NOT. bc_lr_cyc .OR. .NOT. bc_ns_cyc ) THEN |
---|
[1] | 663 | ! |
---|
| 664 | !-- Runs with a non-cyclic lateral wall need perturbations |
---|
| 665 | !-- near the inflow throughout the whole simulation |
---|
| 666 | dist_range = 1 |
---|
[75] | 667 | CALL disturb_field( nzb_u_inner, tend, u ) |
---|
| 668 | CALL disturb_field( nzb_v_inner, tend, v ) |
---|
[1] | 669 | dist_range = 0 |
---|
| 670 | ENDIF |
---|
| 671 | time_disturb = time_disturb - dt_disturb |
---|
| 672 | ENDIF |
---|
| 673 | ENDIF |
---|
| 674 | |
---|
| 675 | ! |
---|
| 676 | !-- Reduce the velocity divergence via the equation for perturbation |
---|
| 677 | !-- pressure. |
---|
[106] | 678 | IF ( intermediate_timestep_count == 1 .OR. & |
---|
| 679 | call_psolver_at_all_substeps ) THEN |
---|
[1] | 680 | CALL pres |
---|
| 681 | ENDIF |
---|
| 682 | |
---|
| 683 | ! |
---|
| 684 | !-- If required, compute liquid water content |
---|
[1015] | 685 | IF ( cloud_physics ) THEN |
---|
| 686 | CALL calc_liquid_water_content |
---|
| 687 | ENDIF |
---|
[1115] | 688 | ! |
---|
| 689 | !-- If required, compute virtual potential temperature |
---|
| 690 | IF ( humidity ) THEN |
---|
| 691 | CALL compute_vpt |
---|
| 692 | ENDIF |
---|
[1585] | 693 | |
---|
[1] | 694 | ! |
---|
| 695 | !-- Compute the diffusion quantities |
---|
| 696 | IF ( .NOT. constant_diffusion ) THEN |
---|
| 697 | |
---|
| 698 | ! |
---|
[1276] | 699 | !-- Determine surface fluxes shf and qsws and surface values |
---|
| 700 | !-- pt_surface and q_surface in dependence on data from external |
---|
| 701 | !-- file LSF_DATA respectively |
---|
| 702 | IF ( ( large_scale_forcing .AND. lsf_surf ) .AND. & |
---|
| 703 | intermediate_timestep_count == intermediate_timestep_count_max )& |
---|
| 704 | THEN |
---|
| 705 | CALL ls_forcing_surf ( simulated_time ) |
---|
| 706 | ENDIF |
---|
| 707 | |
---|
| 708 | ! |
---|
[1691] | 709 | !-- First the vertical fluxes in the surface (constant flux) layer are computed |
---|
| 710 | IF ( constant_flux_layer ) THEN |
---|
| 711 | CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'start' ) |
---|
| 712 | CALL surface_layer_fluxes |
---|
| 713 | CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'stop' ) |
---|
[1] | 714 | ENDIF |
---|
[1241] | 715 | |
---|
[1] | 716 | ! |
---|
[1691] | 717 | !-- If required, solve the energy balance for the surface and run soil |
---|
| 718 | !-- model |
---|
| 719 | IF ( land_surface .AND. simulated_time > skip_time_do_lsm) THEN |
---|
| 720 | |
---|
| 721 | CALL cpu_log( log_point(54), 'land_surface', 'start' ) |
---|
| 722 | CALL lsm_energy_balance |
---|
| 723 | CALL lsm_soil_model |
---|
| 724 | CALL cpu_log( log_point(54), 'land_surface', 'stop' ) |
---|
| 725 | ENDIF |
---|
[2007] | 726 | |
---|
[1691] | 727 | ! |
---|
[2007] | 728 | !-- If required, solve the energy balance for urban surfaces and run |
---|
| 729 | !-- the material heat model |
---|
| 730 | IF (urban_surface) THEN |
---|
| 731 | CALL cpu_log( log_point(74), 'urban_surface', 'start' ) |
---|
| 732 | CALL usm_surface_energy_balance |
---|
| 733 | IF ( usm_material_model ) THEN |
---|
| 734 | CALL usm_material_heat_model |
---|
| 735 | ENDIF |
---|
| 736 | CALL cpu_log( log_point(74), 'urban_surface', 'stop' ) |
---|
| 737 | ENDIF |
---|
| 738 | |
---|
| 739 | ! |
---|
[1] | 740 | !-- Compute the diffusion coefficients |
---|
| 741 | CALL cpu_log( log_point(17), 'diffusivities', 'start' ) |
---|
[75] | 742 | IF ( .NOT. humidity ) THEN |
---|
[97] | 743 | IF ( ocean ) THEN |
---|
[388] | 744 | CALL diffusivities( prho, prho_reference ) |
---|
[97] | 745 | ELSE |
---|
| 746 | CALL diffusivities( pt, pt_reference ) |
---|
| 747 | ENDIF |
---|
[1] | 748 | ELSE |
---|
[97] | 749 | CALL diffusivities( vpt, pt_reference ) |
---|
[1] | 750 | ENDIF |
---|
| 751 | CALL cpu_log( log_point(17), 'diffusivities', 'stop' ) |
---|
| 752 | |
---|
| 753 | ENDIF |
---|
| 754 | |
---|
[1691] | 755 | ! |
---|
| 756 | !-- If required, calculate radiative fluxes and heating rates |
---|
| 757 | IF ( radiation .AND. intermediate_timestep_count & |
---|
| 758 | == intermediate_timestep_count_max .AND. simulated_time > & |
---|
| 759 | skip_time_do_radiation ) THEN |
---|
| 760 | |
---|
| 761 | time_radiation = time_radiation + dt_3d |
---|
| 762 | |
---|
| 763 | IF ( time_radiation >= dt_radiation .OR. force_radiation_call ) & |
---|
| 764 | THEN |
---|
| 765 | |
---|
| 766 | CALL cpu_log( log_point(50), 'radiation', 'start' ) |
---|
| 767 | |
---|
| 768 | IF ( .NOT. force_radiation_call ) THEN |
---|
| 769 | time_radiation = time_radiation - dt_radiation |
---|
| 770 | ENDIF |
---|
| 771 | |
---|
[1976] | 772 | CALL radiation_control |
---|
[1691] | 773 | |
---|
| 774 | CALL cpu_log( log_point(50), 'radiation', 'stop' ) |
---|
[2007] | 775 | |
---|
| 776 | IF (urban_surface) THEN |
---|
| 777 | CALL cpu_log( log_point(75), 'usm_radiation', 'start' ) |
---|
| 778 | CALL usm_radiation |
---|
| 779 | CALL cpu_log( log_point(75), 'usm_radiation', 'stop' ) |
---|
| 780 | ENDIF |
---|
| 781 | |
---|
[1691] | 782 | ENDIF |
---|
| 783 | ENDIF |
---|
| 784 | |
---|
[1] | 785 | ENDDO ! Intermediate step loop |
---|
| 786 | |
---|
| 787 | ! |
---|
| 788 | !-- Increase simulation time and output times |
---|
[1111] | 789 | nr_timesteps_this_run = nr_timesteps_this_run + 1 |
---|
[291] | 790 | current_timestep_number = current_timestep_number + 1 |
---|
| 791 | simulated_time = simulated_time + dt_3d |
---|
| 792 | simulated_time_chr = time_to_string( simulated_time ) |
---|
| 793 | time_since_reference_point = simulated_time - coupling_start_time |
---|
| 794 | |
---|
[1957] | 795 | |
---|
| 796 | |
---|
[1] | 797 | IF ( simulated_time >= skip_time_data_output_av ) THEN |
---|
| 798 | time_do_av = time_do_av + dt_3d |
---|
| 799 | ENDIF |
---|
| 800 | IF ( simulated_time >= skip_time_do2d_xy ) THEN |
---|
| 801 | time_do2d_xy = time_do2d_xy + dt_3d |
---|
| 802 | ENDIF |
---|
| 803 | IF ( simulated_time >= skip_time_do2d_xz ) THEN |
---|
| 804 | time_do2d_xz = time_do2d_xz + dt_3d |
---|
| 805 | ENDIF |
---|
| 806 | IF ( simulated_time >= skip_time_do2d_yz ) THEN |
---|
| 807 | time_do2d_yz = time_do2d_yz + dt_3d |
---|
| 808 | ENDIF |
---|
| 809 | IF ( simulated_time >= skip_time_do3d ) THEN |
---|
| 810 | time_do3d = time_do3d + dt_3d |
---|
| 811 | ENDIF |
---|
[410] | 812 | DO mid = 1, masks |
---|
| 813 | IF ( simulated_time >= skip_time_domask(mid) ) THEN |
---|
| 814 | time_domask(mid)= time_domask(mid) + dt_3d |
---|
| 815 | ENDIF |
---|
| 816 | ENDDO |
---|
[1] | 817 | time_dvrp = time_dvrp + dt_3d |
---|
| 818 | IF ( simulated_time >= skip_time_dosp ) THEN |
---|
| 819 | time_dosp = time_dosp + dt_3d |
---|
| 820 | ENDIF |
---|
| 821 | time_dots = time_dots + dt_3d |
---|
[849] | 822 | IF ( .NOT. first_call_lpm ) THEN |
---|
[1] | 823 | time_dopts = time_dopts + dt_3d |
---|
| 824 | ENDIF |
---|
| 825 | IF ( simulated_time >= skip_time_dopr ) THEN |
---|
| 826 | time_dopr = time_dopr + dt_3d |
---|
| 827 | ENDIF |
---|
| 828 | time_dopr_listing = time_dopr_listing + dt_3d |
---|
| 829 | time_run_control = time_run_control + dt_3d |
---|
| 830 | |
---|
| 831 | ! |
---|
[102] | 832 | !-- Data exchange between coupled models |
---|
[291] | 833 | IF ( coupling_mode /= 'uncoupled' .AND. run_coupled ) THEN |
---|
[102] | 834 | time_coupling = time_coupling + dt_3d |
---|
[343] | 835 | |
---|
[108] | 836 | ! |
---|
| 837 | !-- In case of model termination initiated by the local model |
---|
| 838 | !-- (terminate_coupled > 0), the coupler must be skipped because it would |
---|
| 839 | !-- cause an MPI intercomminucation hang. |
---|
| 840 | !-- If necessary, the coupler will be called at the beginning of the |
---|
| 841 | !-- next restart run. |
---|
| 842 | DO WHILE ( time_coupling >= dt_coupling .AND. terminate_coupled == 0 ) |
---|
[102] | 843 | CALL surface_coupler |
---|
| 844 | time_coupling = time_coupling - dt_coupling |
---|
| 845 | ENDDO |
---|
| 846 | ENDIF |
---|
| 847 | |
---|
| 848 | ! |
---|
[46] | 849 | !-- Execute user-defined actions |
---|
| 850 | CALL user_actions( 'after_integration' ) |
---|
| 851 | |
---|
| 852 | ! |
---|
[1] | 853 | !-- If Galilei transformation is used, determine the distance that the |
---|
| 854 | !-- model has moved so far |
---|
| 855 | IF ( galilei_transformation ) THEN |
---|
| 856 | advected_distance_x = advected_distance_x + u_gtrans * dt_3d |
---|
| 857 | advected_distance_y = advected_distance_y + v_gtrans * dt_3d |
---|
| 858 | ENDIF |
---|
| 859 | |
---|
| 860 | ! |
---|
| 861 | !-- Check, if restart is necessary (because cpu-time is expiring or |
---|
| 862 | !-- because it is forced by user) and set stop flag |
---|
[108] | 863 | !-- This call is skipped if the remote model has already initiated a restart. |
---|
| 864 | IF ( .NOT. terminate_run ) CALL check_for_restart |
---|
[1] | 865 | |
---|
| 866 | ! |
---|
| 867 | !-- Carry out statistical analysis and output at the requested output times. |
---|
| 868 | !-- The MOD function is used for calculating the output time counters (like |
---|
| 869 | !-- time_dopr) in order to regard a possible decrease of the output time |
---|
| 870 | !-- interval in case of restart runs |
---|
| 871 | |
---|
| 872 | ! |
---|
| 873 | !-- Set a flag indicating that so far no statistics have been created |
---|
| 874 | !-- for this time step |
---|
| 875 | flow_statistics_called = .FALSE. |
---|
| 876 | |
---|
| 877 | ! |
---|
| 878 | !-- If required, call flow_statistics for averaging in time |
---|
[1342] | 879 | IF ( averaging_interval_pr /= 0.0_wp .AND. & |
---|
[1] | 880 | ( dt_dopr - time_dopr ) <= averaging_interval_pr .AND. & |
---|
| 881 | simulated_time >= skip_time_dopr ) THEN |
---|
| 882 | time_dopr_av = time_dopr_av + dt_3d |
---|
| 883 | IF ( time_dopr_av >= dt_averaging_input_pr ) THEN |
---|
| 884 | do_sum = .TRUE. |
---|
| 885 | time_dopr_av = MOD( time_dopr_av, & |
---|
| 886 | MAX( dt_averaging_input_pr, dt_3d ) ) |
---|
| 887 | ENDIF |
---|
| 888 | ENDIF |
---|
| 889 | IF ( do_sum ) CALL flow_statistics |
---|
| 890 | |
---|
| 891 | ! |
---|
[410] | 892 | !-- Sum-up 3d-arrays for later output of time-averaged 2d/3d/masked data |
---|
[1342] | 893 | IF ( averaging_interval /= 0.0_wp .AND. & |
---|
[1] | 894 | ( dt_data_output_av - time_do_av ) <= averaging_interval .AND. & |
---|
| 895 | simulated_time >= skip_time_data_output_av ) & |
---|
| 896 | THEN |
---|
| 897 | time_do_sla = time_do_sla + dt_3d |
---|
| 898 | IF ( time_do_sla >= dt_averaging_input ) THEN |
---|
| 899 | CALL sum_up_3d_data |
---|
| 900 | average_count_3d = average_count_3d + 1 |
---|
| 901 | time_do_sla = MOD( time_do_sla, MAX( dt_averaging_input, dt_3d ) ) |
---|
| 902 | ENDIF |
---|
| 903 | ENDIF |
---|
| 904 | |
---|
| 905 | ! |
---|
| 906 | !-- Calculate spectra for time averaging |
---|
[1342] | 907 | IF ( averaging_interval_sp /= 0.0_wp .AND. & |
---|
[1] | 908 | ( dt_dosp - time_dosp ) <= averaging_interval_sp .AND. & |
---|
| 909 | simulated_time >= skip_time_dosp ) THEN |
---|
| 910 | time_dosp_av = time_dosp_av + dt_3d |
---|
| 911 | IF ( time_dosp_av >= dt_averaging_input_pr ) THEN |
---|
| 912 | CALL calc_spectra |
---|
| 913 | time_dosp_av = MOD( time_dosp_av, & |
---|
| 914 | MAX( dt_averaging_input_pr, dt_3d ) ) |
---|
| 915 | ENDIF |
---|
| 916 | ENDIF |
---|
| 917 | |
---|
| 918 | ! |
---|
[1957] | 919 | !-- Call flight module and output data |
---|
| 920 | IF ( virtual_flight ) THEN |
---|
| 921 | CALL flight_measurement |
---|
| 922 | CALL data_output_flight |
---|
| 923 | ENDIF |
---|
| 924 | |
---|
| 925 | ! |
---|
[1] | 926 | !-- Profile output (ASCII) on file |
---|
| 927 | IF ( time_dopr_listing >= dt_dopr_listing ) THEN |
---|
| 928 | CALL print_1d |
---|
| 929 | time_dopr_listing = MOD( time_dopr_listing, MAX( dt_dopr_listing, & |
---|
| 930 | dt_3d ) ) |
---|
| 931 | ENDIF |
---|
| 932 | |
---|
| 933 | ! |
---|
| 934 | !-- Graphic output for PROFIL |
---|
| 935 | IF ( time_dopr >= dt_dopr ) THEN |
---|
| 936 | IF ( dopr_n /= 0 ) CALL data_output_profiles |
---|
| 937 | time_dopr = MOD( time_dopr, MAX( dt_dopr, dt_3d ) ) |
---|
[1342] | 938 | time_dopr_av = 0.0_wp ! due to averaging (see above) |
---|
[1] | 939 | ENDIF |
---|
| 940 | |
---|
| 941 | ! |
---|
| 942 | !-- Graphic output for time series |
---|
| 943 | IF ( time_dots >= dt_dots ) THEN |
---|
[48] | 944 | CALL data_output_tseries |
---|
[1] | 945 | time_dots = MOD( time_dots, MAX( dt_dots, dt_3d ) ) |
---|
| 946 | ENDIF |
---|
| 947 | |
---|
| 948 | ! |
---|
| 949 | !-- Output of spectra (formatted for use with PROFIL), in case of no |
---|
| 950 | !-- time averaging, spectra has to be calculated before |
---|
| 951 | IF ( time_dosp >= dt_dosp ) THEN |
---|
| 952 | IF ( average_count_sp == 0 ) CALL calc_spectra |
---|
| 953 | CALL data_output_spectra |
---|
| 954 | time_dosp = MOD( time_dosp, MAX( dt_dosp, dt_3d ) ) |
---|
| 955 | ENDIF |
---|
| 956 | |
---|
| 957 | ! |
---|
| 958 | !-- 2d-data output (cross-sections) |
---|
| 959 | IF ( time_do2d_xy >= dt_do2d_xy ) THEN |
---|
| 960 | CALL data_output_2d( 'xy', 0 ) |
---|
| 961 | time_do2d_xy = MOD( time_do2d_xy, MAX( dt_do2d_xy, dt_3d ) ) |
---|
| 962 | ENDIF |
---|
| 963 | IF ( time_do2d_xz >= dt_do2d_xz ) THEN |
---|
| 964 | CALL data_output_2d( 'xz', 0 ) |
---|
| 965 | time_do2d_xz = MOD( time_do2d_xz, MAX( dt_do2d_xz, dt_3d ) ) |
---|
| 966 | ENDIF |
---|
| 967 | IF ( time_do2d_yz >= dt_do2d_yz ) THEN |
---|
| 968 | CALL data_output_2d( 'yz', 0 ) |
---|
| 969 | time_do2d_yz = MOD( time_do2d_yz, MAX( dt_do2d_yz, dt_3d ) ) |
---|
| 970 | ENDIF |
---|
| 971 | |
---|
| 972 | ! |
---|
| 973 | !-- 3d-data output (volume data) |
---|
| 974 | IF ( time_do3d >= dt_do3d ) THEN |
---|
| 975 | CALL data_output_3d( 0 ) |
---|
| 976 | time_do3d = MOD( time_do3d, MAX( dt_do3d, dt_3d ) ) |
---|
| 977 | ENDIF |
---|
| 978 | |
---|
| 979 | ! |
---|
[1783] | 980 | !-- Masked data output |
---|
[410] | 981 | DO mid = 1, masks |
---|
| 982 | IF ( time_domask(mid) >= dt_domask(mid) ) THEN |
---|
| 983 | CALL data_output_mask( 0 ) |
---|
| 984 | time_domask(mid) = MOD( time_domask(mid), & |
---|
| 985 | MAX( dt_domask(mid), dt_3d ) ) |
---|
| 986 | ENDIF |
---|
| 987 | ENDDO |
---|
| 988 | |
---|
| 989 | ! |
---|
| 990 | !-- Output of time-averaged 2d/3d/masked data |
---|
[1] | 991 | IF ( time_do_av >= dt_data_output_av ) THEN |
---|
| 992 | CALL average_3d_data |
---|
| 993 | CALL data_output_2d( 'xy', 1 ) |
---|
| 994 | CALL data_output_2d( 'xz', 1 ) |
---|
| 995 | CALL data_output_2d( 'yz', 1 ) |
---|
| 996 | CALL data_output_3d( 1 ) |
---|
[410] | 997 | DO mid = 1, masks |
---|
| 998 | CALL data_output_mask( 1 ) |
---|
| 999 | ENDDO |
---|
[1] | 1000 | time_do_av = MOD( time_do_av, MAX( dt_data_output_av, dt_3d ) ) |
---|
| 1001 | ENDIF |
---|
| 1002 | |
---|
| 1003 | ! |
---|
| 1004 | !-- Output of particle time series |
---|
[253] | 1005 | IF ( particle_advection ) THEN |
---|
| 1006 | IF ( time_dopts >= dt_dopts .OR. & |
---|
| 1007 | ( simulated_time >= particle_advection_start .AND. & |
---|
[849] | 1008 | first_call_lpm ) ) THEN |
---|
[253] | 1009 | CALL data_output_ptseries |
---|
| 1010 | time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) ) |
---|
| 1011 | ENDIF |
---|
[1] | 1012 | ENDIF |
---|
| 1013 | |
---|
| 1014 | ! |
---|
| 1015 | !-- Output of dvrp-graphics (isosurface, particles, slicer) |
---|
| 1016 | #if defined( __dvrp_graphics ) |
---|
| 1017 | CALL DVRP_LOG_EVENT( -2, current_timestep_number-1 ) |
---|
| 1018 | #endif |
---|
| 1019 | IF ( time_dvrp >= dt_dvrp ) THEN |
---|
| 1020 | CALL data_output_dvrp |
---|
| 1021 | time_dvrp = MOD( time_dvrp, MAX( dt_dvrp, dt_3d ) ) |
---|
| 1022 | ENDIF |
---|
| 1023 | #if defined( __dvrp_graphics ) |
---|
| 1024 | CALL DVRP_LOG_EVENT( 2, current_timestep_number ) |
---|
| 1025 | #endif |
---|
| 1026 | |
---|
| 1027 | ! |
---|
| 1028 | !-- If required, set the heat flux for the next time step at a random value |
---|
| 1029 | IF ( constant_heatflux .AND. random_heatflux ) CALL disturb_heatflux |
---|
| 1030 | |
---|
| 1031 | ! |
---|
| 1032 | !-- Execute user-defined actions |
---|
| 1033 | CALL user_actions( 'after_timestep' ) |
---|
| 1034 | |
---|
[1402] | 1035 | ! |
---|
[1918] | 1036 | !-- Determine size of next time step. Save timestep dt_3d because it is |
---|
| 1037 | !-- newly calculated in routine timestep, but required further below for |
---|
| 1038 | !-- steering the run control output interval |
---|
| 1039 | dt_3d_old = dt_3d |
---|
| 1040 | CALL timestep |
---|
| 1041 | |
---|
| 1042 | ! |
---|
[1925] | 1043 | !-- Synchronize the timestep in case of nested run. |
---|
| 1044 | IF ( nested_run ) THEN |
---|
| 1045 | ! |
---|
| 1046 | !-- Synchronize by unifying the time step. |
---|
| 1047 | !-- Global minimum of all time-steps is used for all. |
---|
| 1048 | CALL pmci_synchronize |
---|
| 1049 | ENDIF |
---|
| 1050 | |
---|
| 1051 | ! |
---|
[1918] | 1052 | !-- Computation and output of run control parameters. |
---|
| 1053 | !-- This is also done whenever perturbations have been imposed |
---|
| 1054 | IF ( time_run_control >= dt_run_control .OR. & |
---|
| 1055 | timestep_scheme(1:5) /= 'runge' .OR. disturbance_created ) & |
---|
| 1056 | THEN |
---|
| 1057 | CALL run_control |
---|
| 1058 | IF ( time_run_control >= dt_run_control ) THEN |
---|
| 1059 | time_run_control = MOD( time_run_control, & |
---|
| 1060 | MAX( dt_run_control, dt_3d_old ) ) |
---|
| 1061 | ENDIF |
---|
| 1062 | ENDIF |
---|
| 1063 | |
---|
| 1064 | ! |
---|
[1402] | 1065 | !-- Output elapsed simulated time in form of a progress bar on stdout |
---|
| 1066 | IF ( myid == 0 ) CALL output_progress_bar |
---|
| 1067 | |
---|
[1] | 1068 | CALL cpu_log( log_point_s(10), 'timesteps', 'stop' ) |
---|
| 1069 | |
---|
[667] | 1070 | |
---|
[1] | 1071 | ENDDO ! time loop |
---|
| 1072 | |
---|
[1402] | 1073 | IF ( myid == 0 ) CALL finish_progress_bar |
---|
| 1074 | |
---|
[1] | 1075 | #if defined( __dvrp_graphics ) |
---|
| 1076 | CALL DVRP_LOG_EVENT( -2, current_timestep_number ) |
---|
| 1077 | #endif |
---|
| 1078 | |
---|
[1402] | 1079 | CALL location_message( 'finished time-stepping', .TRUE. ) |
---|
[1384] | 1080 | |
---|
[1] | 1081 | END SUBROUTINE time_integration |
---|