mirror of
https://github.com/LCPQ/quantum_package
synced 2025-01-09 12:44:07 +01:00
Eigenvalues are now printed with the correct sign
This commit is contained in:
parent
47ac033293
commit
61bcd8c103
@ -119,7 +119,7 @@ subroutine damping_SCF
|
|||||||
write(output_hartree_fock,'(A4,X,A16, X, A16, X, A16, X, A4 )'), '====','================','================','================', '===='
|
write(output_hartree_fock,'(A4,X,A16, X, A16, X, A16, X, A4 )'), '====','================','================','================', '===='
|
||||||
write(output_hartree_fock,*)
|
write(output_hartree_fock,*)
|
||||||
|
|
||||||
call mo_as_eigvectors_of_mo_matrix(Fock_matrix_mo,size(Fock_matrix_mo,1),size(Fock_matrix_mo,2),mo_label)
|
call mo_as_eigvectors_of_mo_matrix(Fock_matrix_mo,size(Fock_matrix_mo,1),size(Fock_matrix_mo,2),mo_label,1)
|
||||||
|
|
||||||
call write_double(output_hartree_fock, E_min, 'Hartree-Fock energy')
|
call write_double(output_hartree_fock, E_min, 'Hartree-Fock energy')
|
||||||
call ezfio_set_hartree_fock_energy(E_min)
|
call ezfio_set_hartree_fock_energy(E_min)
|
||||||
|
@ -13,7 +13,7 @@ subroutine huckel_guess
|
|||||||
label = "Guess"
|
label = "Guess"
|
||||||
call mo_as_eigvectors_of_mo_matrix(mo_mono_elec_integral, &
|
call mo_as_eigvectors_of_mo_matrix(mo_mono_elec_integral, &
|
||||||
size(mo_mono_elec_integral,1), &
|
size(mo_mono_elec_integral,1), &
|
||||||
size(mo_mono_elec_integral,2),label)
|
size(mo_mono_elec_integral,2),label,1)
|
||||||
TOUCH mo_coef
|
TOUCH mo_coef
|
||||||
|
|
||||||
c = 0.5d0 * 1.75d0
|
c = 0.5d0 * 1.75d0
|
||||||
|
@ -185,9 +185,9 @@ subroutine set_natural_mos
|
|||||||
allocate(tmp(size(one_body_dm_mo,1),size(one_body_dm_mo,2)))
|
allocate(tmp(size(one_body_dm_mo,1),size(one_body_dm_mo,2)))
|
||||||
|
|
||||||
! Negation to have the occupied MOs first after the diagonalization
|
! Negation to have the occupied MOs first after the diagonalization
|
||||||
tmp = -one_body_dm_mo
|
tmp = one_body_dm_mo
|
||||||
label = "Natural"
|
label = "Natural"
|
||||||
call mo_as_eigvectors_of_mo_matrix(tmp,size(tmp,1),size(tmp,2),label)
|
call mo_as_eigvectors_of_mo_matrix(tmp,size(tmp,1),size(tmp,2),label,-1)
|
||||||
deallocate(tmp)
|
deallocate(tmp)
|
||||||
|
|
||||||
end
|
end
|
||||||
|
@ -10,7 +10,7 @@ program H_CORE_guess
|
|||||||
label = "Guess"
|
label = "Guess"
|
||||||
call mo_as_eigvectors_of_mo_matrix(mo_mono_elec_integral, &
|
call mo_as_eigvectors_of_mo_matrix(mo_mono_elec_integral, &
|
||||||
size(mo_mono_elec_integral,1), &
|
size(mo_mono_elec_integral,1), &
|
||||||
size(mo_mono_elec_integral,2),label)
|
size(mo_mono_elec_integral,2),label,1)
|
||||||
print *, 'save mos'
|
print *, 'save mos'
|
||||||
call save_mos
|
call save_mos
|
||||||
|
|
||||||
|
@ -7,8 +7,8 @@ program guess_mimi
|
|||||||
mo_coef = ao_ortho_lowdin_coef
|
mo_coef = ao_ortho_lowdin_coef
|
||||||
TOUCH mo_coef
|
TOUCH mo_coef
|
||||||
label = "Guess"
|
label = "Guess"
|
||||||
call mo_as_eigvectors_of_mo_matrix(-ao_overlap, &
|
call mo_as_eigvectors_of_mo_matrix(ao_overlap, &
|
||||||
size(ao_overlap,1), &
|
size(ao_overlap,1), &
|
||||||
size(ao_overlap,2),label)
|
size(ao_overlap,2),label,-1)
|
||||||
call save_mos
|
call save_mos
|
||||||
end
|
end
|
||||||
|
@ -9,7 +9,7 @@ subroutine hcore_guess
|
|||||||
label = "Guess"
|
label = "Guess"
|
||||||
call mo_as_eigvectors_of_mo_matrix(mo_mono_elec_integral, &
|
call mo_as_eigvectors_of_mo_matrix(mo_mono_elec_integral, &
|
||||||
size(mo_mono_elec_integral,1), &
|
size(mo_mono_elec_integral,1), &
|
||||||
size(mo_mono_elec_integral,2),label)
|
size(mo_mono_elec_integral,2),label,1)
|
||||||
print *, 'save mos'
|
print *, 'save mos'
|
||||||
call save_mos
|
call save_mos
|
||||||
SOFT_TOUCH mo_coef mo_label
|
SOFT_TOUCH mo_coef mo_label
|
||||||
|
@ -44,13 +44,14 @@ subroutine save_mos_truncated(n)
|
|||||||
|
|
||||||
end
|
end
|
||||||
|
|
||||||
subroutine mo_as_eigvectors_of_mo_matrix(matrix,n,m,label)
|
subroutine mo_as_eigvectors_of_mo_matrix(matrix,n,m,label,sign)
|
||||||
implicit none
|
implicit none
|
||||||
integer,intent(in) :: n,m
|
integer,intent(in) :: n,m, sign
|
||||||
character*(64), intent(in) :: label
|
character*(64), intent(in) :: label
|
||||||
double precision, intent(in) :: matrix(n,m)
|
double precision, intent(in) :: matrix(n,m)
|
||||||
|
|
||||||
double precision, allocatable :: mo_coef_new(:,:), R(:,:),eigvalues(:)
|
integer :: i,j
|
||||||
|
double precision, allocatable :: mo_coef_new(:,:), R(:,:),eigvalues(:), A(:,:)
|
||||||
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: mo_coef_new, R
|
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: mo_coef_new, R
|
||||||
|
|
||||||
call write_time(output_mo_basis)
|
call write_time(output_mo_basis)
|
||||||
@ -58,30 +59,48 @@ subroutine mo_as_eigvectors_of_mo_matrix(matrix,n,m,label)
|
|||||||
print *, irp_here, ': Error : m/= mo_tot_num'
|
print *, irp_here, ': Error : m/= mo_tot_num'
|
||||||
stop 1
|
stop 1
|
||||||
endif
|
endif
|
||||||
allocate(R(n,m),mo_coef_new(ao_num_align,m),eigvalues(m))
|
allocate(A(n,m),R(n,m),mo_coef_new(ao_num_align,m),eigvalues(m))
|
||||||
|
if (sign == -1) then
|
||||||
|
do j=1,m
|
||||||
|
do i=1,n
|
||||||
|
A(i,j) = -matrix(i,j)
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
else
|
||||||
|
do j=1,m
|
||||||
|
do i=1,n
|
||||||
|
A(i,j) = matrix(i,j)
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
endif
|
||||||
mo_coef_new = mo_coef
|
mo_coef_new = mo_coef
|
||||||
|
|
||||||
call lapack_diag(eigvalues,R,matrix,size(matrix,1),size(matrix,2))
|
call lapack_diag(eigvalues,R,A,n,m)
|
||||||
integer :: i
|
|
||||||
write (output_mo_basis,'(A)'), 'MOs are now **'//trim(label)//'**'
|
write (output_mo_basis,'(A)'), 'MOs are now **'//trim(label)//'**'
|
||||||
write (output_mo_basis,'(A)'), ''
|
write (output_mo_basis,'(A)'), ''
|
||||||
write (output_mo_basis,'(A)'), 'Eigenvalues'
|
write (output_mo_basis,'(A)'), 'Eigenvalues'
|
||||||
write (output_mo_basis,'(A)'), '-----------'
|
write (output_mo_basis,'(A)'), '-----------'
|
||||||
write (output_mo_basis,'(A)'), ''
|
write (output_mo_basis,'(A)'), ''
|
||||||
write (output_mo_basis,'(A)'), '======== ================'
|
write (output_mo_basis,'(A)'), '======== ================'
|
||||||
do i = 1, m
|
if (sign == -1) then
|
||||||
write (output_mo_basis,'(I8,X,F16.10)'), i,eigvalues(i)
|
do i=1,m
|
||||||
|
eigvalues(i) = -eigvalues(i)
|
||||||
|
enddo
|
||||||
|
endif
|
||||||
|
do i=1,m
|
||||||
|
write (output_mo_basis,'(I8,X,F16.10)'), i,eigvalues(i)
|
||||||
enddo
|
enddo
|
||||||
write (output_mo_basis,'(A)'), '======== ================'
|
write (output_mo_basis,'(A)'), '======== ================'
|
||||||
write (output_mo_basis,'(A)'), ''
|
write (output_mo_basis,'(A)'), ''
|
||||||
|
|
||||||
call dgemm('N','N',ao_num,m,m,1.d0,mo_coef_new,size(mo_coef_new,1),R,size(R,1),0.d0,mo_coef,size(mo_coef,1))
|
call dgemm('N','N',ao_num,m,m,1.d0,mo_coef_new,size(mo_coef_new,1),R,size(R,1),0.d0,mo_coef,size(mo_coef,1))
|
||||||
deallocate(mo_coef_new,R,eigvalues)
|
deallocate(A,mo_coef_new,R,eigvalues)
|
||||||
call write_time(output_mo_basis)
|
call write_time(output_mo_basis)
|
||||||
|
|
||||||
mo_label = label
|
mo_label = label
|
||||||
SOFT_TOUCH mo_coef mo_label
|
SOFT_TOUCH mo_coef mo_label
|
||||||
end
|
end
|
||||||
|
|
||||||
subroutine mo_as_eigvectors_of_mo_matrix_sort_by_observable(matrix,observable,n,m,label)
|
subroutine mo_as_eigvectors_of_mo_matrix_sort_by_observable(matrix,observable,n,m,label)
|
||||||
implicit none
|
implicit none
|
||||||
integer,intent(in) :: n,m
|
integer,intent(in) :: n,m
|
||||||
|
Loading…
Reference in New Issue
Block a user