10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-06-27 07:32:21 +02:00
QuantumPackage/src/utils_cc/org/phase.org

4.7 KiB

program run
  implicit none

  integer :: n(2), degree1, degree2, exc(0:2,2,2)
  integer, allocatable :: list_anni(:,:), list_crea(:,:)
  double precision :: phase1, phase2
  integer :: h1,h2,p1,p2,s1,s2,i,j

  allocate(list_anni(N_int*bit_kind_size,2))
  allocate(list_crea(N_int*bit_kind_size,2))

  do i = 1, N_det-1
    do j = i+1, N_det
      !call print_det(psi_det(1,1,j),N_int)
      call get_excitation(psi_det(1,1,i),psi_det(1,1,j),exc,degree1,phase1,N_int)
      call decode_exc(exc,degree1,h1,p1,h2,p2,s1,s2)
      !print*,'old',degree1,phase1
      !print*,'h1:',h1,'h2:',h2,'s1:',s1,'s2:',s2
      !print*,'p1:',p1,'p2:',p2
      call get_excitation_degree(psi_det(1,1,i),psi_det(1,1,j),degree1,N_int)
      call get_excitation_general(psi_det(1,1,i),psi_det(1,1,j),degree2,n,list_anni,list_crea,phase2,N_int)
      !print*,'new',degree2,phase2
      !print*,'ha:',list_anni(1:n(1),1),'hb',list_anni(1:n(2),2)
      !print*,'pa:',list_crea(1:n(1),1),'pb',list_crea(1:n(2),2)
      !print*,''
      if (degree1 /= degree2) then
        print*,'Error degree:',degree1,degree2
        call abort
      endif
      if (degree1 <= 2 .and. phase1 /= phase2) then
        print*,'Error phase',phase1,phase2
        call abort
      endif
    enddo
  enddo
  
end

phase

subroutine get_phase_general(det1,det2,phase,degree,Nint)
  implicit none

  integer, intent(in)           :: Nint
  integer(bit_kind), intent(in) :: det1(Nint,2), det2(Nint,2)
  double precision, intent(out) :: phase
  integer, intent(out)          :: degree
  integer :: n(2)
  integer, allocatable :: list_anni(:,:), list_crea(:,:)

  allocate(list_anni(N_int*bit_kind_size,2))
  allocate(list_crea(N_int*bit_kind_size,2))

  call get_excitation_general(det1,det2,degree,n,list_anni,list_crea,phase,Nint)
end

Get excitation general

subroutine get_excitation_general(det1,det2,degree,n,list_anni,list_crea,phase,Nint)

  use bitmasks

  implicit none

  integer, intent(in) :: Nint
  integer(bit_kind), intent(in)  :: det1(Nint,2), det2(Nint,2)
  double precision, intent(out)  :: phase
  integer, intent(out)           :: list_crea(Nint*bit_kind_size,2)
  integer, intent(out)           :: list_anni(Nint*bit_kind_size,2)
  integer, intent(out)           :: degree, n(2)
  
  integer, allocatable           :: l1(:,:), l2(:,:) 
  integer(bit_kind), allocatable :: det_crea(:,:), det_anni(:,:)
  integer, allocatable           :: pos_anni(:,:), pos_crea(:,:)

  integer :: n1(2),n2(2),n_crea(2),n_anni(2),i,j,k,d

  allocate(l1(Nint*bit_kind_size,2))
  allocate(l2(Nint*bit_kind_size,2))
  allocate(det_crea(Nint,2),det_anni(Nint,2))

  ! 1      111010
  ! 2      110101
  !
  !not 1-> 000101
  !    2   110101
  !and     000101 -> crea
  !
  !    1   111010
  !not 2-> 001010
  !        001010 -> anni

  do j = 1, 2
    do i = 1, Nint
      det_crea(i,j) = iand(not(det1(i,j)),det2(i,j))
    enddo
  enddo
  
  do j = 1, 2
    do i = 1, Nint
      det_anni(i,j) = iand(det1(i,j),not(det2(i,j)))
    enddo
  enddo

  call bitstring_to_list_ab(det1,l1,n1,Nint)
  call bitstring_to_list_ab(det2,l2,n2,Nint)
  call bitstring_to_list_ab(det_crea,list_crea,n_crea,Nint)
  call bitstring_to_list_ab(det_anni,list_anni,n_anni,Nint)

  do i = 1, 2
    if (n_crea(i) /= n_anni(i)) then
      print*,'Well, it seems we have a problem here...'
      call abort
    endif
  enddo

  !1    11110011001  1 2 3 4 7 8  11
  !pos               1 2 3 4 5 6  7 
  !2    11100101011  1 2 3 6 8 10 11
  !anni 00010010000  4 7
  !pos               4 5
  !crea 00000100010  6 10
  !pos               4 6
  !4 -> 6  pos(4 -> 4)
  !7 -> 10 pos(5 -> 6)

  n = n_anni
  degree = n_anni(1) + n_anni(2)

  allocate(pos_anni(max(n(1),n(2)),2))
  allocate(pos_crea(max(n(1),n(2)),2))
  
  ! Search pos anni
  do j = 1, 2
    k = 1
    do i = 1, n1(j)
       if (k > n_anni(j)) exit
       if (l1(i,j) /= list_anni(k,j)) cycle
       pos_anni(k,j) = i
       k = k + 1
    enddo
  enddo

  ! Search pos crea
  do j = 1, 2
    k = 1
    do i = 1, n2(j)
       if (k > n_crea(j)) exit
       if (l2(i,j) /= list_crea(k,j)) cycle
       pos_crea(k,j) = i
       k = k + 1
    enddo
  enddo

  ! Distance between the ith anni and the ith crea op
  ! By doing so there is no crossing between the different pairs of anni/crea
  ! and the phase is determined by the sum of the distances
  ! -> (-1)^{sum of the distances}
  d = 0
  do j = 1, 2
    do i = 1, n(j)
      d = d + abs(pos_anni(i,j) - pos_crea(i,j))
    enddo
  enddo
  
  phase = dble((-1)**d)

  ! Debug
  !print*,l2(1:n2(1),1)
  !print*,l2(1:n2(2),2)
  !!call print_det(det1,Nint)
  !!call print_det(det2,Nint)
  !print*,phase
  !print*,''
end