mirror of
https://github.com/QuantumPackage/qp2.git
synced 2024-12-23 21:03:49 +01:00
modified reorder_core_orb for periodic
This commit is contained in:
parent
cc840cdbc1
commit
0722e12882
@ -30,47 +30,81 @@ subroutine initialize_mo_coef_begin_iteration
|
|||||||
end
|
end
|
||||||
|
|
||||||
subroutine reorder_core_orb
|
subroutine reorder_core_orb
|
||||||
implicit none
|
implicit none
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! TODO: modify for complex
|
! TODO: test for complex
|
||||||
! routines that takes the current :c:data:`mo_coef` and reorder the core orbitals (see :c:data:`list_core` and :c:data:`n_core_orb`) according to the overlap with :c:data:`mo_coef_begin_iteration`
|
! routines that takes the current :c:data:`mo_coef` and reorder the core orbitals (see :c:data:`list_core` and :c:data:`n_core_orb`) according to the overlap with :c:data:`mo_coef_begin_iteration`
|
||||||
END_DOC
|
END_DOC
|
||||||
integer :: i,j,iorb
|
integer :: i,j,iorb
|
||||||
integer :: k,l
|
integer :: k,l
|
||||||
double precision, allocatable :: accu(:)
|
integer, allocatable :: index_core_orb(:),iorder(:)
|
||||||
integer, allocatable :: index_core_orb(:),iorder(:)
|
double precision, allocatable :: accu(:)
|
||||||
double precision, allocatable :: mo_coef_tmp(:,:)
|
integer :: i1,i2
|
||||||
allocate(accu(mo_num),index_core_orb(n_core_orb),iorder(mo_num))
|
if (is_periodic) then
|
||||||
allocate(mo_coef_tmp(ao_num,mo_num))
|
complex*16, allocatable :: accu_c(:)
|
||||||
|
allocate(accu(mo_num),accu_c(mo_num),index_core_orb(n_core_orb),iorder(mo_num))
|
||||||
do i = 1, n_core_orb
|
do i = 1, n_core_orb
|
||||||
iorb = list_core(i)
|
iorb = list_core(i)
|
||||||
do j = 1, mo_num
|
do j = 1, mo_num
|
||||||
accu(j) = 0.d0
|
accu(j) = 0.d0
|
||||||
iorder(j) = j
|
accu_c(j) = (0.d0,0.d0)
|
||||||
do k = 1, ao_num
|
iorder(j) = j
|
||||||
do l = 1, ao_num
|
do k = 1, ao_num
|
||||||
accu(j) += mo_coef_begin_iteration(k,iorb) * mo_coef(l,j) * ao_overlap(k,l)
|
do l = 1, ao_num
|
||||||
|
accu_c(j) += dconjg(mo_coef_begin_iteration_complex(k,iorb)) * &
|
||||||
|
mo_coef_complex(l,j) * ao_overlap_complex(k,l)
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
accu(j) = -cdabs(accu_c(j))
|
||||||
|
enddo
|
||||||
|
call dsort(accu,iorder,mo_num)
|
||||||
|
index_core_orb(i) = iorder(1)
|
||||||
enddo
|
enddo
|
||||||
enddo
|
|
||||||
accu(j) = -dabs(accu(j))
|
|
||||||
enddo
|
|
||||||
call dsort(accu,iorder,mo_num)
|
|
||||||
index_core_orb(i) = iorder(1)
|
|
||||||
enddo
|
|
||||||
|
|
||||||
double precision :: x
|
complex*16 :: x_c
|
||||||
integer :: i1,i2
|
do j = 1, n_core_orb
|
||||||
do j = 1, n_core_orb
|
i1 = list_core(j)
|
||||||
i1 = list_core(j)
|
i2 = index_core_orb(j)
|
||||||
i2 = index_core_orb(j)
|
do i=1,ao_num
|
||||||
do i=1,ao_num
|
x_c = mo_coef_complex(i,i1)
|
||||||
x = mo_coef(i,i1)
|
mo_coef_complex(i,i1) = mo_coef_complex(i,i2)
|
||||||
mo_coef(i,i1) = mo_coef(i,i2)
|
mo_coef_complex(i,i2) = x_c
|
||||||
mo_coef(i,i2) = x
|
enddo
|
||||||
enddo
|
enddo
|
||||||
enddo
|
!call loc_cele_routine
|
||||||
!call loc_cele_routine
|
|
||||||
|
deallocate(accu,accu_c,index_core_orb, iorder)
|
||||||
deallocate(accu,index_core_orb, iorder)
|
else
|
||||||
|
allocate(accu(mo_num),index_core_orb(n_core_orb),iorder(mo_num))
|
||||||
|
|
||||||
|
do i = 1, n_core_orb
|
||||||
|
iorb = list_core(i)
|
||||||
|
do j = 1, mo_num
|
||||||
|
accu(j) = 0.d0
|
||||||
|
iorder(j) = j
|
||||||
|
do k = 1, ao_num
|
||||||
|
do l = 1, ao_num
|
||||||
|
accu(j) += mo_coef_begin_iteration(k,iorb) * mo_coef(l,j) * ao_overlap(k,l)
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
accu(j) = -dabs(accu(j))
|
||||||
|
enddo
|
||||||
|
call dsort(accu,iorder,mo_num)
|
||||||
|
index_core_orb(i) = iorder(1)
|
||||||
|
enddo
|
||||||
|
|
||||||
|
double precision :: x
|
||||||
|
do j = 1, n_core_orb
|
||||||
|
i1 = list_core(j)
|
||||||
|
i2 = index_core_orb(j)
|
||||||
|
do i=1,ao_num
|
||||||
|
x = mo_coef(i,i1)
|
||||||
|
mo_coef(i,i1) = mo_coef(i,i2)
|
||||||
|
mo_coef(i,i2) = x
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
!call loc_cele_routine
|
||||||
|
|
||||||
|
deallocate(accu,index_core_orb, iorder)
|
||||||
|
endif
|
||||||
end
|
end
|
||||||
|
@ -14,9 +14,6 @@ ao_ints
|
|||||||
[ double precision, ao_two_e_integral_schwartz,(ao_num,ao_num) ]
|
[ double precision, ao_two_e_integral_schwartz,(ao_num,ao_num) ]
|
||||||
compute_ao_integrals_jl
|
compute_ao_integrals_jl
|
||||||
|
|
||||||
mo_basis
|
|
||||||
reorder_core_orb: implement for periodic
|
|
||||||
|
|
||||||
scf
|
scf
|
||||||
finish ao_two_e_integral_{alpha,beta}_complex (need reverse index?)
|
finish ao_two_e_integral_{alpha,beta}_complex (need reverse index?)
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user