[1] | 1 | SUBROUTINE timestep |
---|
| 2 | |
---|
[1036] | 3 | !--------------------------------------------------------------------------------! |
---|
| 4 | ! This file is part of PALM. |
---|
| 5 | ! |
---|
| 6 | ! PALM is free software: you can redistribute it and/or modify it under the terms |
---|
| 7 | ! of the GNU General Public License as published by the Free Software Foundation, |
---|
| 8 | ! either version 3 of the License, or (at your option) any later version. |
---|
| 9 | ! |
---|
| 10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 13 | ! |
---|
| 14 | ! You should have received a copy of the GNU General Public License along with |
---|
| 15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 16 | ! |
---|
| 17 | ! Copyright 1997-2012 Leibniz University Hannover |
---|
| 18 | !--------------------------------------------------------------------------------! |
---|
| 19 | ! |
---|
[258] | 20 | ! Current revisions: |
---|
[866] | 21 | ! ------------------ |
---|
[1257] | 22 | ! openacc porting |
---|
| 23 | ! bugfix for calculation of advective timestep in case of vertically stretched |
---|
| 24 | ! grids |
---|
[316] | 25 | ! |
---|
[1] | 26 | ! Former revisions: |
---|
| 27 | ! ----------------- |
---|
[3] | 28 | ! $Id: timestep.f90 1257 2013-11-08 15:18:40Z raasch $ |
---|
[110] | 29 | ! |
---|
[1093] | 30 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
| 31 | ! unused variables removed |
---|
| 32 | ! |
---|
[1054] | 33 | ! 1053 2012-11-13 17:11:03Z hoffmann |
---|
| 34 | ! timestep is reduced in two-moment cloud scheme according to the maximum |
---|
| 35 | ! terminal velocity of rain drops |
---|
| 36 | ! |
---|
[1037] | 37 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 38 | ! code put under GPL (PALM 3.9) |
---|
| 39 | ! |
---|
[1002] | 40 | ! 1001 2012-09-13 14:08:46Z raasch |
---|
| 41 | ! all actions concerning leapfrog scheme removed |
---|
| 42 | ! |
---|
[979] | 43 | ! 978 2012-08-09 08:28:32Z fricke |
---|
| 44 | ! restriction of the outflow damping layer in the diffusion criterion removed |
---|
| 45 | ! |
---|
[867] | 46 | ! 866 2012-03-28 06:44:41Z raasch |
---|
| 47 | ! bugfix for timestep calculation in case of Galilei transformation, |
---|
| 48 | ! special treatment in case of mirror velocity boundary condition removed |
---|
| 49 | ! |
---|
[708] | 50 | ! 707 2011-03-29 11:39:40Z raasch |
---|
| 51 | ! bc_lr/ns replaced by bc_lr/ns_cyc |
---|
| 52 | ! |
---|
[668] | 53 | ! 667 2010-12-23 12:06:00Z suehring/gryschka |
---|
| 54 | ! Exchange of terminate_coupled between ocean and atmosphere via PE0 |
---|
| 55 | ! Minimum grid spacing dxyz2_min(k) is now calculated using dzw instead of dzu |
---|
| 56 | ! |
---|
[623] | 57 | ! 622 2010-12-10 08:08:13Z raasch |
---|
| 58 | ! optional barriers included in order to speed up collective operations |
---|
| 59 | ! |
---|
[392] | 60 | ! 343 2009-06-24 12:59:09Z maronga |
---|
| 61 | ! Additional timestep criterion in case of simulations with plant canopy |
---|
| 62 | ! Output of messages replaced by message handling routine. |
---|
| 63 | ! |
---|
[226] | 64 | ! 222 2009-01-12 16:04:16Z letzel |
---|
| 65 | ! Implementation of a MPI-1 Coupling: replaced myid with target_id |
---|
| 66 | ! Bugfix for nonparallel execution |
---|
| 67 | ! |
---|
[110] | 68 | ! 108 2007-08-24 15:10:38Z letzel |
---|
| 69 | ! modifications to terminate coupled runs |
---|
| 70 | ! |
---|
[3] | 71 | ! RCS Log replace by Id keyword, revision history cleaned up |
---|
| 72 | ! |
---|
[1] | 73 | ! Revision 1.21 2006/02/23 12:59:44 raasch |
---|
| 74 | ! nt_anz renamed current_timestep_number |
---|
| 75 | ! |
---|
| 76 | ! Revision 1.1 1997/08/11 06:26:19 raasch |
---|
| 77 | ! Initial revision |
---|
| 78 | ! |
---|
| 79 | ! |
---|
| 80 | ! Description: |
---|
| 81 | ! ------------ |
---|
| 82 | ! Compute the time step under consideration of the FCL and diffusion criterion. |
---|
| 83 | !------------------------------------------------------------------------------! |
---|
| 84 | |
---|
| 85 | USE arrays_3d |
---|
[1053] | 86 | USE cloud_parameters |
---|
[1] | 87 | USE control_parameters |
---|
| 88 | USE cpulog |
---|
| 89 | USE grid_variables |
---|
| 90 | USE indices |
---|
| 91 | USE interfaces |
---|
| 92 | USE pegrid |
---|
| 93 | USE statistics |
---|
| 94 | |
---|
| 95 | IMPLICIT NONE |
---|
| 96 | |
---|
[1257] | 97 | INTEGER :: i, j, k |
---|
[1] | 98 | |
---|
[1257] | 99 | REAL :: div, dt_diff, dt_diff_l, dt_plant_canopy, dt_plant_canopy_l, & |
---|
| 100 | dt_plant_canopy_u, dt_plant_canopy_v, dt_plant_canopy_w, & |
---|
| 101 | dt_u, dt_u_l, dt_v, dt_v_l, dt_w, dt_w_l, u_gtrans_l, u_max_l, & |
---|
| 102 | u_min_l, value, v_gtrans_l, v_max_l, v_min_l, w_max_l, w_min_l |
---|
[1] | 103 | |
---|
| 104 | REAL, DIMENSION(2) :: uv_gtrans, uv_gtrans_l |
---|
[1257] | 105 | REAL, DIMENSION(3) :: reduce, reduce_l |
---|
[1] | 106 | REAL, DIMENSION(nzb+1:nzt) :: dxyz2_min |
---|
| 107 | |
---|
[667] | 108 | |
---|
| 109 | |
---|
[1] | 110 | CALL cpu_log( log_point(12), 'calculate_timestep', 'start' ) |
---|
| 111 | |
---|
| 112 | ! |
---|
| 113 | !-- In case of Galilei-transform not using the geostrophic wind as translation |
---|
| 114 | !-- velocity, compute the volume-averaged horizontal velocity components, which |
---|
| 115 | !-- will then be subtracted from the horizontal wind for the time step and |
---|
| 116 | !-- horizontal advection routines. |
---|
| 117 | IF ( galilei_transformation .AND. .NOT. use_ug_for_galilei_tr ) THEN |
---|
| 118 | IF ( flow_statistics_called ) THEN |
---|
| 119 | ! |
---|
| 120 | !-- Horizontal averages already existent, just need to average them |
---|
| 121 | !-- vertically. |
---|
| 122 | u_gtrans = 0.0 |
---|
| 123 | v_gtrans = 0.0 |
---|
| 124 | DO k = nzb+1, nzt |
---|
| 125 | u_gtrans = u_gtrans + hom(k,1,1,0) |
---|
| 126 | v_gtrans = v_gtrans + hom(k,1,2,0) |
---|
| 127 | ENDDO |
---|
| 128 | u_gtrans = u_gtrans / REAL( nzt - nzb ) |
---|
| 129 | v_gtrans = v_gtrans / REAL( nzt - nzb ) |
---|
| 130 | ELSE |
---|
| 131 | ! |
---|
| 132 | !-- Averaging over the entire model domain. |
---|
[1257] | 133 | u_gtrans_l = 0.0 |
---|
| 134 | v_gtrans_l = 0.0 |
---|
| 135 | !$acc parallel present( u, v ) |
---|
[1] | 136 | DO i = nxl, nxr |
---|
| 137 | DO j = nys, nyn |
---|
| 138 | DO k = nzb+1, nzt |
---|
[1257] | 139 | u_gtrans_l = u_gtrans_l + u(k,j,i) |
---|
| 140 | v_gtrans_l = v_gtrans_l + v(k,j,i) |
---|
[1] | 141 | ENDDO |
---|
| 142 | ENDDO |
---|
| 143 | ENDDO |
---|
[1257] | 144 | !$acc end parallel |
---|
| 145 | uv_gtrans_l(1) = u_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) ) |
---|
| 146 | uv_gtrans_l(2) = v_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) ) |
---|
[1] | 147 | #if defined( __parallel ) |
---|
[622] | 148 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[1] | 149 | CALL MPI_ALLREDUCE( uv_gtrans_l, uv_gtrans, 2, MPI_REAL, MPI_SUM, & |
---|
| 150 | comm2d, ierr ) |
---|
| 151 | u_gtrans = uv_gtrans(1) / REAL( numprocs ) |
---|
| 152 | v_gtrans = uv_gtrans(2) / REAL( numprocs ) |
---|
| 153 | #else |
---|
| 154 | u_gtrans = uv_gtrans_l(1) |
---|
| 155 | v_gtrans = uv_gtrans_l(2) |
---|
| 156 | #endif |
---|
| 157 | ENDIF |
---|
| 158 | ENDIF |
---|
| 159 | |
---|
[866] | 160 | ! |
---|
[1257] | 161 | !-- Determine the maxima of the velocity components, including their |
---|
| 162 | !-- grid index positions. |
---|
| 163 | #if defined( __openacc ) |
---|
| 164 | IF ( dt_fixed ) THEN ! otherwise do it further below for better cache usage |
---|
| 165 | u_max_l = -999999.9 |
---|
| 166 | u_min_l = 999999.9 |
---|
| 167 | v_max_l = -999999.9 |
---|
| 168 | v_min_l = 999999.9 |
---|
| 169 | w_max_l = -999999.9 |
---|
| 170 | w_min_l = 999999.9 |
---|
| 171 | !$acc parallel present( u, v, w ) |
---|
| 172 | DO i = nxl, nxr |
---|
| 173 | DO j = nys, nyn |
---|
| 174 | DO k = nzb+1, nzt |
---|
| 175 | u_max_l = MAX( u_max_l, u(k,j,i) ) |
---|
| 176 | u_min_l = MIN( u_min_l, u(k,j,i) ) |
---|
| 177 | v_max_l = MAX( v_max_l, v(k,j,i) ) |
---|
| 178 | v_min_l = MIN( v_min_l, v(k,j,i) ) |
---|
| 179 | w_max_l = MAX( w_max_l, w(k,j,i) ) |
---|
| 180 | w_min_l = MIN( w_min_l, w(k,j,i) ) |
---|
| 181 | ENDDO |
---|
| 182 | ENDDO |
---|
| 183 | ENDDO |
---|
| 184 | !$acc end parallel |
---|
| 185 | #if defined( __parallel ) |
---|
| 186 | reduce_l(1) = u_max_l |
---|
| 187 | reduce_l(2) = v_max_l |
---|
| 188 | reduce_l(3) = w_max_l |
---|
| 189 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 190 | CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MAX, comm2d, ierr ) |
---|
| 191 | u_max = reduce(1) |
---|
| 192 | v_max = reduce(2) |
---|
| 193 | w_max = reduce(3) |
---|
| 194 | reduce_l(1) = u_min_l |
---|
| 195 | reduce_l(2) = v_min_l |
---|
| 196 | reduce_l(3) = w_min_l |
---|
| 197 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 198 | CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) |
---|
| 199 | IF ( ABS( reduce(1) ) > u_max ) u_max = reduce(1) |
---|
| 200 | IF ( ABS( reduce(2) ) > v_max ) v_max = reduce(2) |
---|
| 201 | IF ( ABS( reduce(3) ) > w_max ) w_max = reduce(3) |
---|
| 202 | #else |
---|
| 203 | IF ( ABS( u_min_l ) > u_max_l ) THEN |
---|
| 204 | u_max = u_min_l |
---|
| 205 | ELSE |
---|
| 206 | u_max = u_max_l |
---|
| 207 | ENDIF |
---|
| 208 | IF ( ABS( v_min_l ) > v_max_l ) THEN |
---|
| 209 | v_max = v_min_l |
---|
| 210 | ELSE |
---|
| 211 | v_max = v_max_l |
---|
| 212 | ENDIF |
---|
| 213 | IF ( ABS( w_min_l ) > w_max_l ) THEN |
---|
| 214 | w_max = w_min_l |
---|
| 215 | ELSE |
---|
| 216 | w_max = w_max_l |
---|
| 217 | ENDIF |
---|
| 218 | #endif |
---|
| 219 | ENDIF |
---|
| 220 | #else |
---|
[866] | 221 | CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, u, 'abs', 0.0, & |
---|
| 222 | u_max, u_max_ijk ) |
---|
| 223 | CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, v, 'abs', 0.0, & |
---|
| 224 | v_max, v_max_ijk ) |
---|
| 225 | CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, w, 'abs', 0.0, & |
---|
| 226 | w_max, w_max_ijk ) |
---|
[1257] | 227 | #endif |
---|
[866] | 228 | |
---|
[1257] | 229 | IF ( .NOT. dt_fixed ) THEN |
---|
| 230 | #if defined( __openacc ) |
---|
[866] | 231 | ! |
---|
[1257] | 232 | !-- Variable time step: |
---|
| 233 | !-- Calculate the maximum time step according to the CFL-criterion, |
---|
| 234 | !-- individually for each velocity component |
---|
| 235 | dt_u_l = 999999.9 |
---|
| 236 | dt_v_l = 999999.9 |
---|
| 237 | dt_w_l = 999999.9 |
---|
| 238 | u_max_l = -999999.9 |
---|
| 239 | u_min_l = 999999.9 |
---|
| 240 | v_max_l = -999999.9 |
---|
| 241 | v_min_l = 999999.9 |
---|
| 242 | w_max_l = -999999.9 |
---|
| 243 | w_min_l = 999999.9 |
---|
| 244 | !$acc parallel loop collapse(3) present( u, v, w ) |
---|
| 245 | DO i = nxl, nxr |
---|
| 246 | DO j = nys, nyn |
---|
| 247 | DO k = nzb+1, nzt |
---|
| 248 | dt_u_l = MIN( dt_u_l, ( dx / ( ABS( u(k,j,i) - u_gtrans ) + 1.0E-10 ) ) ) |
---|
| 249 | dt_v_l = MIN( dt_v_l, ( dy / ( ABS( v(k,j,i) - v_gtrans ) + 1.0E-10 ) ) ) |
---|
| 250 | dt_w_l = MIN( dt_w_l, ( dzu(k) / ( ABS( w(k,j,i) ) + 1.0E-10 ) ) ) |
---|
| 251 | u_max_l = MAX( u_max_l, u(k,j,i) ) |
---|
| 252 | u_min_l = MIN( u_min_l, u(k,j,i) ) |
---|
| 253 | v_max_l = MAX( v_max_l, v(k,j,i) ) |
---|
| 254 | v_min_l = MIN( v_min_l, v(k,j,i) ) |
---|
| 255 | w_max_l = MAX( w_max_l, w(k,j,i) ) |
---|
| 256 | w_min_l = MIN( w_min_l, w(k,j,i) ) |
---|
| 257 | ENDDO |
---|
| 258 | ENDDO |
---|
| 259 | ENDDO |
---|
| 260 | !$acc end parallel |
---|
[866] | 261 | |
---|
[1257] | 262 | #if defined( __parallel ) |
---|
| 263 | reduce_l(1) = dt_u_l |
---|
| 264 | reduce_l(2) = dt_v_l |
---|
| 265 | reduce_l(3) = dt_w_l |
---|
| 266 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 267 | CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) |
---|
| 268 | dt_u = reduce(1) |
---|
| 269 | dt_v = reduce(2) |
---|
| 270 | dt_w = reduce(3) |
---|
[866] | 271 | |
---|
[1257] | 272 | reduce_l(1) = u_max_l |
---|
| 273 | reduce_l(2) = v_max_l |
---|
| 274 | reduce_l(3) = w_max_l |
---|
| 275 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 276 | CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MAX, comm2d, ierr ) |
---|
| 277 | u_max = reduce(1) |
---|
| 278 | v_max = reduce(2) |
---|
| 279 | w_max = reduce(3) |
---|
| 280 | reduce_l(1) = u_min_l |
---|
| 281 | reduce_l(2) = v_min_l |
---|
| 282 | reduce_l(3) = w_min_l |
---|
| 283 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 284 | CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) |
---|
| 285 | IF ( ABS( reduce(1) ) > u_max ) u_max = reduce(1) |
---|
| 286 | IF ( ABS( reduce(2) ) > v_max ) v_max = reduce(2) |
---|
| 287 | IF ( ABS( reduce(3) ) > w_max ) w_max = reduce(3) |
---|
| 288 | #else |
---|
| 289 | dt_u = dt_u_l |
---|
| 290 | dt_v = dt_v_l |
---|
| 291 | dt_w = dt_w_l |
---|
| 292 | |
---|
| 293 | IF ( ABS( u_min_l ) > u_max_l ) THEN |
---|
| 294 | u_max = u_min_l |
---|
| 295 | ELSE |
---|
| 296 | u_max = u_max_l |
---|
| 297 | ENDIF |
---|
| 298 | IF ( ABS( v_min_l ) > v_max_l ) THEN |
---|
| 299 | v_max = v_min_l |
---|
| 300 | ELSE |
---|
| 301 | v_max = v_max_l |
---|
| 302 | ENDIF |
---|
| 303 | IF ( ABS( w_min_l ) > w_max_l ) THEN |
---|
| 304 | w_max = w_min_l |
---|
| 305 | ELSE |
---|
| 306 | w_max = w_max_l |
---|
| 307 | ENDIF |
---|
| 308 | #endif |
---|
| 309 | |
---|
| 310 | #else |
---|
[1] | 311 | ! |
---|
| 312 | !-- Variable time step: |
---|
[1257] | 313 | !-- Calculate the maximum time step according to the CFL-criterion, |
---|
| 314 | !-- individually for each velocity component |
---|
| 315 | dt_u_l = 999999.9 |
---|
| 316 | dt_v_l = 999999.9 |
---|
| 317 | dt_w_l = 999999.9 |
---|
| 318 | DO i = nxl, nxr |
---|
| 319 | DO j = nys, nyn |
---|
| 320 | DO k = nzb+1, nzt |
---|
| 321 | dt_u_l = MIN( dt_u_l, ( dx / ( ABS( u(k,j,i) - u_gtrans ) + 1.0E-10 ) ) ) |
---|
| 322 | dt_v_l = MIN( dt_v_l, ( dy / ( ABS( v(k,j,i) - v_gtrans ) + 1.0E-10 ) ) ) |
---|
| 323 | dt_w_l = MIN( dt_w_l, ( dzu(k) / ( ABS( w(k,j,i) ) + 1.0E-10 ) ) ) |
---|
| 324 | ENDDO |
---|
| 325 | ENDDO |
---|
| 326 | ENDDO |
---|
[1] | 327 | |
---|
[1257] | 328 | #if defined( __parallel ) |
---|
| 329 | reduce_l(1) = dt_u_l |
---|
| 330 | reduce_l(2) = dt_v_l |
---|
| 331 | reduce_l(3) = dt_w_l |
---|
| 332 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 333 | CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) |
---|
| 334 | dt_u = reduce(1) |
---|
| 335 | dt_v = reduce(2) |
---|
| 336 | dt_w = reduce(3) |
---|
| 337 | #else |
---|
| 338 | dt_u = dt_u_l |
---|
| 339 | dt_v = dt_v_l |
---|
| 340 | dt_w = dt_w_l |
---|
| 341 | #endif |
---|
| 342 | |
---|
| 343 | #endif |
---|
| 344 | |
---|
[1] | 345 | ! |
---|
| 346 | !-- Compute time step according to the diffusion criterion. |
---|
| 347 | !-- First calculate minimum grid spacing which only depends on index k |
---|
| 348 | !-- Note: also at k=nzb+1 a full grid length is being assumed, although |
---|
| 349 | !-- in the Prandtl-layer friction term only dz/2 is used. |
---|
| 350 | !-- Experience from the old model seems to justify this. |
---|
| 351 | dt_diff_l = 999999.0 |
---|
| 352 | |
---|
| 353 | DO k = nzb+1, nzt |
---|
[667] | 354 | dxyz2_min(k) = MIN( dx2, dy2, dzw(k)*dzw(k) ) * 0.125 |
---|
[1] | 355 | ENDDO |
---|
| 356 | |
---|
| 357 | !$OMP PARALLEL private(i,j,k,value) reduction(MIN: dt_diff_l) |
---|
| 358 | !$OMP DO |
---|
[1257] | 359 | !$acc parallel loop collapse(3) present( kh, km ) |
---|
[1] | 360 | DO i = nxl, nxr |
---|
| 361 | DO j = nys, nyn |
---|
| 362 | DO k = nzb+1, nzt |
---|
[1257] | 363 | dt_diff_l = MIN( dt_diff_l, dxyz2_min(k) / & |
---|
| 364 | ( MAX( kh(k,j,i), km(k,j,i) ) + 1E-20 ) ) |
---|
[1] | 365 | ENDDO |
---|
| 366 | ENDDO |
---|
| 367 | ENDDO |
---|
[1257] | 368 | !$acc end parallel |
---|
[1] | 369 | !$OMP END PARALLEL |
---|
| 370 | #if defined( __parallel ) |
---|
[622] | 371 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[1] | 372 | CALL MPI_ALLREDUCE( dt_diff_l, dt_diff, 1, MPI_REAL, MPI_MIN, comm2d, & |
---|
| 373 | ierr ) |
---|
| 374 | #else |
---|
| 375 | dt_diff = dt_diff_l |
---|
| 376 | #endif |
---|
| 377 | |
---|
| 378 | ! |
---|
[316] | 379 | !-- Additional timestep criterion with plant canopies: |
---|
| 380 | !-- it is not allowed to extract more than the available momentum |
---|
| 381 | IF ( plant_canopy ) THEN |
---|
[318] | 382 | |
---|
| 383 | dt_plant_canopy_l = 0.0 |
---|
| 384 | DO i = nxl, nxr |
---|
| 385 | DO j = nys, nyn |
---|
| 386 | DO k = nzb+1, nzt |
---|
| 387 | dt_plant_canopy_u = cdc(k,j,i) * lad_u(k,j,i) * & |
---|
| 388 | SQRT( u(k,j,i)**2 + & |
---|
| 389 | ( ( v(k,j,i-1) + & |
---|
| 390 | v(k,j,i) + & |
---|
| 391 | v(k,j+1,i) + & |
---|
| 392 | v(k,j+1,i-1) ) & |
---|
| 393 | / 4.0 )**2 + & |
---|
| 394 | ( ( w(k-1,j,i-1) + & |
---|
| 395 | w(k-1,j,i) + & |
---|
| 396 | w(k,j,i-1) + & |
---|
| 397 | w(k,j,i) ) & |
---|
| 398 | / 4.0 )**2 ) |
---|
| 399 | IF ( dt_plant_canopy_u > dt_plant_canopy_l ) THEN |
---|
| 400 | dt_plant_canopy_l = dt_plant_canopy_u |
---|
| 401 | ENDIF |
---|
| 402 | dt_plant_canopy_v = cdc(k,j,i) * lad_v(k,j,i) * & |
---|
| 403 | SQRT( ( ( u(k,j-1,i) + & |
---|
| 404 | u(k,j-1,i+1) + & |
---|
| 405 | u(k,j,i) + & |
---|
| 406 | u(k,j,i+1) ) & |
---|
| 407 | / 4.0 )**2 + & |
---|
| 408 | v(k,j,i)**2 + & |
---|
| 409 | ( ( w(k-1,j-1,i) + & |
---|
| 410 | w(k-1,j,i) + & |
---|
| 411 | w(k,j-1,i) + & |
---|
| 412 | w(k,j,i) ) & |
---|
| 413 | / 4.0 )**2 ) |
---|
| 414 | IF ( dt_plant_canopy_v > dt_plant_canopy_l ) THEN |
---|
| 415 | dt_plant_canopy_l = dt_plant_canopy_v |
---|
| 416 | ENDIF |
---|
| 417 | dt_plant_canopy_w = cdc(k,j,i) * lad_w(k,j,i) * & |
---|
| 418 | SQRT( ( ( u(k,j,i) + & |
---|
| 419 | u(k,j,i+1) + & |
---|
| 420 | u(k+1,j,i) + & |
---|
| 421 | u(k+1,j,i+1) ) & |
---|
| 422 | / 4.0 )**2 + & |
---|
| 423 | ( ( v(k,j,i) + & |
---|
| 424 | v(k,j+1,i) + & |
---|
| 425 | v(k+1,j,i) + & |
---|
| 426 | v(k+1,j+1,i) ) & |
---|
| 427 | / 4.0 )**2 + & |
---|
| 428 | w(k,j,i)**2 ) |
---|
| 429 | IF ( dt_plant_canopy_w > dt_plant_canopy_l ) THEN |
---|
| 430 | dt_plant_canopy_l = dt_plant_canopy_w |
---|
| 431 | ENDIF |
---|
| 432 | ENDDO |
---|
| 433 | ENDDO |
---|
| 434 | ENDDO |
---|
| 435 | |
---|
| 436 | IF ( dt_plant_canopy_l > 0.0 ) THEN |
---|
[320] | 437 | ! |
---|
| 438 | !-- Invert dt_plant_canopy_l and apply a security timestep factor 0.1 |
---|
[318] | 439 | dt_plant_canopy_l = 0.1 / dt_plant_canopy_l |
---|
[320] | 440 | ELSE |
---|
| 441 | ! |
---|
| 442 | !-- In case of inhomogeneous plant canopy, some processors may have no |
---|
| 443 | !-- canopy at all. Then use dt_max as dummy instead. |
---|
| 444 | dt_plant_canopy_l = dt_max |
---|
[318] | 445 | ENDIF |
---|
[320] | 446 | |
---|
[316] | 447 | ! |
---|
[318] | 448 | !-- Determine the global minumum |
---|
| 449 | #if defined( __parallel ) |
---|
[622] | 450 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[866] | 451 | CALL MPI_ALLREDUCE( dt_plant_canopy_l, dt_plant_canopy, 1, MPI_REAL, & |
---|
[318] | 452 | MPI_MIN, comm2d, ierr ) |
---|
| 453 | #else |
---|
| 454 | dt_plant_canopy = dt_plant_canopy_l |
---|
| 455 | #endif |
---|
[316] | 456 | |
---|
| 457 | ELSE |
---|
| 458 | ! |
---|
| 459 | !-- Use dt_diff as dummy value to avoid extra IF branches further below |
---|
| 460 | dt_plant_canopy = dt_diff |
---|
| 461 | |
---|
| 462 | ENDIF |
---|
| 463 | |
---|
| 464 | ! |
---|
| 465 | !-- The time step is the minimum of the 3-4 components and the diffusion time |
---|
[1001] | 466 | !-- step minus a reduction (cfl_factor) to be on the safe side. |
---|
[1] | 467 | !-- The time step must not exceed the maximum allowed value. |
---|
[1053] | 468 | dt_3d = cfl_factor * MIN( dt_diff, dt_plant_canopy, dt_u, dt_v, dt_w, & |
---|
| 469 | dt_precipitation ) |
---|
[1] | 470 | dt_3d = MIN( dt_3d, dt_max ) |
---|
| 471 | |
---|
| 472 | ! |
---|
| 473 | !-- Remember the restricting time step criterion for later output. |
---|
[316] | 474 | IF ( MIN( dt_u, dt_v, dt_w ) < MIN( dt_diff, dt_plant_canopy ) ) THEN |
---|
[1] | 475 | timestep_reason = 'A' |
---|
[316] | 476 | ELSEIF ( dt_plant_canopy < dt_diff ) THEN |
---|
| 477 | timestep_reason = 'C' |
---|
[1] | 478 | ELSE |
---|
| 479 | timestep_reason = 'D' |
---|
| 480 | ENDIF |
---|
| 481 | |
---|
| 482 | ! |
---|
| 483 | !-- Set flag if the time step becomes too small. |
---|
| 484 | IF ( dt_3d < ( 0.00001 * dt_max ) ) THEN |
---|
| 485 | stop_dt = .TRUE. |
---|
[108] | 486 | |
---|
[320] | 487 | WRITE( message_string, * ) 'Time step has reached minimum limit.', & |
---|
| 488 | '&dt = ', dt_3d, ' s Simulation is terminated.', & |
---|
| 489 | '&old_dt = ', old_dt, ' s', & |
---|
| 490 | '&dt_u = ', dt_u, ' s', & |
---|
| 491 | '&dt_v = ', dt_v, ' s', & |
---|
| 492 | '&dt_w = ', dt_w, ' s', & |
---|
| 493 | '&dt_diff = ', dt_diff, ' s', & |
---|
| 494 | '&dt_plant_canopy = ', dt_plant_canopy, ' s', & |
---|
[1257] | 495 | '&u_max = ', u_max, ' m/s k=', u_max_ijk(1), & |
---|
[320] | 496 | ' j=', u_max_ijk(2), ' i=', u_max_ijk(3), & |
---|
[1257] | 497 | '&v_max = ', v_max, ' m/s k=', v_max_ijk(1), & |
---|
[320] | 498 | ' j=', v_max_ijk(2), ' i=', v_max_ijk(3), & |
---|
[1257] | 499 | '&w_max = ', w_max, ' m/s k=', w_max_ijk(1), & |
---|
[320] | 500 | ' j=', w_max_ijk(2), ' i=', w_max_ijk(3) |
---|
[258] | 501 | CALL message( 'timestep', 'PA0312', 0, 1, 0, 6, 0 ) |
---|
[108] | 502 | ! |
---|
| 503 | !-- In case of coupled runs inform the remote model of the termination |
---|
| 504 | !-- and its reason, provided the remote model has not already been |
---|
| 505 | !-- informed of another termination reason (terminate_coupled > 0) before. |
---|
[222] | 506 | #if defined( __parallel ) |
---|
[108] | 507 | IF ( coupling_mode /= 'uncoupled' .AND. terminate_coupled == 0 ) THEN |
---|
| 508 | terminate_coupled = 2 |
---|
[667] | 509 | IF ( myid == 0 ) THEN |
---|
| 510 | CALL MPI_SENDRECV( & |
---|
| 511 | terminate_coupled, 1, MPI_INTEGER, target_id, 0, & |
---|
| 512 | terminate_coupled_remote, 1, MPI_INTEGER, target_id, 0, & |
---|
| 513 | comm_inter, status, ierr ) |
---|
| 514 | ENDIF |
---|
| 515 | CALL MPI_BCAST( terminate_coupled_remote, 1, MPI_INTEGER, 0, comm2d, ierr) |
---|
[108] | 516 | ENDIF |
---|
[222] | 517 | #endif |
---|
[1] | 518 | ENDIF |
---|
| 519 | |
---|
| 520 | ! |
---|
[1001] | 521 | !-- Ensure a smooth value (two significant digits) of the timestep. |
---|
| 522 | div = 1000.0 |
---|
| 523 | DO WHILE ( dt_3d < div ) |
---|
| 524 | div = div / 10.0 |
---|
| 525 | ENDDO |
---|
| 526 | dt_3d = NINT( dt_3d * 100.0 / div ) * div / 100.0 |
---|
[1] | 527 | |
---|
| 528 | ! |
---|
[1001] | 529 | !-- Adjust the time step |
---|
| 530 | old_dt = dt_3d |
---|
[1] | 531 | |
---|
[1001] | 532 | ENDIF |
---|
[1] | 533 | |
---|
| 534 | CALL cpu_log( log_point(12), 'calculate_timestep', 'stop' ) |
---|
| 535 | |
---|
| 536 | END SUBROUTINE timestep |
---|