SUBROUTINE pres !------------------------------------------------------------------------------! ! Actual revisions: ! ----------------- ! Volume flow conservation added for the remaining three outflow boundaries ! ! Former revisions: ! ----------------- ! $Id: pres.f90 106 2007-08-16 14:30:26Z letzel $ ! ! 85 2007-05-11 09:35:14Z raasch ! Division through dt_3d replaced by multiplication of the inverse. ! For performance optimisation, this is done in the loop calculating the ! divergence instead of using a seperate loop. ! ! 75 2007-03-22 09:54:05Z raasch ! Volume flow control for non-cyclic boundary conditions added (currently only ! for the north boundary!!), 2nd+3rd argument removed from exchange horiz, ! mean vertical velocity is removed in case of Neumann boundary conditions ! both at the bottom and the top ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.25 2006/04/26 13:26:12 raasch ! OpenMP optimization (+localsum, threadsum) ! ! Revision 1.1 1997/07/24 11:24:44 raasch ! Initial revision ! ! ! Description: ! ------------ ! Compute the divergence of the provisional velocity field. Solve the Poisson ! equation for the perturbation pressure. Compute the final velocities using ! this perturbation pressure. Compute the remaining divergence. !------------------------------------------------------------------------------! USE arrays_3d USE constants USE control_parameters USE cpulog USE grid_variables USE indices USE interfaces USE pegrid USE poisfft_mod USE poisfft_hybrid_mod USE statistics IMPLICIT NONE INTEGER :: i, j, k, sr REAL :: ddt_3d, localsum, threadsum REAL, DIMENSION(1:2) :: volume_flow_l, volume_flow_offset REAL, DIMENSION(1:nzt) :: w_l, w_l_l CALL cpu_log( log_point(8), 'pres', 'start' ) ddt_3d = 1.0 / dt_3d ! !-- Multigrid method needs additional grid points for the divergence array IF ( psolver == 'multigrid' ) THEN DEALLOCATE( d ) ALLOCATE( d(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ) ENDIF ! !-- Conserve the volume flow at the outflow in case of non-cyclic lateral !-- boundary conditions !-- WARNING: so far, this conservation does not work at the left/south !-- boundary if the topography at the inflow differs from that at the !-- outflow! For this case, volume_flow_area needs adjustment! ! !-- Left/right IF ( conserve_volume_flow .AND. ( outflow_l .OR. outflow_r ) ) THEN volume_flow(1) = 0.0 volume_flow_l(1) = 0.0 IF ( outflow_l ) THEN i = 0 ELSEIF ( outflow_r ) THEN i = nx+1 ENDIF DO j = nys, nyn ! !-- Sum up the volume flow through the south/north boundary DO k = nzb_2d(j,i) + 1, nzt volume_flow_l(1) = volume_flow_l(1) + u(k,j,i) * dzu(k) ENDDO ENDDO #if defined( __parallel ) CALL MPI_ALLREDUCE( volume_flow_l(1), volume_flow(1), 1, MPI_REAL, & MPI_SUM, comm1dy, ierr ) #else volume_flow = volume_flow_l #endif volume_flow_offset(1) = ( volume_flow_initial(1) - volume_flow(1) ) & / volume_flow_area(1) DO j = nys, nyn DO k = nzb_v_inner(j,i) + 1, nzt u(k,j,i) = u(k,j,i) + volume_flow_offset(1) ENDDO ENDDO CALL exchange_horiz( u ) ENDIF ! !-- South/north IF ( conserve_volume_flow .AND. ( outflow_n .OR. outflow_s ) ) THEN volume_flow(2) = 0.0 volume_flow_l(2) = 0.0 IF ( outflow_s ) THEN j = 0 ELSEIF ( outflow_n ) THEN j = ny+1 ENDIF DO i = nxl, nxr ! !-- Sum up the volume flow through the south/north boundary DO k = nzb_2d(j,i) + 1, nzt volume_flow_l(2) = volume_flow_l(2) + v(k,j,i) * dzu(k) ENDDO ENDDO #if defined( __parallel ) CALL MPI_ALLREDUCE( volume_flow_l(2), volume_flow(2), 1, MPI_REAL, & MPI_SUM, comm1dx, ierr ) #else volume_flow = volume_flow_l #endif volume_flow_offset(2) = ( volume_flow_initial(2) - volume_flow(2) ) & / volume_flow_area(2) DO i = nxl, nxr DO k = nzb_v_inner(j,i) + 1, nzt v(k,j,i) = v(k,j,i) + volume_flow_offset(2) ENDDO ENDDO CALL exchange_horiz( v ) ENDIF ! !-- Remove mean vertical velocity IF ( ibc_p_b == 1 .AND. ibc_p_t == 1 ) THEN IF ( simulated_time > 0.0 ) THEN ! otherwise nzb_w_inner is not yet known w_l = 0.0; w_l_l = 0.0 DO i = nxl, nxr DO j = nys, nyn DO k = nzb_w_inner(j,i)+1, nzt w_l_l(k) = w_l_l(k) + w(k,j,i) ENDDO ENDDO ENDDO #if defined( __parallel ) CALL MPI_ALLREDUCE( w_l_l(1), w_l(1), nzt, MPI_REAL, MPI_SUM, comm2d, & ierr ) #else w_l = w_l_l #endif DO k = 1, nzt w_l(k) = w_l(k) / ngp_2dh_outer(k,0) ENDDO DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 DO k = nzb_w_inner(j,i)+1, nzt w(k,j,i) = w(k,j,i) - w_l(k) ENDDO ENDDO ENDDO ENDIF ENDIF ! !-- Compute the divergence of the provisional velocity field. CALL cpu_log( log_point_s(1), 'divergence', 'start' ) IF ( psolver == 'multigrid' ) THEN !$OMP PARALLEL DO SCHEDULE( STATIC ) DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 DO k = nzb, nzt+1 d(k,j,i) = 0.0 ENDDO ENDDO ENDDO ELSE !$OMP PARALLEL DO SCHEDULE( STATIC ) DO i = nxl, nxra DO j = nys, nyna DO k = nzb+1, nzta d(k,j,i) = 0.0 ENDDO ENDDO ENDDO ENDIF localsum = 0.0 threadsum = 0.0 #if defined( __ibm ) !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum) !$OMP DO SCHEDULE( STATIC ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * ddx + & ( v(k,j+1,i) - v(k,j,i) ) * ddy + & ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) ) * ddt_3d ENDDO ! !-- Additional pressure boundary condition at the bottom boundary for !-- inhomogeneous Prandtl layer heat fluxes and temperatures, respectively !-- dp/dz = -(dtau13/dx + dtau23/dy) + g*pt'/pt0. !-- This condition must not be applied at the start of a run, because then !-- flow_statistics has not yet been called and thus sums = 0. IF ( ibc_p_b == 2 .AND. sums(nzb+1,4) /= 0.0 ) THEN k = nzb_s_inner(j,i) d(k+1,j,i) = d(k+1,j,i) + ( & ( usws(j,i+1) - usws(j,i) ) * ddx & + ( vsws(j+1,i) - vsws(j,i) ) * ddy & - g * ( pt(k+1,j,i) - sums(k+1,4) ) / & sums(k+1,4) & ) * ddzw(k+1) * ddt_3d ENDIF ! !-- Compute possible PE-sum of divergences for flow_statistics DO k = nzb_s_inner(j,i)+1, nzt threadsum = threadsum + ABS( d(k,j,i) ) ENDDO ENDDO ENDDO localsum = ( localsum + threadsum ) * dt_3d !$OMP END PARALLEL #else IF ( ibc_p_b == 2 .AND. sums(nzb+1,4) /= 0.0 ) THEN !$OMP PARALLEL PRIVATE (i,j,k) !$OMP DO SCHEDULE( STATIC ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * ddx + & ( v(k,j+1,i) - v(k,j,i) ) * ddy + & ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) ) * ddt_3d ENDDO ENDDO ! !-- Additional pressure boundary condition at the bottom boundary for !-- inhomogeneous Prandtl layer heat fluxes and temperatures, respectively !-- dp/dz = -(dtau13/dx + dtau23/dy) + g*pt'/pt0. !-- This condition must not be applied at the start of a run, because then !-- flow_statistics has not yet been called and thus sums = 0. DO j = nys, nyn k = nzb_s_inner(j,i) d(k+1,j,i) = d(k+1,j,i) + ( & ( usws(j,i+1) - usws(j,i) ) * ddx & + ( vsws(j+1,i) - vsws(j,i) ) * ddy & - g * ( pt(k+1,j,i) - sums(k+1,4) ) / & sums(k+1,4) & ) * ddzw(k+1) * ddt_3d ENDDO ENDDO !$OMP END PARALLEL ELSE !$OMP PARALLEL PRIVATE (i,j,k) !$OMP DO SCHEDULE( STATIC ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * ddx + & ( v(k,j+1,i) - v(k,j,i) ) * ddy + & ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) ) * ddt_3d ENDDO ENDDO ENDDO !$OMP END PARALLEL ENDIF ! !-- Compute possible PE-sum of divergences for flow_statistics !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum) !$OMP DO SCHEDULE( STATIC ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt threadsum = threadsum + ABS( d(k,j,i) ) ENDDO ENDDO ENDDO localsum = ( localsum + threadsum ) * dt_3d !$OMP END PARALLEL #endif ! !-- For completeness, set the divergence sum of all statistic regions to those !-- of the total domain sums_divold_l(0:statistic_regions) = localsum ! !-- Determine absolute minimum/maximum (only for test cases, therefore as !-- comment line) ! CALL global_min_max( nzb+1, nzt, nys, nyn, nxl, nxr, d, 'abs', divmax, & ! divmax_ijk ) CALL cpu_log( log_point_s(1), 'divergence', 'stop' ) ! !-- Compute the pressure perturbation solving the Poisson equation IF ( psolver(1:7) == 'poisfft' ) THEN ! !-- Enlarge the size of tend, used as a working array for the transpositions IF ( nxra > nxr .OR. nyna > nyn .OR. nza > nz ) THEN DEALLOCATE( tend ) ALLOCATE( tend(1:nza,nys:nyna,nxl:nxra) ) ENDIF ! !-- Solve Poisson equation via FFT and solution of tridiagonal matrices IF ( psolver == 'poisfft' ) THEN ! !-- Solver for 2d-decomposition CALL poisfft( d, tend ) ELSEIF ( psolver == 'poisfft_hybrid' ) THEN ! !-- Solver for 1d-decomposition (using MPI and OpenMP). !-- The old hybrid-solver is still included here, as long as there !-- are some optimization problems in poisfft CALL poisfft_hybrid( d ) ENDIF ! !-- Resize tend to its normal size IF ( nxra > nxr .OR. nyna > nyn .OR. nza > nz ) THEN DEALLOCATE( tend ) ALLOCATE( tend(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ) ENDIF ! !-- Store computed perturbation pressure and set boundary condition in !-- z-direction !$OMP PARALLEL DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt tend(k,j,i) = d(k,j,i) ENDDO ENDDO ENDDO ! !-- Bottom boundary: !-- This condition is only required for internal output. The pressure !-- gradient (dp(nzb+1)-dp(nzb))/dz is not used anywhere else. IF ( ibc_p_b == 1 ) THEN ! !-- Neumann (dp/dz = 0) !$OMP PARALLEL DO DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 tend(nzb_s_inner(j,i),j,i) = tend(nzb_s_inner(j,i)+1,j,i) ENDDO ENDDO ELSEIF ( ibc_p_b == 2 ) THEN ! !-- Neumann condition for inhomogeneous surfaces, !-- here currently still in the form of a zero gradient. Actually !-- dp/dz = -(dtau13/dx + dtau23/dy) + g*pt'/pt0 would have to be used for !-- the computation (cf. above: computation of divergences). !$OMP PARALLEL DO DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 tend(nzb_s_inner(j,i),j,i) = tend(nzb_s_inner(j,i)+1,j,i) ENDDO ENDDO ELSE ! !-- Dirichlet !$OMP PARALLEL DO DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 tend(nzb_s_inner(j,i),j,i) = 0.0 ENDDO ENDDO ENDIF ! !-- Top boundary IF ( ibc_p_t == 1 ) THEN ! !-- Neumann !$OMP PARALLEL DO DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 tend(nzt+1,j,i) = tend(nzt,j,i) ENDDO ENDDO ELSE ! !-- Dirichlet !$OMP PARALLEL DO DO i = nxl-1, nxr+1 DO j = nys-1, nyn+1 tend(nzt+1,j,i) = 0.0 ENDDO ENDDO ENDIF ! !-- Exchange boundaries for p CALL exchange_horiz( tend ) ELSEIF ( psolver == 'sor' ) THEN ! !-- Solve Poisson equation for perturbation pressure using SOR-Red/Black !-- scheme CALL sor( d, ddzu, ddzw, p ) tend = p ELSEIF ( psolver == 'multigrid' ) THEN ! !-- Solve Poisson equation for perturbation pressure using Multigrid scheme, !-- array tend is used to store the residuals CALL poismg( tend ) ! !-- Restore perturbation pressure on tend because this array is used !-- further below to correct the velocity fields tend = p ENDIF ! !-- Store perturbation pressure on array p, used in the momentum equations IF ( psolver(1:7) == 'poisfft' ) THEN ! !-- Here, only the values from the left and right boundaries are copied !-- The remaining values are copied in the following loop due to speed !-- optimization !$OMP PARALLEL DO DO j = nys-1, nyn+1 DO k = nzb, nzt+1 p(k,j,nxl-1) = tend(k,j,nxl-1) p(k,j,nxr+1) = tend(k,j,nxr+1) ENDDO ENDDO ENDIF ! !-- Correction of the provisional velocities with the current perturbation !-- pressure just computed IF ( conserve_volume_flow .AND. & ( bc_lr == 'cyclic' .OR. bc_ns == 'cyclic' ) ) THEN volume_flow_l(1) = 0.0 volume_flow_l(2) = 0.0 ENDIF !$OMP PARALLEL PRIVATE (i,j,k) !$OMP DO DO i = nxl, nxr IF ( psolver(1:7) == 'poisfft' ) THEN DO j = nys-1, nyn+1 DO k = nzb, nzt+1 p(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDIF DO j = nys, nyn DO k = nzb_w_inner(j,i)+1, nzt w(k,j,i) = w(k,j,i) - dt_3d * & ( tend(k+1,j,i) - tend(k,j,i) ) * ddzu(k+1) ENDDO DO k = nzb_u_inner(j,i)+1, nzt u(k,j,i) = u(k,j,i) - dt_3d * ( tend(k,j,i) - tend(k,j,i-1) ) * ddx ENDDO DO k = nzb_v_inner(j,i)+1, nzt v(k,j,i) = v(k,j,i) - dt_3d * ( tend(k,j,i) - tend(k,j-1,i) ) * ddy ENDDO ! !-- Sum up the volume flow through the right and north boundary IF ( conserve_volume_flow .AND. bc_lr == 'cyclic' .AND. & i == nx ) THEN !$OMP CRITICAL DO k = nzb_2d(j,i) + 1, nzt volume_flow_l(1) = volume_flow_l(1) + u(k,j,i) * dzu(k) ENDDO !$OMP END CRITICAL ENDIF IF ( conserve_volume_flow .AND. bc_ns == 'cyclic' .AND. & j == ny ) THEN !$OMP CRITICAL DO k = nzb_2d(j,i) + 1, nzt volume_flow_l(2) = volume_flow_l(2) + v(k,j,i) * dzu(k) ENDDO !$OMP END CRITICAL ENDIF ENDDO ENDDO !$OMP END PARALLEL ! !-- Conserve the volume flow IF ( conserve_volume_flow .AND. & ( bc_lr == 'cyclic' .OR. bc_ns == 'cyclic' ) ) THEN #if defined( __parallel ) CALL MPI_ALLREDUCE( volume_flow_l(1), volume_flow(1), 2, MPI_REAL, & MPI_SUM, comm2d, ierr ) #else volume_flow = volume_flow_l #endif volume_flow_offset = ( volume_flow_initial - volume_flow ) / & volume_flow_area !$OMP PARALLEL PRIVATE (i,j,k) !$OMP DO DO i = nxl, nxr DO j = nys, nyn IF ( bc_lr == 'cyclic' ) THEN DO k = nzb_u_inner(j,i) + 1, nzt u(k,j,i) = u(k,j,i) + volume_flow_offset(1) ENDDO ENDIF IF ( bc_ns == 'cyclic' ) THEN DO k = nzb_v_inner(j,i) + 1, nzt v(k,j,i) = v(k,j,i) + volume_flow_offset(2) ENDDO ENDIF ENDDO ENDDO !$OMP END PARALLEL ENDIF ! !-- Exchange of boundaries for the velocities CALL exchange_horiz( u ) CALL exchange_horiz( v ) CALL exchange_horiz( w ) ! !-- Compute the divergence of the corrected velocity field, !-- a possible PE-sum is computed in flow_statistics CALL cpu_log( log_point_s(1), 'divergence', 'start' ) sums_divnew_l = 0.0 ! !-- d must be reset to zero because it can contain nonzero values below the !-- topography IF ( topography /= 'flat' ) d = 0.0 localsum = 0.0 threadsum = 0.0 !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum) !$OMP DO SCHEDULE( STATIC ) #if defined( __ibm ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt d(k,j,i) = ( u(k,j,i+1) - u(k,j,i) ) * ddx + & ( v(k,j+1,i) - v(k,j,i) ) * ddy + & ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) ENDDO DO k = nzb+1, nzt threadsum = threadsum + ABS( d(k,j,i) ) ENDDO ENDDO ENDDO #else DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt d(k,j,i) = ( u(k,j,i+1) - u(k,j,i) ) * ddx + & ( v(k,j+1,i) - v(k,j,i) ) * ddy + & ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) threadsum = threadsum + ABS( d(k,j,i) ) ENDDO ENDDO ENDDO #endif localsum = localsum + threadsum !$OMP END PARALLEL ! !-- For completeness, set the divergence sum of all statistic regions to those !-- of the total domain sums_divnew_l(0:statistic_regions) = localsum CALL cpu_log( log_point_s(1), 'divergence', 'stop' ) CALL cpu_log( log_point(8), 'pres', 'stop' ) END SUBROUTINE pres