SUBROUTINE timestep !--------------------------------------------------------------------------------! ! This file is part of PALM. ! ! PALM is free software: you can redistribute it and/or modify it under the terms ! of the GNU General Public License as published by the Free Software Foundation, ! either version 3 of the License, or (at your option) any later version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 1997-2012 Leibniz University Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ------------------ ! ! ! Former revisions: ! ----------------- ! $Id: timestep.f90 1258 2013-11-08 16:09:09Z witha $ ! ! 1257 2013-11-08 15:18:40Z raasch ! openacc porting ! bugfix for calculation of advective timestep in case of vertically stretched ! grids ! ! 1092 2013-02-02 11:24:22Z raasch ! unused variables removed ! ! 1053 2012-11-13 17:11:03Z hoffmann ! timestep is reduced in two-moment cloud scheme according to the maximum ! terminal velocity of rain drops ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1001 2012-09-13 14:08:46Z raasch ! all actions concerning leapfrog scheme removed ! ! 978 2012-08-09 08:28:32Z fricke ! restriction of the outflow damping layer in the diffusion criterion removed ! ! 866 2012-03-28 06:44:41Z raasch ! bugfix for timestep calculation in case of Galilei transformation, ! special treatment in case of mirror velocity boundary condition removed ! ! 707 2011-03-29 11:39:40Z raasch ! bc_lr/ns replaced by bc_lr/ns_cyc ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! Exchange of terminate_coupled between ocean and atmosphere via PE0 ! Minimum grid spacing dxyz2_min(k) is now calculated using dzw instead of dzu ! ! 622 2010-12-10 08:08:13Z raasch ! optional barriers included in order to speed up collective operations ! ! 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 cloud_parameters 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_u_l, dt_v, dt_v_l, dt_w, dt_w_l, u_gtrans_l, u_max_l, & u_min_l, value, v_gtrans_l, v_max_l, v_min_l, w_max_l, w_min_l REAL, DIMENSION(2) :: uv_gtrans, uv_gtrans_l REAL, DIMENSION(3) :: reduce, reduce_l REAL, DIMENSION(nzb+1:nzt) :: dxyz2_min CALL cpu_log( log_point(12), 'calculate_timestep', 'start' ) ! !-- 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. u_gtrans_l = 0.0 v_gtrans_l = 0.0 !$acc parallel present( u, v ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt u_gtrans_l = u_gtrans_l + u(k,j,i) v_gtrans_l = v_gtrans_l + v(k,j,i) ENDDO ENDDO ENDDO !$acc end parallel uv_gtrans_l(1) = u_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) ) uv_gtrans_l(2) = v_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) ) #if defined( __parallel ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) 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 ! !-- Determine the maxima of the velocity components, including their !-- grid index positions. #if defined( __openacc ) IF ( dt_fixed ) THEN ! otherwise do it further below for better cache usage u_max_l = -999999.9 u_min_l = 999999.9 v_max_l = -999999.9 v_min_l = 999999.9 w_max_l = -999999.9 w_min_l = 999999.9 !$acc parallel present( u, v, w ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt u_max_l = MAX( u_max_l, u(k,j,i) ) u_min_l = MIN( u_min_l, u(k,j,i) ) v_max_l = MAX( v_max_l, v(k,j,i) ) v_min_l = MIN( v_min_l, v(k,j,i) ) w_max_l = MAX( w_max_l, w(k,j,i) ) w_min_l = MIN( w_min_l, w(k,j,i) ) ENDDO ENDDO ENDDO !$acc end parallel #if defined( __parallel ) reduce_l(1) = u_max_l reduce_l(2) = v_max_l reduce_l(3) = w_max_l IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MAX, comm2d, ierr ) u_max = reduce(1) v_max = reduce(2) w_max = reduce(3) reduce_l(1) = u_min_l reduce_l(2) = v_min_l reduce_l(3) = w_min_l IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) IF ( ABS( reduce(1) ) > u_max ) u_max = reduce(1) IF ( ABS( reduce(2) ) > v_max ) v_max = reduce(2) IF ( ABS( reduce(3) ) > w_max ) w_max = reduce(3) #else IF ( ABS( u_min_l ) > u_max_l ) THEN u_max = u_min_l ELSE u_max = u_max_l ENDIF IF ( ABS( v_min_l ) > v_max_l ) THEN v_max = v_min_l ELSE v_max = v_max_l ENDIF IF ( ABS( w_min_l ) > w_max_l ) THEN w_max = w_min_l ELSE w_max = w_max_l ENDIF #endif ENDIF #else CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, u, 'abs', 0.0, & u_max, u_max_ijk ) CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, v, 'abs', 0.0, & v_max, v_max_ijk ) CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, w, 'abs', 0.0, & w_max, w_max_ijk ) #endif IF ( .NOT. dt_fixed ) THEN #if defined( __openacc ) ! !-- Variable time step: !-- Calculate the maximum time step according to the CFL-criterion, !-- individually for each velocity component dt_u_l = 999999.9 dt_v_l = 999999.9 dt_w_l = 999999.9 u_max_l = -999999.9 u_min_l = 999999.9 v_max_l = -999999.9 v_min_l = 999999.9 w_max_l = -999999.9 w_min_l = 999999.9 !$acc parallel loop collapse(3) present( u, v, w ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt dt_u_l = MIN( dt_u_l, ( dx / ( ABS( u(k,j,i) - u_gtrans ) + 1.0E-10 ) ) ) dt_v_l = MIN( dt_v_l, ( dy / ( ABS( v(k,j,i) - v_gtrans ) + 1.0E-10 ) ) ) dt_w_l = MIN( dt_w_l, ( dzu(k) / ( ABS( w(k,j,i) ) + 1.0E-10 ) ) ) u_max_l = MAX( u_max_l, u(k,j,i) ) u_min_l = MIN( u_min_l, u(k,j,i) ) v_max_l = MAX( v_max_l, v(k,j,i) ) v_min_l = MIN( v_min_l, v(k,j,i) ) w_max_l = MAX( w_max_l, w(k,j,i) ) w_min_l = MIN( w_min_l, w(k,j,i) ) ENDDO ENDDO ENDDO !$acc end parallel #if defined( __parallel ) reduce_l(1) = dt_u_l reduce_l(2) = dt_v_l reduce_l(3) = dt_w_l IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) dt_u = reduce(1) dt_v = reduce(2) dt_w = reduce(3) reduce_l(1) = u_max_l reduce_l(2) = v_max_l reduce_l(3) = w_max_l IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MAX, comm2d, ierr ) u_max = reduce(1) v_max = reduce(2) w_max = reduce(3) reduce_l(1) = u_min_l reduce_l(2) = v_min_l reduce_l(3) = w_min_l IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) IF ( ABS( reduce(1) ) > u_max ) u_max = reduce(1) IF ( ABS( reduce(2) ) > v_max ) v_max = reduce(2) IF ( ABS( reduce(3) ) > w_max ) w_max = reduce(3) #else dt_u = dt_u_l dt_v = dt_v_l dt_w = dt_w_l IF ( ABS( u_min_l ) > u_max_l ) THEN u_max = u_min_l ELSE u_max = u_max_l ENDIF IF ( ABS( v_min_l ) > v_max_l ) THEN v_max = v_min_l ELSE v_max = v_max_l ENDIF IF ( ABS( w_min_l ) > w_max_l ) THEN w_max = w_min_l ELSE w_max = w_max_l ENDIF #endif #else ! !-- Variable time step: !-- Calculate the maximum time step according to the CFL-criterion, !-- individually for each velocity component dt_u_l = 999999.9 dt_v_l = 999999.9 dt_w_l = 999999.9 DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt dt_u_l = MIN( dt_u_l, ( dx / ( ABS( u(k,j,i) - u_gtrans ) + 1.0E-10 ) ) ) dt_v_l = MIN( dt_v_l, ( dy / ( ABS( v(k,j,i) - v_gtrans ) + 1.0E-10 ) ) ) dt_w_l = MIN( dt_w_l, ( dzu(k) / ( ABS( w(k,j,i) ) + 1.0E-10 ) ) ) ENDDO ENDDO ENDDO #if defined( __parallel ) reduce_l(1) = dt_u_l reduce_l(2) = dt_v_l reduce_l(3) = dt_w_l IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr ) dt_u = reduce(1) dt_v = reduce(2) dt_w = reduce(3) #else dt_u = dt_u_l dt_v = dt_v_l dt_w = dt_w_l #endif #endif ! !-- 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, dzw(k)*dzw(k) ) * 0.125 ENDDO !$OMP PARALLEL private(i,j,k,value) reduction(MIN: dt_diff_l) !$OMP DO !$acc parallel loop collapse(3) present( kh, km ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt dt_diff_l = MIN( dt_diff_l, dxyz2_min(k) / & ( MAX( kh(k,j,i), km(k,j,i) ) + 1E-20 ) ) ENDDO ENDDO ENDDO !$acc end parallel !$OMP END PARALLEL #if defined( __parallel ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( dt_diff_l, dt_diff, 1, MPI_REAL, MPI_MIN, comm2d, & ierr ) #else dt_diff = dt_diff_l #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 ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) 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 (cfl_factor) to be on the safe side. !-- The time step must not exceed the maximum allowed value. dt_3d = cfl_factor * MIN( dt_diff, dt_plant_canopy, dt_u, dt_v, dt_w, & dt_precipitation ) 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 IF ( myid == 0 ) THEN 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 CALL MPI_BCAST( terminate_coupled_remote, 1, MPI_INTEGER, 0, comm2d, ierr) ENDIF #endif ENDIF ! !-- Ensure a smooth value (two significant digits) of the timestep. div = 1000.0 DO WHILE ( dt_3d < div ) div = div / 10.0 ENDDO dt_3d = NINT( dt_3d * 100.0 / div ) * div / 100.0 ! !-- Adjust the time step old_dt = dt_3d ENDIF CALL cpu_log( log_point(12), 'calculate_timestep', 'stop' ) END SUBROUTINE timestep