9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2025-04-25 17:54:44 +02:00

Introduced qp_bug.irp.f

This commit is contained in:
Anthony Scemama 2025-02-13 10:23:52 +01:00
parent 4a8551be15
commit eb6e1d4339
2 changed files with 39 additions and 6 deletions

View File

@ -197,7 +197,9 @@ integer function get_index_in_psi_det_alpha_unique(key,Nint)
enddo
i += 1
ASSERT (i <= N_det_alpha_unique)
if (i> N_det_alpha_unique) then
call qp_bug(irp_here, i, 'i> N_det_alpha_unique')
endif
!DIR$ FORCEINLINE
do while (spin_det_search_key(psi_det_alpha_unique(1,i),Nint) == det_ref)
@ -219,12 +221,15 @@ integer function get_index_in_psi_det_alpha_unique(key,Nint)
endif
i += 1
if (i > N_det_alpha_unique) then
ASSERT (get_index_in_psi_det_alpha_unique > 0)
return
exit
endif
enddo
if (get_index_in_psi_det_alpha_unique <= 0) then
call qp_bug(irp_here, get_index_in_psi_det_alpha_unique, 'get_index_in_psi_det_alpha_unique <= 0')
endif
end
integer function get_index_in_psi_det_beta_unique(key,Nint)
@ -277,7 +282,9 @@ integer function get_index_in_psi_det_beta_unique(key,Nint)
enddo
i += 1
ASSERT (i <= N_det_beta_unique)
if (i > N_det_beta_unique) then
call qp_bug(irp_here, i, 'i> N_det_beta_unique')
endif
!DIR$ FORCEINLINE
do while (spin_det_search_key(psi_det_beta_unique(1,i),Nint) == det_ref)
@ -299,12 +306,15 @@ integer function get_index_in_psi_det_beta_unique(key,Nint)
endif
i += 1
if (i > N_det_beta_unique) then
ASSERT (get_index_in_psi_det_beta_unique > 0)
return
exit
endif
enddo
if (get_index_in_psi_det_beta_unique <= 0) then
call qp_bug(irp_here, i, 'get_index_in_psi_det_beta_unique <= 0')
endif
end

23
src/utils/bug.irp.f Normal file
View File

@ -0,0 +1,23 @@
subroutine qp_bug(from, code, message)
implicit none
BEGIN_DOC
! This routine prints a bug report
END_DOC
character*(*) :: from
integer :: code
character*(*) :: message
print *, ''
print *, '======================='
print *, 'Bug in Quantum Package!'
print *, '======================='
print *, ''
print *, ' from: ', trim(from)
print *, ' code: ', code
print *, ' info: ', trim(message)
print *, ''
print *, 'Please report this bug at https://github.com/QuantumPackage/qp2/issues'
print *, 'with your output file attached.'
print *, ''
stop -1
end subroutine qp_bug