EVOLUTION: Parallelize evolve_ssp324() using OpenMP.

Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
Grzegorz Kowal 2021-12-08 18:27:43 -03:00
parent 537f8693bc
commit ee89bc42e6

View File

@ -2682,8 +2682,8 @@ module evolution
! !
subroutine evolve_ssp324() subroutine evolve_ssp324()
use blocks , only : block_data, list_data use blocks , only : block_data, list_data, data_blocks, get_dblocks
use boundaries , only : boundary_fluxes use boundaries, only : boundary_fluxes
use equations , only : errors, ibp, cmax use equations , only : errors, ibp, cmax
use helpers , only : print_message use helpers , only : print_message
use sources , only : update_sources use sources , only : update_sources
@ -2693,7 +2693,7 @@ module evolution
type(block_data), pointer :: pdata type(block_data), pointer :: pdata
logical :: test, cond logical :: test, cond
integer :: nrej, i, status integer :: i, l, n, nrej, status
real(kind=8) :: tm, dtm, dh, fc real(kind=8) :: tm, dtm, dh, fc
real(kind=8), parameter :: onethird = 1.0d+00 / 3.0d+00, & real(kind=8), parameter :: onethird = 1.0d+00 / 3.0d+00, &
@ -2706,6 +2706,7 @@ module evolution
test = .true. test = .true.
cond = .true. cond = .true.
nrej = 0 nrej = 0
n = get_dblocks()
! at the entry point we assume the previous solution of conserved variables U(n) ! at the entry point we assume the previous solution of conserved variables U(n)
! is stored in pdata%u(:,:,:,:,1) and the primitive variables are already ! is stored in pdata%u(:,:,:,:,1) and the primitive variables are already
@ -2723,15 +2724,15 @@ module evolution
!= preparation step: U(1) = U(n) != preparation step: U(1) = U(n)
! !
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%uu(:,:,:,:,2) = pdata%uu(:,:,:,:,1) pdata%uu(:,:,:,:,2) = pdata%uu(:,:,:,:,1)
pdata%u => pdata%uu(:,:,:,:,2) pdata%u => pdata%uu(:,:,:,:,2)
pdata => pdata%next
end do end do
!$omp end parallel do
!= 1st to 3rd steps: U(1) = U(1) + ½ dt F[U(1)], for i = 1...3 != 1st to 3rd steps: U(1) = U(1) + ½ dt F[U(1)], for i = 1...3
! !
@ -2740,24 +2741,24 @@ module evolution
tm = time + (i - 1) * dh tm = time + (i - 1) * dh
dtm = dh dtm = dh
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
call update_increment(pdata) call update_increment(pdata)
call update_sources(pdata, tm, dtm, pdata%du(:,:,:,:)) call update_sources(pdata, tm, dtm, pdata%du(:,:,:,:))
pdata => pdata%next
end do end do
!$omp end parallel do
call boundary_fluxes() call boundary_fluxes()
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%uu(:,:,:,:,2) = pdata%uu(:,:,:,:,2) + dh * pdata%du(:,:,:,:) pdata%uu(:,:,:,:,2) = pdata%uu(:,:,:,:,2) + dh * pdata%du(:,:,:,:)
pdata => pdata%next
end do end do
!$omp end parallel do
call update_variables(tm, dtm, status) call update_variables(tm, dtm, status)
if (status /= 0) go to 100 if (status /= 0) go to 100
@ -2766,16 +2767,16 @@ module evolution
!= 4th step: U(1) = U(n) + U(1), U(2) = U(n) + U(1) != 4th step: U(1) = U(n) + U(1), U(2) = U(n) + U(1)
! !
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%uu(:,:,:,:,3) = onethird * pdata%uu(:,:,:,:,1) & pdata%uu(:,:,:,:,3) = onethird * pdata%uu(:,:,:,:,1) &
+ twothird * pdata%uu(:,:,:,:,2) + twothird * pdata%uu(:,:,:,:,2)
pdata%uu(:,:,:,:,2) = twothird * pdata%uu(:,:,:,:,1) & pdata%uu(:,:,:,:,2) = twothird * pdata%uu(:,:,:,:,1) &
+ onethird * pdata%uu(:,:,:,:,2) + onethird * pdata%uu(:,:,:,:,2)
pdata => pdata%next
end do end do
!$omp end parallel do
call update_variables(tm, dtm, status) call update_variables(tm, dtm, status)
if (status /= 0) go to 100 if (status /= 0) go to 100
@ -2785,54 +2786,54 @@ module evolution
tm = time + dh tm = time + dh
dtm = dh dtm = dh
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
call update_increment(pdata) call update_increment(pdata)
call update_sources(pdata, tm, dtm, pdata%du(:,:,:,:)) call update_sources(pdata, tm, dtm, pdata%du(:,:,:,:))
pdata => pdata%next
end do end do
!$omp end parallel do
call boundary_fluxes() call boundary_fluxes()
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%uu(:,:,:,:,2) = pdata%uu(:,:,:,:,2) + dh * pdata%du(:,:,:,:) pdata%uu(:,:,:,:,2) = pdata%uu(:,:,:,:,2) + dh * pdata%du(:,:,:,:)
pdata => pdata%next
end do end do
!$omp end parallel do
call update_variables(tm, dtm, status) call update_variables(tm, dtm, status)
if (status /= 0) go to 100 if (status /= 0) go to 100
!= 6th step: U(2) = ½ (U(1) + U(2)) <- 2-order approximation != 6th step: U(2) = ½ (U(1) + U(2)) <- 2-order approximation
! !
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%uu(:,:,:,:,3) = 5.0d-01 * (pdata%uu(:,:,:,:,3) & pdata%uu(:,:,:,:,3) = 5.0d-01 * (pdata%uu(:,:,:,:,3) &
+ pdata%uu(:,:,:,:,2)) + pdata%uu(:,:,:,:,2))
pdata => pdata%next
end do end do
!$omp end parallel do
!= 7th step: decay the magnetic divergence potential of both solutions != 7th step: decay the magnetic divergence potential of both solutions
! !
if (ibp > 0) then if (ibp > 0) then
adecay(:) = exp(aglm(:) * cmax * dt) adecay(:) = exp(aglm(:) * cmax * dt)
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%uu(ibp,:,:,:,2) = adecay(pdata%meta%level) & pdata%uu(ibp,:,:,:,2) = adecay(pdata%meta%level) &
* pdata%uu(ibp,:,:,:,2) * pdata%uu(ibp,:,:,:,2)
pdata%uu(ibp,:,:,:,3) = adecay(pdata%meta%level) & pdata%uu(ibp,:,:,:,3) = adecay(pdata%meta%level) &
* pdata%uu(ibp,:,:,:,3) * pdata%uu(ibp,:,:,:,3)
pdata => pdata%next
end do end do
!$omp end parallel do
end if end if
! update the vector of errors ! update the vector of errors
@ -2869,13 +2870,13 @@ module evolution
! since the solution was rejected, we have to revert the primitive variables ! since the solution was rejected, we have to revert the primitive variables
! to the previous state ! to the previous state
! !
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%u => pdata%uu(:,:,:,:,1) pdata%u => pdata%uu(:,:,:,:,1)
pdata => pdata%next
end do end do
!$omp end parallel do
call update_variables(tm, dtm, status) call update_variables(tm, dtm, status)
if (status /= 0) & if (status /= 0) &
@ -2889,15 +2890,15 @@ module evolution
!= final step: U(n+1) = U(1) - update the accepted solution != final step: U(n+1) = U(1) - update the accepted solution
! !
pdata => list_data !$omp parallel do default(shared) private(pdata)
do while (associated(pdata)) do l = 1, n
pdata => data_blocks(l)%ptr
pdata%uu(:,:,:,:,1) = pdata%uu(:,:,:,:,2) pdata%uu(:,:,:,:,1) = pdata%uu(:,:,:,:,2)
pdata%u => pdata%uu(:,:,:,:,1) pdata%u => pdata%uu(:,:,:,:,1)
pdata => pdata%next
end do end do
!$omp end parallel do
!------------------------------------------------------------------------------- !-------------------------------------------------------------------------------
! !