Skip to content

Commit

Permalink
fix: add omp directives for race condition in tridiagonal (#1109)
Browse files Browse the repository at this point in the history
  • Loading branch information
ganganoaa authored Mar 9, 2023
1 parent 9339b88 commit 0ff254e
Showing 1 changed file with 12 additions and 6 deletions.
18 changes: 12 additions & 6 deletions tridiagonal/tridiagonal.F90
Original file line number Diff line number Diff line change
Expand Up @@ -89,8 +89,10 @@ subroutine tri_invert(x,d,a,b,c)
integer :: k

if(present(a)) then
init_tridiagonal = .true.

!< Check if module variables are allocated
!$OMP SINGLE
init_tridiagonal = .true.
if(allocated(e)) deallocate(e)
if(allocated(g)) deallocate(g)
if(allocated(bb)) deallocate(bb)
Expand All @@ -99,6 +101,7 @@ subroutine tri_invert(x,d,a,b,c)
allocate(g (size(x,1),size(x,2),size(x,3)))
allocate(bb(size(x,1),size(x,2)))
allocate(cc(size(x,1),size(x,2),size(x,3)))
!$OMP END SINGLE !< There is an implicit barrier.

e(:,:,1) = - a(:,:,1)/b(:,:,1)
a(:,:,size(x,3)) = 0.0
Expand Down Expand Up @@ -132,12 +135,15 @@ end subroutine tri_invert
!> @brief Releases memory used by the solver
subroutine close_tridiagonal

implicit none
implicit none

deallocate(e)
deallocate(g)
deallocate(bb)
deallocate(cc)
!< Check if module variables are allocated
!$OMP SINGLE
if(allocated(e)) deallocate(e)
if(allocated(g)) deallocate(g)
if(allocated(bb)) deallocate(bb)
if(allocated(cc)) deallocate(cc)
!$OMP END SINGLE !< There is an implicit barrier.

return
end subroutine close_tridiagonal
Expand Down

0 comments on commit 0ff254e

Please sign in to comment.