Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

modifications for safer thread handling #289

Merged
merged 5 commits into from
Apr 21, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion .github/workflows/build.yml
Original file line number Diff line number Diff line change
Expand Up @@ -146,7 +146,7 @@ jobs:
run: >-
github-release ${{ env.RELEASE_COMMAND }}
--tag ${{ env.RELEASE_TAG }}
--name "Continous release version"
--name "Continuous release version"
--description "$DESCRIPTION"
--pre-release
env:
Expand Down
28 changes: 23 additions & 5 deletions src/calculator/calculator.F90
Original file line number Diff line number Diff line change
Expand Up @@ -138,8 +138,10 @@ subroutine engrad_mol(mol,calc,energy,gradient,iostatus)
iostatus = 0
dum1 = 1.0_wp
dum2 = 1.0_wp
if(n>0)then
calc%etmp = 0.0_wp
calc%grdtmp = 0.0_wp
endif
!$omp end critical

!**************************************
Expand Down Expand Up @@ -635,17 +637,33 @@ subroutine constrhess(nat,at,xyz,calc,phess)
real(wp) :: energy,el,er
real(wp),allocatable :: hess(:,:)
logical :: consgeo
integer :: i,j,k,n3,io
integer :: i,j,k,n3,io,ncons

if (calc%nconstraints <= 0) return
ncons = calc%nconstraints
if (ncons <= 0) return
!>--- skip if only nonadiabatic constraints

dummycalc = calc !> new dummy calculation
dummycalc%id = -1000 !> set to something arbitrary so that only constraints are considered
!> new dummy calculation
!> set to something arbitrary so that only constraints are considered
!> and copy the neccesities
dummycalc%id = -1000
dummycalc%ncalculations = 0
dummycalc%pr_energies = .false.

!> copy the constraints
dummycalc%nfreeze = calc%nfreeze
dummycalc%nconstraints = ncons
!$omp critical
allocate(dummycalc%cons( ncons ))
do i=1,ncons
dummycalc%cons(i) = calc%cons(i)
enddo
if(calc%nfreeze > 0)then
dummycalc%freezelist = calc%freezelist
endif
n3 = nat*3
allocate (hess(n3,n3),source=0.0_wp)
!$omp end critical

call numhess1(nat,at,xyz,dummycalc,hess,io)

Expand All @@ -656,7 +674,7 @@ subroutine constrhess(nat,at,xyz,calc,phess)
phess(k) = phess(k)+0.5_wp*(hess(j,i)+hess(i,j))
end do
end do

deallocate (hess)
return
end subroutine constrhess
Expand Down
26 changes: 13 additions & 13 deletions src/ompmklset.F90
Original file line number Diff line number Diff line change
Expand Up @@ -17,9 +17,9 @@
! along with crest. If not, see <https://www.gnu.org/licenses/>.
!================================================================================!

!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!
!c OMP and MKL parallelization settings
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!

subroutine ompmklset(threads)
use omp_lib
Expand All @@ -32,9 +32,9 @@ subroutine ompmklset(threads)
#endif
end subroutine ompmklset

!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!
!c OMP and MKL parallelization settings (short routine)
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!

subroutine ompenvset(omp)
use iomod
Expand All @@ -47,9 +47,10 @@ subroutine ompenvset(omp)

end subroutine ompenvset

!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!
!c OMP and MKL autoset switchcase routine
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!

subroutine new_ompautoset(env,modus,maxjobs,parallel_jobs,cores_per_job)
use omp_lib
use crest_data
Expand Down Expand Up @@ -82,12 +83,10 @@ subroutine new_ompautoset(env,modus,maxjobs,parallel_jobs,cores_per_job)
Tfloor = floor(Tfrac)
cores_per_job = max(nint(Tfloor),1)
end if
if (index(modus,'_nested') .ne. 0) then
if (index(modus,'_nested') .ne. 0 .and. cores_per_job > 1) then
if (env%omp_allow_nested) then
!> We should never need more than two active nested layers
call omp_set_max_active_levels(2)
!else
! parallel_jobs = T
end if
end if

Expand Down Expand Up @@ -115,9 +114,10 @@ subroutine new_ompautoset(env,modus,maxjobs,parallel_jobs,cores_per_job)

end subroutine new_ompautoset

!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!
!c get omp/mkl automatically from the global variables
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!

subroutine ompgetauto(threads,omp,maxrun)
use omp_lib
use iomod
Expand All @@ -138,9 +138,9 @@ subroutine ompgetauto(threads,omp,maxrun)

end subroutine ompgetauto

!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!
!c print omp/mkl threads that are used at the moment
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc!
subroutine ompprint_intern()
use omp_lib
implicit none
Expand Down
Loading