SUBROUTINE timestep !------------------------------------------------------------------------------! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: timestep.f90 392 2009-09-24 10:39:14Z maronga $ ! ! 343 2009-06-24 12:59:09Z maronga ! Additional timestep criterion in case of simulations with plant canopy ! Output of messages replaced by message handling routine. ! ! 222 2009-01-12 16:04:16Z letzel ! Implementation of a MPI-1 Coupling: replaced myid with target_id ! Bugfix for nonparallel execution ! ! 108 2007-08-24 15:10:38Z letzel ! modifications to terminate coupled runs ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.21 2006/02/23 12:59:44 raasch ! nt_anz renamed current_timestep_number ! ! Revision 1.1 1997/08/11 06:26:19 raasch ! Initial revision ! ! ! Description: ! ------------ ! Compute the time step under consideration of the FCL and diffusion criterion. !------------------------------------------------------------------------------! USE arrays_3d USE control_parameters USE cpulog USE grid_variables USE indices USE interfaces USE pegrid USE statistics IMPLICIT NONE INTEGER :: i, j, k REAL :: div, dt_diff, dt_diff_l, dt_plant_canopy, & dt_plant_canopy_l, & dt_plant_canopy_u, dt_plant_canopy_v, dt_plant_canopy_w, & dt_u, dt_v, dt_w, lad_max, percent_change, & u_gtrans_l, vabs_max, value, v_gtrans_l REAL, DIMENSION(2) :: uv_gtrans, uv_gtrans_l REAL, DIMENSION(nzb+1:nzt) :: dxyz2_min CALL cpu_log( log_point(12), 'calculate_timestep', 'start' ) ! !-- Determine the maxima of the velocity components. CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, u, 'abs', & u_max, u_max_ijk ) CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, v, 'abs', & v_max, v_max_ijk ) CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, w, 'abs', & w_max, w_max_ijk ) ! !-- In case maxima of the horizontal velocity components have been found at the !-- bottom boundary (k=nzb), the corresponding maximum at level k=1 is chosen !-- if the Dirichlet-boundary condition ('mirror') has been set. This is !-- necessary, because otherwise in case of Galilei-transform a far too large !-- velocity (having the respective opposite sign) would be used for the time !-- step determination (almost double the mean flow velocity). IF ( ibc_uv_b == 0 ) THEN IF ( u_max_ijk(1) == nzb ) THEN u_max = -u_max u_max_ijk(1) = nzb + 1 ENDIF IF ( v_max_ijk(1) == nzb ) THEN v_max = -v_max v_max_ijk(1) = nzb + 1 ENDIF ENDIF ! !-- In case of Galilei-transform not using the geostrophic wind as translation !-- velocity, compute the volume-averaged horizontal velocity components, which !-- will then be subtracted from the horizontal wind for the time step and !-- horizontal advection routines. IF ( galilei_transformation .AND. .NOT. use_ug_for_galilei_tr ) THEN IF ( flow_statistics_called ) THEN ! !-- Horizontal averages already existent, just need to average them !-- vertically. u_gtrans = 0.0 v_gtrans = 0.0 DO k = nzb+1, nzt u_gtrans = u_gtrans + hom(k,1,1,0) v_gtrans = v_gtrans + hom(k,1,2,0) ENDDO u_gtrans = u_gtrans / REAL( nzt - nzb ) v_gtrans = v_gtrans / REAL( nzt - nzb ) ELSE ! !-- Averaging over the entire model domain. uv_gtrans_l = 0.0 DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt uv_gtrans_l(1) = uv_gtrans_l(1) + u(k,j,i) uv_gtrans_l(2) = uv_gtrans_l(2) + v(k,j,i) ENDDO ENDDO ENDDO uv_gtrans_l = uv_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) ) #if defined( __parallel ) CALL MPI_ALLREDUCE( uv_gtrans_l, uv_gtrans, 2, MPI_REAL, MPI_SUM, & comm2d, ierr ) u_gtrans = uv_gtrans(1) / REAL( numprocs ) v_gtrans = uv_gtrans(2) / REAL( numprocs ) #else u_gtrans = uv_gtrans_l(1) v_gtrans = uv_gtrans_l(2) #endif ENDIF ENDIF IF ( .NOT. dt_fixed ) THEN ! !-- Variable time step: ! !-- For each component, compute the maximum time step according to the !-- FCL-criterion. dt_u = dx / ( ABS( u_max - u_gtrans ) + 1.0E-10 ) dt_v = dy / ( ABS( v_max - v_gtrans ) + 1.0E-10 ) dt_w = dzu(MAX( 1, w_max_ijk(1) )) / ( ABS( w_max ) + 1.0E-10 ) ! !-- Compute time step according to the diffusion criterion. !-- First calculate minimum grid spacing which only depends on index k !-- Note: also at k=nzb+1 a full grid length is being assumed, although !-- in the Prandtl-layer friction term only dz/2 is used. !-- Experience from the old model seems to justify this. dt_diff_l = 999999.0 DO k = nzb+1, nzt dxyz2_min(k) = MIN( dx2, dy2, dzu(k)*dzu(k) ) * 0.125 ENDDO !$OMP PARALLEL private(i,j,k,value) reduction(MIN: dt_diff_l) !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt value = dxyz2_min(k) / ( MAX( kh(k,j,i), km(k,j,i) ) + 1E-20 ) dt_diff_l = MIN( value, dt_diff_l ) ENDDO ENDDO ENDDO !$OMP END PARALLEL #if defined( __parallel ) CALL MPI_ALLREDUCE( dt_diff_l, dt_diff, 1, MPI_REAL, MPI_MIN, comm2d, & ierr ) #else dt_diff = dt_diff_l #endif ! !-- In case of non-cyclic lateral boundaries, the diffusion time step !-- may be further restricted by the lateral damping layer (damping only !-- along x and y) IF ( bc_lr /= 'cyclic' ) THEN dt_diff = MIN( dt_diff, 0.125 * dx2 / ( km_damp_max + 1E-20 ) ) ELSEIF ( bc_ns /= 'cyclic' ) THEN dt_diff = MIN( dt_diff, 0.125 * dy2 / ( km_damp_max + 1E-20 ) ) ENDIF ! !-- Additional timestep criterion with plant canopies: !-- it is not allowed to extract more than the available momentum IF ( plant_canopy ) THEN dt_plant_canopy_l = 0.0 DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt dt_plant_canopy_u = cdc(k,j,i) * lad_u(k,j,i) * & SQRT( u(k,j,i)**2 + & ( ( v(k,j,i-1) + & v(k,j,i) + & v(k,j+1,i) + & v(k,j+1,i-1) ) & / 4.0 )**2 + & ( ( w(k-1,j,i-1) + & w(k-1,j,i) + & w(k,j,i-1) + & w(k,j,i) ) & / 4.0 )**2 ) IF ( dt_plant_canopy_u > dt_plant_canopy_l ) THEN dt_plant_canopy_l = dt_plant_canopy_u ENDIF dt_plant_canopy_v = cdc(k,j,i) * lad_v(k,j,i) * & SQRT( ( ( u(k,j-1,i) + & u(k,j-1,i+1) + & u(k,j,i) + & u(k,j,i+1) ) & / 4.0 )**2 + & v(k,j,i)**2 + & ( ( w(k-1,j-1,i) + & w(k-1,j,i) + & w(k,j-1,i) + & w(k,j,i) ) & / 4.0 )**2 ) IF ( dt_plant_canopy_v > dt_plant_canopy_l ) THEN dt_plant_canopy_l = dt_plant_canopy_v ENDIF dt_plant_canopy_w = cdc(k,j,i) * lad_w(k,j,i) * & SQRT( ( ( u(k,j,i) + & u(k,j,i+1) + & u(k+1,j,i) + & u(k+1,j,i+1) ) & / 4.0 )**2 + & ( ( v(k,j,i) + & v(k,j+1,i) + & v(k+1,j,i) + & v(k+1,j+1,i) ) & / 4.0 )**2 + & w(k,j,i)**2 ) IF ( dt_plant_canopy_w > dt_plant_canopy_l ) THEN dt_plant_canopy_l = dt_plant_canopy_w ENDIF ENDDO ENDDO ENDDO IF ( dt_plant_canopy_l > 0.0 ) THEN ! !-- Invert dt_plant_canopy_l and apply a security timestep factor 0.1 dt_plant_canopy_l = 0.1 / dt_plant_canopy_l ELSE ! !-- In case of inhomogeneous plant canopy, some processors may have no !-- canopy at all. Then use dt_max as dummy instead. dt_plant_canopy_l = dt_max ENDIF ! !-- Determine the global minumum #if defined( __parallel ) CALL MPI_ALLREDUCE( dt_plant_canopy_l, dt_plant_canopy, 1, MPI_REAL, & MPI_MIN, comm2d, ierr ) #else dt_plant_canopy = dt_plant_canopy_l #endif ELSE ! !-- Use dt_diff as dummy value to avoid extra IF branches further below dt_plant_canopy = dt_diff ENDIF ! !-- The time step is the minimum of the 3-4 components and the diffusion time !-- step minus a reduction to be on the safe side. Factor 0.5 is necessary !-- since the leap-frog scheme always progresses by 2 * delta t. !-- The user has to set the cfl_factor small enough to ensure that the !-- divergences do not become too large. !-- The time step must not exceed the maximum allowed value. IF ( timestep_scheme(1:5) == 'runge' ) THEN dt_3d = cfl_factor * MIN( dt_diff, dt_plant_canopy, dt_u, dt_v, dt_w ) ELSE dt_3d = cfl_factor * 0.5 * & MIN( dt_diff, dt_plant_canopy, dt_u, dt_v, dt_w ) ENDIF dt_3d = MIN( dt_3d, dt_max ) ! !-- Remember the restricting time step criterion for later output. IF ( MIN( dt_u, dt_v, dt_w ) < MIN( dt_diff, dt_plant_canopy ) ) THEN timestep_reason = 'A' ELSEIF ( dt_plant_canopy < dt_diff ) THEN timestep_reason = 'C' ELSE timestep_reason = 'D' ENDIF ! !-- Set flag if the time step becomes too small. IF ( dt_3d < ( 0.00001 * dt_max ) ) THEN stop_dt = .TRUE. WRITE( message_string, * ) 'Time step has reached minimum limit.', & '&dt = ', dt_3d, ' s Simulation is terminated.', & '&old_dt = ', old_dt, ' s', & '&dt_u = ', dt_u, ' s', & '&dt_v = ', dt_v, ' s', & '&dt_w = ', dt_w, ' s', & '&dt_diff = ', dt_diff, ' s', & '&dt_plant_canopy = ', dt_plant_canopy, ' s', & '&u_max = ', u_max, ' m/s k=', u_max_ijk(1), & ' j=', u_max_ijk(2), ' i=', u_max_ijk(3), & '&v_max = ', v_max, ' m/s k=', v_max_ijk(1), & ' j=', v_max_ijk(2), ' i=', v_max_ijk(3), & '&w_max = ', w_max, ' m/s k=', w_max_ijk(1), & ' j=', w_max_ijk(2), ' i=', w_max_ijk(3) CALL message( 'timestep', 'PA0312', 0, 1, 0, 6, 0 ) ! !-- In case of coupled runs inform the remote model of the termination !-- and its reason, provided the remote model has not already been !-- informed of another termination reason (terminate_coupled > 0) before. #if defined( __parallel ) IF ( coupling_mode /= 'uncoupled' .AND. terminate_coupled == 0 ) THEN terminate_coupled = 2 CALL MPI_SENDRECV( & terminate_coupled, 1, MPI_INTEGER, target_id, 0, & terminate_coupled_remote, 1, MPI_INTEGER, target_id, 0, & comm_inter, status, ierr ) ENDIF #endif ENDIF ! !-- Ensure a smooth value (two significant digits) of the timestep. For !-- other schemes than Runge-Kutta, the following restrictions appear: !-- The current timestep is only then changed, if the change relative to !-- its previous value exceeds +5 % or -2 %. In case of a timestep !-- reduction, at least 30 iterations have to be performed before a timestep !-- enlargement is permitted again. percent_change = dt_3d / old_dt - 1.0 IF ( percent_change > 0.05 .OR. percent_change < -0.02 .OR. & timestep_scheme(1:5) == 'runge' ) THEN ! !-- Time step enlargement by no more than 2 %. IF ( percent_change > 0.0 .AND. simulated_time /= 0.0 .AND. & timestep_scheme(1:5) /= 'runge' ) THEN dt_3d = 1.02 * old_dt ENDIF ! !-- A relatively smooth value of the time step is ensured by taking !-- only the first two significant digits. div = 1000.0 DO WHILE ( dt_3d < div ) div = div / 10.0 ENDDO dt_3d = NINT( dt_3d * 100.0 / div ) * div / 100.0 ! !-- Now the time step can be adjusted. IF ( percent_change < 0.0 .OR. timestep_scheme(1:5) == 'runge' ) & THEN ! !-- Time step reduction. old_dt = dt_3d dt_changed = .TRUE. ELSE ! !-- For other timestep schemes , the time step is only enlarged !-- after at least 30 iterations since the previous time step !-- change or, of course, after model initialization. IF ( current_timestep_number >= last_dt_change + 30 .OR. & simulated_time == 0.0 ) THEN old_dt = dt_3d dt_changed = .TRUE. ELSE dt_3d = old_dt dt_changed = .FALSE. ENDIF ENDIF ELSE ! !-- No time step change since the difference is too small. dt_3d = old_dt dt_changed = .FALSE. ENDIF IF ( dt_changed ) last_dt_change = current_timestep_number ENDIF CALL cpu_log( log_point(12), 'calculate_timestep', 'stop' ) END SUBROUTINE timestep