SUBROUTINE swap_timelevel !--------------------------------------------------------------------------------! ! 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-2014 Leibniz Universitaet Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: swap_timelevel.f90 1321 2014-03-20 09:40:40Z witha $ ! ! 1320 2014-03-20 08:40:49Z raasch ! ONLY-attribute added to USE-statements, ! revision history before 2012 removed, ! 1318 2014-03-17 13:35:16Z raasch ! module interfaces removed ! ! 1115 2013-03-26 18:16:16Z hoffmann ! calculation of qr and nr is restricted to precipitation ! ! 1111 2013-03-08 23:54:10Z raasch ! openACC directives added ! ! 1053 2012-11-13 17:11:03Z hoffmann ! swap of timelevels for nr, qr added ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1032 2012-10-21 13:03:21Z letzel ! save memory by not allocating pt_2 in case of neutral = .T. ! ! 1010 2012-09-20 07:59:54Z raasch ! cpp switch __nopointer added for pointer free version ! ! 1001 2012-09-13 14:08:46Z raasch ! all actions concerning leapfrog scheme removed ! ! Revision 1.1 2000/01/10 10:08:58 10:08:58 raasch (Siegfried Raasch) ! Initial revision ! ! ! Description: ! ------------ ! Swap of timelevels of variables after each timestep !------------------------------------------------------------------------------! USE arrays_3d, & ONLY: e, e_1, e_2, e_p, nr, nr_1, nr_2, nr_p, pt, pt_1, pt_2, pt_p, q,& q_1, q_2, q_p, qr, qr_1, qr_2, qr_p, sa, sa_1, sa_2, sa_p, u, & u_1, u_2, u_p, v, v_1, v_2, v_p, w, w_1, w_2, w_p USE cpulog, & ONLY: cpu_log, log_point USE control_parameters, & ONLY: cloud_physics, constant_diffusion, humidity, icloud_scheme, & neutral, ocean, passive_scalar, precipitation, timestep_count IMPLICIT NONE ! !-- Incrementing timestep counter timestep_count = timestep_count + 1 ! !-- Swap of variables #if defined( __nopointer ) CALL cpu_log( log_point(28), 'swap_timelevel (nop)', 'start' ) !$acc kernels present( pt, pt_p, u, u_p, v, v_p, w, w_p ) u = u_p v = v_p w = w_p pt = pt_p !$acc end kernels IF ( .NOT. constant_diffusion ) THEN !$acc kernels present( e, e_p ) e = e_p !$acc end kernels ENDIF IF ( ocean ) THEN sa = sa_p ENDIF IF ( humidity .OR. passive_scalar ) THEN q = q_p IF ( cloud_physics .AND. icloud_scheme == 0 ) THEN qr = qr_p nr = nr_p ENDIF ENDIF CALL cpu_log( log_point(28), 'swap_timelevel (nop)', 'stop' ) #else CALL cpu_log( log_point(28), 'swap_timelevel', 'start' ) SELECT CASE ( MOD( timestep_count, 2 ) ) CASE ( 0 ) u => u_1; u_p => u_2 v => v_1; v_p => v_2 w => w_1; w_p => w_2 IF ( .NOT. neutral ) THEN pt => pt_1; pt_p => pt_2 ENDIF IF ( .NOT. constant_diffusion ) THEN e => e_1; e_p => e_2 ENDIF IF ( ocean ) THEN sa => sa_1; sa_p => sa_2 ENDIF IF ( humidity .OR. passive_scalar ) THEN q => q_1; q_p => q_2 IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & precipitation ) THEN qr => qr_1; qr_p => qr_2 nr => nr_1; nr_p => nr_2 ENDIF ENDIF CASE ( 1 ) u => u_2; u_p => u_1 v => v_2; v_p => v_1 w => w_2; w_p => w_1 IF ( .NOT. neutral ) THEN pt => pt_2; pt_p => pt_1 ENDIF IF ( .NOT. constant_diffusion ) THEN e => e_2; e_p => e_1 ENDIF IF ( ocean ) THEN sa => sa_2; sa_p => sa_1 ENDIF IF ( humidity .OR. passive_scalar ) THEN q => q_2; q_p => q_1 IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & precipitation) THEN qr => qr_2; qr_p => qr_1 nr => nr_2; nr_p => nr_1 ENDIF ENDIF END SELECT CALL cpu_log( log_point(28), 'swap_timelevel', 'stop' ) #endif END SUBROUTINE swap_timelevel