mirror of
https://github.com/LCPQ/quantum_package
synced 2024-11-03 20:54:00 +01:00
Add the possibility to abort cleanly a running Davidson or integrals
This commit is contained in:
parent
4fca291344
commit
9af8b74f4a
@ -196,7 +196,7 @@ BEGIN_PROVIDER [ logical, ao_bielec_integrals_in_map ]
|
|||||||
|
|
||||||
integer :: i,j,k,l
|
integer :: i,j,k,l
|
||||||
double precision :: ao_bielec_integral,cpu_1,cpu_2, wall_1, wall_2
|
double precision :: ao_bielec_integral,cpu_1,cpu_2, wall_1, wall_2
|
||||||
double precision :: integral
|
double precision :: integral, wall_0
|
||||||
double precision :: thresh
|
double precision :: thresh
|
||||||
thresh = ao_integrals_threshold
|
thresh = ao_integrals_threshold
|
||||||
|
|
||||||
@ -206,7 +206,7 @@ BEGIN_PROVIDER [ logical, ao_bielec_integrals_in_map ]
|
|||||||
real(integral_kind),allocatable :: buffer_value(:)
|
real(integral_kind),allocatable :: buffer_value(:)
|
||||||
integer(omp_lock_kind) :: lock
|
integer(omp_lock_kind) :: lock
|
||||||
|
|
||||||
integer :: n_integrals, n_centers
|
integer :: n_integrals, n_centers, thread_num
|
||||||
integer :: jl_pairs(2,ao_num*(ao_num+1)/2), kk, m, j1, i1, lmax
|
integer :: jl_pairs(2,ao_num*(ao_num+1)/2), kk, m, j1, i1, lmax
|
||||||
|
|
||||||
PROVIDE gauleg_t2 ao_integrals_map all_utils
|
PROVIDE gauleg_t2 ao_integrals_map all_utils
|
||||||
@ -235,22 +235,28 @@ BEGIN_PROVIDER [ logical, ao_bielec_integrals_in_map ]
|
|||||||
call omp_init_lock(lock)
|
call omp_init_lock(lock)
|
||||||
lmax = ao_num*(ao_num+1)/2
|
lmax = ao_num*(ao_num+1)/2
|
||||||
write(output_BiInts,*) 'providing the AO integrals'
|
write(output_BiInts,*) 'providing the AO integrals'
|
||||||
|
call wall_time(wall_0)
|
||||||
call wall_time(wall_1)
|
call wall_time(wall_1)
|
||||||
call cpu_time(cpu_1)
|
call cpu_time(cpu_1)
|
||||||
!$OMP PARALLEL PRIVATE(i,j,k,l,kk, &
|
!$OMP PARALLEL PRIVATE(i,j,k,l,kk, &
|
||||||
!$OMP integral,buffer_i,buffer_value,n_integrals, &
|
!$OMP integral,buffer_i,buffer_value,n_integrals, &
|
||||||
!$OMP cpu_2,wall_2,i1,j1) &
|
!$OMP cpu_2,wall_2,i1,j1,thread_num) &
|
||||||
!$OMP DEFAULT(NONE) &
|
!$OMP DEFAULT(NONE) &
|
||||||
!$OMP SHARED (ao_num, jl_pairs, ao_integrals_map,thresh, &
|
!$OMP SHARED (ao_num, jl_pairs, ao_integrals_map,thresh, &
|
||||||
!$OMP cpu_1,wall_1,lock, lmax,n_centers,ao_nucl, &
|
!$OMP cpu_1,wall_1,lock, lmax,n_centers,ao_nucl, &
|
||||||
!$OMP ao_overlap_abs,ao_overlap,output_BiInts)
|
!$OMP ao_overlap_abs,ao_overlap,output_BiInts,abort_here, &
|
||||||
|
!$OMP wall_0)
|
||||||
|
|
||||||
allocate(buffer_i(size_buffer))
|
allocate(buffer_i(size_buffer))
|
||||||
allocate(buffer_value(size_buffer))
|
allocate(buffer_value(size_buffer))
|
||||||
n_integrals = 0
|
n_integrals = 0
|
||||||
|
!$ thread_num = omp_get_thread_num()
|
||||||
|
|
||||||
!$OMP DO SCHEDULE(dynamic)
|
!$OMP DO SCHEDULE(dynamic)
|
||||||
do kk=1,lmax
|
do kk=1,lmax
|
||||||
|
if (abort_here) then
|
||||||
|
cycle
|
||||||
|
endif
|
||||||
j = jl_pairs(1,kk)
|
j = jl_pairs(1,kk)
|
||||||
l = jl_pairs(2,kk)
|
l = jl_pairs(2,kk)
|
||||||
j1 = j+ishft(l*l-l,-1)
|
j1 = j+ishft(l*l-l,-1)
|
||||||
@ -293,7 +299,14 @@ BEGIN_PROVIDER [ logical, ao_bielec_integrals_in_map ]
|
|||||||
enddo
|
enddo
|
||||||
enddo
|
enddo
|
||||||
call wall_time(wall_2)
|
call wall_time(wall_2)
|
||||||
write(output_BiInts,*) 100.*float(kk)/float(lmax), '% in ', wall_2-wall_1, 's'
|
|
||||||
|
if (thread_num == 0) then
|
||||||
|
if (wall_2 - wall_0 > 1.d0) then
|
||||||
|
wall_0 = wall_2
|
||||||
|
write(output_BiInts,*) 100.*float(kk)/float(lmax), '% in ', &
|
||||||
|
wall_2-wall_1, 's', map_mb(ao_integrals_map) ,'MB'
|
||||||
|
endif
|
||||||
|
endif
|
||||||
enddo
|
enddo
|
||||||
!$OMP END DO NOWAIT
|
!$OMP END DO NOWAIT
|
||||||
call insert_into_ao_integrals_map(n_integrals,buffer_i,buffer_value)
|
call insert_into_ao_integrals_map(n_integrals,buffer_i,buffer_value)
|
||||||
@ -301,6 +314,9 @@ BEGIN_PROVIDER [ logical, ao_bielec_integrals_in_map ]
|
|||||||
deallocate(buffer_value)
|
deallocate(buffer_value)
|
||||||
!$OMP END PARALLEL
|
!$OMP END PARALLEL
|
||||||
call omp_destroy_lock(lock)
|
call omp_destroy_lock(lock)
|
||||||
|
if (abort_here) then
|
||||||
|
stop 'Aborting in AO integrals calculation'
|
||||||
|
endif
|
||||||
write(output_BiInts,*) 'Sorting the map'
|
write(output_BiInts,*) 'Sorting the map'
|
||||||
call map_sort(ao_integrals_map)
|
call map_sort(ao_integrals_map)
|
||||||
call cpu_time(cpu_2)
|
call cpu_time(cpu_2)
|
||||||
|
@ -49,7 +49,7 @@ subroutine add_integrals_to_map(mask_ijkl)
|
|||||||
|
|
||||||
integer :: i,j,k,l
|
integer :: i,j,k,l
|
||||||
integer :: i0,j0,k0,l0
|
integer :: i0,j0,k0,l0
|
||||||
double precision :: c, cpu_1, cpu_2, wall_1, wall_2
|
double precision :: c, cpu_1, cpu_2, wall_1, wall_2, wall_0
|
||||||
|
|
||||||
integer, allocatable :: list_ijkl(:,:)
|
integer, allocatable :: list_ijkl(:,:)
|
||||||
integer :: n_i, n_j, n_k, n_l
|
integer :: n_i, n_j, n_k, n_l
|
||||||
@ -66,7 +66,7 @@ subroutine add_integrals_to_map(mask_ijkl)
|
|||||||
real(integral_kind),allocatable :: buffer_value(:)
|
real(integral_kind),allocatable :: buffer_value(:)
|
||||||
real :: map_mb
|
real :: map_mb
|
||||||
|
|
||||||
integer :: i1,j1,k1,l1, ii1, kmax, l1_global
|
integer :: i1,j1,k1,l1, ii1, kmax, thread_num
|
||||||
integer :: i2,i3,i4
|
integer :: i2,i3,i4
|
||||||
double precision,parameter :: thr_coef = 0.d0
|
double precision,parameter :: thr_coef = 0.d0
|
||||||
|
|
||||||
@ -81,7 +81,6 @@ subroutine add_integrals_to_map(mask_ijkl)
|
|||||||
call bitstring_to_list( mask_ijkl(1,3), list_ijkl(1,3), n_k, N_int )
|
call bitstring_to_list( mask_ijkl(1,3), list_ijkl(1,3), n_k, N_int )
|
||||||
call bitstring_to_list( mask_ijkl(1,4), list_ijkl(1,4), n_l, N_int )
|
call bitstring_to_list( mask_ijkl(1,4), list_ijkl(1,4), n_l, N_int )
|
||||||
|
|
||||||
l1_global = 0
|
|
||||||
size_buffer = min(ao_num*ao_num*ao_num,16000000)
|
size_buffer = min(ao_num*ao_num*ao_num,16000000)
|
||||||
write(output_BiInts,*) 'Providing the molecular integrals '
|
write(output_BiInts,*) 'Providing the molecular integrals '
|
||||||
write(output_BiInts,*) 'Buffers : ', 8.*(mo_tot_num_align*(n_j)*(n_k+1) + mo_tot_num_align +&
|
write(output_BiInts,*) 'Buffers : ', 8.*(mo_tot_num_align*(n_j)*(n_k+1) + mo_tot_num_align +&
|
||||||
@ -93,13 +92,14 @@ subroutine add_integrals_to_map(mask_ijkl)
|
|||||||
|
|
||||||
!$OMP PARALLEL PRIVATE(l1,k1,j1,i1,i2,i3,i4,i,j,k,l,c, ii1,kmax, &
|
!$OMP PARALLEL PRIVATE(l1,k1,j1,i1,i2,i3,i4,i,j,k,l,c, ii1,kmax, &
|
||||||
!$OMP bielec_tmp_0_idx, bielec_tmp_0, bielec_tmp_1,bielec_tmp_2,bielec_tmp_3,&
|
!$OMP bielec_tmp_0_idx, bielec_tmp_0, bielec_tmp_1,bielec_tmp_2,bielec_tmp_3,&
|
||||||
!$OMP buffer_i,buffer_value,n_integrals,wall_2,i0,j0,k0,l0) &
|
!$OMP buffer_i,buffer_value,n_integrals,wall_2,i0,j0,k0,l0, &
|
||||||
|
!$OMP wall_0,thread_num) &
|
||||||
!$OMP DEFAULT(NONE) &
|
!$OMP DEFAULT(NONE) &
|
||||||
!$OMP SHARED(size_buffer,ao_num,mo_tot_num,n_i,n_j,n_k,n_l,mo_tot_num_align,&
|
!$OMP SHARED(size_buffer,ao_num,mo_tot_num,n_i,n_j,n_k,n_l,mo_tot_num_align,&
|
||||||
!$OMP mo_coef_transp,output_BiInts, &
|
!$OMP mo_coef_transp,output_BiInts, &
|
||||||
!$OMP mo_coef_transp_is_built, list_ijkl, &
|
!$OMP mo_coef_transp_is_built, list_ijkl, &
|
||||||
!$OMP mo_coef_is_built, wall_1, &
|
!$OMP mo_coef_is_built, wall_1, abort_here, &
|
||||||
!$OMP mo_coef,mo_integrals_threshold,l1_global,ao_integrals_map,mo_integrals_map)
|
!$OMP mo_coef,mo_integrals_threshold,ao_integrals_map,mo_integrals_map)
|
||||||
n_integrals = 0
|
n_integrals = 0
|
||||||
allocate(bielec_tmp_3(mo_tot_num_align, n_j, n_k), &
|
allocate(bielec_tmp_3(mo_tot_num_align, n_j, n_k), &
|
||||||
bielec_tmp_1(mo_tot_num_align), &
|
bielec_tmp_1(mo_tot_num_align), &
|
||||||
@ -109,8 +109,12 @@ subroutine add_integrals_to_map(mask_ijkl)
|
|||||||
buffer_i(size_buffer), &
|
buffer_i(size_buffer), &
|
||||||
buffer_value(size_buffer) )
|
buffer_value(size_buffer) )
|
||||||
|
|
||||||
|
!$ thread_num = omp_get_thread_num()
|
||||||
!$OMP DO SCHEDULE(guided)
|
!$OMP DO SCHEDULE(guided)
|
||||||
do l1 = 1,ao_num
|
do l1 = 1,ao_num
|
||||||
|
if (abort_here) then
|
||||||
|
cycle
|
||||||
|
endif
|
||||||
!DEC$ VECTOR ALIGNED
|
!DEC$ VECTOR ALIGNED
|
||||||
bielec_tmp_3 = 0.d0
|
bielec_tmp_3 = 0.d0
|
||||||
do k1 = 1,ao_num
|
do k1 = 1,ao_num
|
||||||
@ -253,11 +257,14 @@ subroutine add_integrals_to_map(mask_ijkl)
|
|||||||
enddo
|
enddo
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
!$OMP ATOMIC
|
|
||||||
l1_global +=1
|
|
||||||
call wall_time(wall_2)
|
call wall_time(wall_2)
|
||||||
write(output_BiInts,*) 100.*float(l1_global)/float(ao_num), '% in ',&
|
if (thread_num == 0) then
|
||||||
wall_2-wall_1, 's', map_mb(mo_integrals_map) ,'MB'
|
if (wall_2 - wall_0 > 1.d0) then
|
||||||
|
wall_0 = wall_2
|
||||||
|
write(output_BiInts,*) 100.*float(l1)/float(ao_num), '% in ', &
|
||||||
|
wall_2-wall_1, 's', map_mb(ao_integrals_map) ,'MB'
|
||||||
|
endif
|
||||||
|
endif
|
||||||
enddo
|
enddo
|
||||||
!$OMP END DO NOWAIT
|
!$OMP END DO NOWAIT
|
||||||
deallocate (bielec_tmp_1,bielec_tmp_2,bielec_tmp_3)
|
deallocate (bielec_tmp_1,bielec_tmp_2,bielec_tmp_3)
|
||||||
@ -266,6 +273,9 @@ subroutine add_integrals_to_map(mask_ijkl)
|
|||||||
real(mo_integrals_threshold,integral_kind))
|
real(mo_integrals_threshold,integral_kind))
|
||||||
deallocate(buffer_i, buffer_value)
|
deallocate(buffer_i, buffer_value)
|
||||||
!$OMP END PARALLEL
|
!$OMP END PARALLEL
|
||||||
|
if (abort_here) then
|
||||||
|
stop 'Aborting in MO integrals calculation'
|
||||||
|
endif
|
||||||
call map_unique(mo_integrals_map)
|
call map_unique(mo_integrals_map)
|
||||||
|
|
||||||
call wall_time(wall_2)
|
call wall_time(wall_2)
|
||||||
@ -327,7 +337,7 @@ end
|
|||||||
!$OMP PRIVATE (i,j,p,q,r,s,integral,c,n,pp,int_value,int_idx, &
|
!$OMP PRIVATE (i,j,p,q,r,s,integral,c,n,pp,int_value,int_idx, &
|
||||||
!$OMP iqrs, iqsr,iqri,iqis) &
|
!$OMP iqrs, iqsr,iqri,iqis) &
|
||||||
!$OMP SHARED(mo_tot_num,mo_coef_transp,mo_tot_num_align,ao_num,&
|
!$OMP SHARED(mo_tot_num,mo_coef_transp,mo_tot_num_align,ao_num,&
|
||||||
!$OMP ao_integrals_threshold,do_direct_integrals) &
|
!$OMP ao_integrals_threshold,do_direct_integrals,abort_here) &
|
||||||
!$OMP REDUCTION(+:mo_bielec_integral_jj,mo_bielec_integral_jj_exchange)
|
!$OMP REDUCTION(+:mo_bielec_integral_jj,mo_bielec_integral_jj_exchange)
|
||||||
|
|
||||||
allocate( int_value(ao_num), int_idx(ao_num), &
|
allocate( int_value(ao_num), int_idx(ao_num), &
|
||||||
@ -336,6 +346,9 @@ end
|
|||||||
|
|
||||||
!$OMP DO SCHEDULE (guided)
|
!$OMP DO SCHEDULE (guided)
|
||||||
do s=1,ao_num
|
do s=1,ao_num
|
||||||
|
if (abort_here) then
|
||||||
|
cycle
|
||||||
|
endif
|
||||||
do q=1,ao_num
|
do q=1,ao_num
|
||||||
|
|
||||||
do j=1,ao_num
|
do j=1,ao_num
|
||||||
@ -421,6 +434,9 @@ end
|
|||||||
!$OMP END DO NOWAIT
|
!$OMP END DO NOWAIT
|
||||||
deallocate(iqrs,iqsr,int_value,int_idx)
|
deallocate(iqrs,iqsr,int_value,int_idx)
|
||||||
!$OMP END PARALLEL
|
!$OMP END PARALLEL
|
||||||
|
if (abort_here) then
|
||||||
|
stop 'Aborting in MO integrals calculation'
|
||||||
|
endif
|
||||||
|
|
||||||
mo_bielec_integral_jj_anti = mo_bielec_integral_jj - mo_bielec_integral_jj_exchange
|
mo_bielec_integral_jj_anti = mo_bielec_integral_jj - mo_bielec_integral_jj_exchange
|
||||||
|
|
||||||
|
@ -89,8 +89,14 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2 $parame
|
|||||||
accu = 0.d0
|
accu = 0.d0
|
||||||
do ispin=1,2
|
do ispin=1,2
|
||||||
other_spin = iand(ispin,1)+1
|
other_spin = iand(ispin,1)+1
|
||||||
|
if (abort_here) then
|
||||||
|
exit
|
||||||
|
endif
|
||||||
$omp_do
|
$omp_do
|
||||||
do ii=1,ia_ja_pairs(1,0,ispin)
|
do ii=1,ia_ja_pairs(1,0,ispin)
|
||||||
|
if (abort_here) then
|
||||||
|
cycle
|
||||||
|
endif
|
||||||
i_a = ia_ja_pairs(1,ii,ispin)
|
i_a = ia_ja_pairs(1,ii,ispin)
|
||||||
ASSERT (i_a > 0)
|
ASSERT (i_a > 0)
|
||||||
ASSERT (i_a <= mo_tot_num)
|
ASSERT (i_a <= mo_tot_num)
|
||||||
@ -152,6 +158,9 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2 $parame
|
|||||||
endif
|
endif
|
||||||
! endif
|
! endif
|
||||||
enddo
|
enddo
|
||||||
|
if (abort_here) then
|
||||||
|
exit
|
||||||
|
endif
|
||||||
enddo
|
enddo
|
||||||
endif
|
endif
|
||||||
! does all the mono excitations of the same spin
|
! does all the mono excitations of the same spin
|
||||||
@ -186,7 +195,9 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2 $parame
|
|||||||
$keys_work
|
$keys_work
|
||||||
key_idx = 0
|
key_idx = 0
|
||||||
endif
|
endif
|
||||||
! endif
|
if (abort_here) then
|
||||||
|
exit
|
||||||
|
endif
|
||||||
enddo
|
enddo
|
||||||
enddo! kk
|
enddo! kk
|
||||||
enddo ! ii
|
enddo ! ii
|
||||||
@ -196,7 +207,7 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2 $parame
|
|||||||
deallocate (keys_out,ia_ja_pairs)
|
deallocate (keys_out,ia_ja_pairs)
|
||||||
$omp_end_parallel
|
$omp_end_parallel
|
||||||
$finalization
|
$finalization
|
||||||
|
abort_here = abort_all
|
||||||
end
|
end
|
||||||
|
|
||||||
subroutine $subroutine_monoexc(key_in, hole_1,particl_1 $parameters )
|
subroutine $subroutine_monoexc(key_in, hole_1,particl_1 $parameters )
|
||||||
@ -330,10 +341,14 @@ subroutine $subroutine($params_main)
|
|||||||
generators_bitmask(1,1,d_hole2,i_bitmask_gen), &
|
generators_bitmask(1,1,d_hole2,i_bitmask_gen), &
|
||||||
generators_bitmask(1,1,d_part2,i_bitmask_gen) &
|
generators_bitmask(1,1,d_part2,i_bitmask_gen) &
|
||||||
$params_post)
|
$params_post)
|
||||||
|
if (abort_here) then
|
||||||
|
exit
|
||||||
|
endif
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
$copy_buffer
|
$copy_buffer
|
||||||
$generate_psi_guess
|
$generate_psi_guess
|
||||||
|
abort_here = abort_all
|
||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
|
@ -1,4 +1,4 @@
|
|||||||
BEGIN_PROVIDER [ integer, davidson_iter_max]
|
BEGIN_PROVIDER [ integer, davidson_iter_max ]
|
||||||
implicit none
|
implicit none
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! Max number of Davidson iterations
|
! Max number of Davidson iterations
|
||||||
@ -6,7 +6,7 @@ BEGIN_PROVIDER [ integer, davidson_iter_max]
|
|||||||
davidson_iter_max = 100
|
davidson_iter_max = 100
|
||||||
END_PROVIDER
|
END_PROVIDER
|
||||||
|
|
||||||
BEGIN_PROVIDER [ integer, davidson_sze_max]
|
BEGIN_PROVIDER [ integer, davidson_sze_max ]
|
||||||
implicit none
|
implicit none
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! Max number of Davidson sizes
|
! Max number of Davidson sizes
|
||||||
@ -322,6 +322,7 @@ subroutine davidson_diag_hjj(dets_in,u_in,H_jj,energies,dim_in,sze,N_st,Nint,iun
|
|||||||
!END DEBUG
|
!END DEBUG
|
||||||
|
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
if (.not.converged) then
|
if (.not.converged) then
|
||||||
iter = davidson_sze_max-1
|
iter = davidson_sze_max-1
|
||||||
endif
|
endif
|
||||||
@ -360,6 +361,7 @@ subroutine davidson_diag_hjj(dets_in,u_in,H_jj,energies,dim_in,sze,N_st,Nint,iun
|
|||||||
y, &
|
y, &
|
||||||
lambda &
|
lambda &
|
||||||
)
|
)
|
||||||
|
abort_here = abort_all
|
||||||
end
|
end
|
||||||
|
|
||||||
BEGIN_PROVIDER [ character(64), davidson_criterion ]
|
BEGIN_PROVIDER [ character(64), davidson_criterion ]
|
||||||
@ -397,4 +399,5 @@ subroutine davidson_converged(energy,residual,N_st,converged)
|
|||||||
else if (davidson_criterion == 'both') then
|
else if (davidson_criterion == 'both') then
|
||||||
converged = dabs(maxval(residual(1:N_st))) + dabs(maxval(E(1:N_st)) ) < davidson_threshold
|
converged = dabs(maxval(residual(1:N_st))) + dabs(maxval(E(1:N_st)) ) < davidson_threshold
|
||||||
endif
|
endif
|
||||||
|
converged = converged.or.abort_here
|
||||||
end
|
end
|
||||||
|
@ -19,6 +19,9 @@ program cisd
|
|||||||
print *, 'PT2 = ', pt2
|
print *, 'PT2 = ', pt2
|
||||||
print *, 'E = ', E_old
|
print *, 'E = ', E_old
|
||||||
print *, 'E+PT2 = ', E_old+pt2
|
print *, 'E+PT2 = ', E_old+pt2
|
||||||
|
if (abort_all) then
|
||||||
|
exit
|
||||||
|
endif
|
||||||
enddo
|
enddo
|
||||||
deallocate(pt2,norm_pert)
|
deallocate(pt2,norm_pert)
|
||||||
end
|
end
|
||||||
|
@ -124,53 +124,53 @@ Documentation
|
|||||||
\int_{-infty}^{+infty} (x-A_center)^(power_A) * (x-B_center)^power_B * exp(-alpha(x-A_center)^2) * exp(-beta(x-B_center)^2) dx
|
\int_{-infty}^{+infty} (x-A_center)^(power_A) * (x-B_center)^power_B * exp(-alpha(x-A_center)^2) * exp(-beta(x-B_center)^2) dx
|
||||||
.br
|
.br
|
||||||
|
|
||||||
`align_double <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L65>`_
|
`align_double <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L70>`_
|
||||||
Compute 1st dimension such that it is aligned for vectorization.
|
Compute 1st dimension such that it is aligned for vectorization.
|
||||||
|
|
||||||
`all_utils <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L1>`_
|
`all_utils <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L1>`_
|
||||||
Dummy provider to provide all utils
|
Dummy provider to provide all utils
|
||||||
|
|
||||||
`binom <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L47>`_
|
`binom <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L52>`_
|
||||||
Binomial coefficients
|
Binomial coefficients
|
||||||
|
|
||||||
`binom_func <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L16>`_
|
`binom_func <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L21>`_
|
||||||
.. math ::
|
.. math ::
|
||||||
.br
|
.br
|
||||||
\frac{i!}{j!(i-j)!}
|
\frac{i!}{j!(i-j)!}
|
||||||
.br
|
.br
|
||||||
|
|
||||||
`binom_transp <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L48>`_
|
`binom_transp <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L53>`_
|
||||||
Binomial coefficients
|
Binomial coefficients
|
||||||
|
|
||||||
`dble_fact <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L124>`_
|
`dble_fact <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L129>`_
|
||||||
n!!
|
n!!
|
||||||
|
|
||||||
`fact <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L80>`_
|
`fact <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L85>`_
|
||||||
n!
|
n!
|
||||||
|
|
||||||
`fact_inv <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L112>`_
|
`fact_inv <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L117>`_
|
||||||
1/n!
|
1/n!
|
||||||
|
|
||||||
`inv_int <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L171>`_
|
`inv_int <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L176>`_
|
||||||
1/i
|
1/i
|
||||||
|
|
||||||
`normalize <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L270>`_
|
`normalize <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L275>`_
|
||||||
Normalizes vector u
|
Normalizes vector u
|
||||||
u is expected to be aligned in memory.
|
u is expected to be aligned in memory.
|
||||||
|
|
||||||
`nproc <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L197>`_
|
`nproc <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L202>`_
|
||||||
Number of current OpenMP threads
|
Number of current OpenMP threads
|
||||||
|
|
||||||
`u_dot_u <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L239>`_
|
`u_dot_u <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L244>`_
|
||||||
Compute <u|u>
|
Compute <u|u>
|
||||||
|
|
||||||
`u_dot_v <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L213>`_
|
`u_dot_v <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L218>`_
|
||||||
Compute <u|v>
|
Compute <u|v>
|
||||||
|
|
||||||
`wall_time <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L182>`_
|
`wall_time <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L187>`_
|
||||||
The equivalent of cpu_time, but for the wall time.
|
The equivalent of cpu_time, but for the wall time.
|
||||||
|
|
||||||
`write_git_log <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L157>`_
|
`write_git_log <http://github.com/LCPQ/quantum_package/tree/master/src/Utils/util.irp.f#L162>`_
|
||||||
Write the last git commit in file iunit.
|
Write the last git commit in file iunit.
|
||||||
|
|
||||||
|
|
||||||
|
49
src/Utils/abort.irp.f
Normal file
49
src/Utils/abort.irp.f
Normal file
@ -0,0 +1,49 @@
|
|||||||
|
BEGIN_PROVIDER [ logical, abort_all ]
|
||||||
|
implicit none
|
||||||
|
BEGIN_DOC
|
||||||
|
! If True, all the calculation is aborted
|
||||||
|
END_DOC
|
||||||
|
abort_all = .False.
|
||||||
|
|
||||||
|
END_PROVIDER
|
||||||
|
|
||||||
|
BEGIN_PROVIDER [ logical, abort_here ]
|
||||||
|
implicit none
|
||||||
|
BEGIN_DOC
|
||||||
|
! If True, all the calculation is aborted
|
||||||
|
END_DOC
|
||||||
|
abort_here = abort_all
|
||||||
|
END_PROVIDER
|
||||||
|
|
||||||
|
subroutine trap_signals
|
||||||
|
use ifport
|
||||||
|
implicit none
|
||||||
|
BEGIN_DOC
|
||||||
|
! What to do when a signal is caught. Here, trap Ctrl-C and call the control_C subroutine.
|
||||||
|
END_DOC
|
||||||
|
integer, external :: control_C
|
||||||
|
integer :: err, flag
|
||||||
|
flag = -1
|
||||||
|
err = signal (sigint, control_C, flag)
|
||||||
|
PROVIDE abort_all
|
||||||
|
PROVIDE abort_here
|
||||||
|
end subroutine trap_signals
|
||||||
|
|
||||||
|
integer function control_C(signum)
|
||||||
|
implicit none
|
||||||
|
integer, intent(in) :: signum
|
||||||
|
BEGIN_DOC
|
||||||
|
! What to do on Ctrl-C. If two Ctrl-C are pressed within 1 sec, the calculation if aborted.
|
||||||
|
END_DOC
|
||||||
|
double precision, save :: last_time
|
||||||
|
double precision :: this_time
|
||||||
|
control_C = 0
|
||||||
|
call wall_time(this_time)
|
||||||
|
if (this_time - last_time < 1.d0) then
|
||||||
|
print *, 'Caught Ctrl-C'
|
||||||
|
abort_all = .True.
|
||||||
|
endif
|
||||||
|
last_time = this_time
|
||||||
|
abort_here = .True.
|
||||||
|
end subroutine control_C
|
||||||
|
|
@ -5,11 +5,16 @@ BEGIN_PROVIDER [ logical, all_utils ]
|
|||||||
END_DOC
|
END_DOC
|
||||||
! Do not move this : it greps itself
|
! Do not move this : it greps itself
|
||||||
BEGIN_SHELL [ /bin/bash ]
|
BEGIN_SHELL [ /bin/bash ]
|
||||||
for i in $(grep "BEGIN_PROVIDER" $QPACKAGE_ROOT/src/Utils/*.irp.f | cut -d ',' -f 2 | cut -d ']' -f 1 | tail --lines=+3 )
|
for i in $(grep "BEGIN_PROVIDER" $QPACKAGE_ROOT/src/Utils/*.irp.f \
|
||||||
|
| grep ',' | cut -d ',' -f 2 | cut -d ']' -f 1 | tail --lines=+3 )
|
||||||
do
|
do
|
||||||
|
if [[ ! -z $i ]]
|
||||||
|
then
|
||||||
echo PROVIDE $i
|
echo PROVIDE $i
|
||||||
|
fi
|
||||||
done
|
done
|
||||||
END_SHELL
|
END_SHELL
|
||||||
|
call trap_signals
|
||||||
|
|
||||||
END_PROVIDER
|
END_PROVIDER
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user