subroutine solve_ks(this)
implicit none
type(ks_alg), intent(inout) :: this
complex(DP), allocatable :: Q(:,:)
complex(DP) :: ctmp
real(DP) :: rtmp
integer :: i,j,jj
integer :: k,n,m,is
select case(this%inext)
case (0)
k = this%k
n = this%n
m = this%m
if (this%debug) write(*,'("%Krylov-Schur iteration:",I5)')this%iter
rtmp = norm(this%V(:,k+1))
rtmp = 1.0_DP/rtmp
!$OMP PARALLEL DO PRIVATE(is)
do is=1,this%n
this%V(is,k+1) = this%V(is,k+1)*rtmp
enddo
!=========================================================================
! Compute Arnoldi factorization, A \ V_m = V_m+1 H_m
! invert method
!=========================================================================
call new_arnoldi(this%arnoldi,n,m,k)
this%inext = this%inext + 1
this%istat = ESLV_OP_NOP
this%arnoldi%V => this%V
this%arnoldi%H => this%H
return
case (1)
call run_arnoldi(this%arnoldi)
select case(this%arnoldi%istat)
case (ESLV_OP_NOP)
this%inext = 1
this%istat = ESLV_OP_NOP
return
case (ESLV_OP_DO_MATVEC)
this%src_vec => this%arnoldi%src_vec
this%dst_vec => this%arnoldi%dst_vec
this%inext = 1
this%istat = ESLV_OP_DO_MATVEC
if (this%debug) write(*,'("Arnoldi iteration.",I4)')this%arnoldi%j
return
case (ESLV_OP_END_ARNOLDI)
this%inext = 2
this%istat = ESLV_OP_NOP
if (this%debug) write(*,'("Arnoldi iteration ENDED.")')
return
end select
case (2)
m = this%m
if (this%debug) write(*,'("%Schur transformation")')
!=========================================================================
! Compute Schur decomposition of H_m = Z_m * T_m * Z_m'
!
! T : upper triangular (sorted)
! Z : orthogonal Schur basis
!=========================================================================
#ifdef _GGGDEBUG
if (this%icnv >0) then
this%T(1:m,1:m) = this%H(1:m,1:m)
this%Z(1:m,1:m) = Z0
call get_schur(this%icnv,m,this%T,this%eval,this%Z,ESLV_MODE_HERMITE)
do j=1,m
do i=1,m
write(*,'("HHH",2I3,2E14.6)')i,j,this%T(i,j)
enddo
enddo
endif
#endif
this%T(1:m,1:m) = this%H(1:m,1:m)
this%Z(1:m,1:m) = Z0
call get_schur(this%icnv,m,this%T,this%eval,this%Z,this%mode)
#ifdef _GGGDEBUG
if (this%icnv >0) then
do j=1,m
do i=1,m
write(*,'("GGG",2I3,2E14.6)')i,j,this%T(i,j)
enddo
enddo
stop
endif
#endif
!=========================================================================
! Krylov-Schur transformation
! _ _
! H <= H Z
!=========================================================================
this%H(1:m,1:m)=this%T(1:m,1:m)
rtmp = this%H(m+1,m)
do i=1,m
this%H(m+1,i) = rtmp*this%Z(m,i)
enddo
!=========================================================================
! Construct Schur vector, V <= V*Z
!=========================================================================
n = this%n
allocate(Q(n,m))
do i=this%icnv+1,m
!$OMP PARALLEL DO PRIVATE(jj,j)
do jj=1,this%n
Q(jj,i) = this%V(jj,this%icnv+1)*this%Z(this%icnv+1,i)
do j=this%icnv+2,m
Q(jj,i) = Q(jj,i) + this%V(jj,j)*this%Z(j,i)
enddo
enddo
!$OMP END PARALLEL DO
enddo
do i=this%icnv+1,m
!$OMP PARALLEL DO PRIVATE(is)
do is=1,this%n
this%V(is,i)=Q(is,i)
enddo
enddo
deallocate(Q)
this%T(1:m,1:m) = this%H(1:m,1:m)
this%Z(1:m,1:m) = Z0
call get_eig(this%icnv,m,this%T,this%eval,this%evec,this%mode)
this%tag(:) = EV_UNCONV
do i=1,m
ctmp = Z0
do j=1,m
ctmp = ctmp + conjg(this%H(m+1,j))*this%evec(j,i)
enddo
this%res(i) = abs(ctmp)/abs(this%eval(i))
if (this%res(i) < this%tol) then
this%tag(i) = EV_CONVED
if (this%debug) write(*,'("%E",I4,3E24.15," T")')i,this%eval(i),this%res(i)
else
if (this%debug) write(*,'("%E",I4,3E24.15," F")')i,this%eval(i),this%res(i)
endif
enddo
this%icnv=0
do i=1,m
if (this%tag(i) == EV_CONVED) then
this%icnv=this%icnv+1
else
exit
endif
enddo
this%tag(:) = EV_UNCONV
do i=1,this%icnv
this%tag(i) = EV_CONVED
this%H(m+1,i) = Z0
enddo
if (this%icnv >= this%nev) then
this%istat = ESLV_OP_NOP
this%inext = 4
return
endif
if (this%debug) write(*,'("%Purge")')
!=========================================================================
! Purge Schur decomposition
!=========================================================================
k= MIN(this%nev+this%icnv,m-1)
this%k = k
do i=1,k
this%H(k+1,i)=this%H(m+1,i)
!
! converged subspace are decoupled.
!
if (this%tag(i) == EV_CONVED) this%H(k+1,i) = Z0
if (this%debug) then
write(*,'("%NEV_ICNV_K",3I4)')this%nev,this%icnv,k
write(*,'("%T",I3,2E15.7)')i,this%H(k+1,i)
endif
enddo
!$OMP PARALLEL DO PRIVATE(is)
do is=1,this%n
this%V(is,k+1) = this%V(is,m+1)
enddo
this%H( 1:m+1,k+1:m) = Z0
this%H(k+2:m+1,1:k+1) = Z0
this%icnv0 = this%icnv
this%iter = this%iter + 1
this%inext = 0
this%istat = ESLV_OP_NOP
if (this%iter > this%maxiter) then
this%inext = 0
this%istat = ESLV_OP_FAIL
return
endif
return
case (4)
call sort_eig(this%icnv,this%m,this%eval,this%n,this%V)
this%istat = ESLV_OP_CONVERGED
this%inext = 0
return
end select
return
end subroutine