2005-03-05 01:46:02 +08:00
|
|
|
!
|
|
|
|
! Copyright (C) 2001 PWSCF group
|
|
|
|
! This file is distributed under the terms of the
|
|
|
|
! GNU General Public License. See the file `License'
|
|
|
|
! in the root directory of the present distribution,
|
|
|
|
! or http://www.gnu.org/copyleft/gpl.txt .
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
subroutine hdiag( max_iter, avg_iter, xk_, et_ )
|
|
|
|
!
|
|
|
|
! Diagonalizes the unperturbed Hamiltonian in a non-selfconsistent way
|
|
|
|
! by Conjugate Gradient (band-by-band)
|
|
|
|
!
|
|
|
|
#include "f_defs.h"
|
2007-11-23 18:03:07 +08:00
|
|
|
USE kinds, ONLY : DP
|
|
|
|
USE cell_base, ONLY: tpiba2
|
|
|
|
USE gvect, ONLY: g, ecfixed, qcutz, q2sigma, gstart
|
|
|
|
USE wvfct, ONLY: g2kin, igk, nbnd, npwx, npw
|
|
|
|
USE uspp, ONLY: vkb, okvan
|
|
|
|
USE noncollin_module, ONLY: npol
|
|
|
|
USE wavefunctions_module,ONLY: evc
|
|
|
|
USE ramanm, ONLY: eth_ns
|
2005-03-05 01:46:02 +08:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
!
|
|
|
|
! I/O variables:
|
|
|
|
!
|
|
|
|
integer :: max_iter
|
|
|
|
! maximum number of iterations
|
2005-08-28 22:09:42 +08:00
|
|
|
real(DP) :: avg_iter, xk_(3), et_(nbnd)
|
2005-03-05 01:46:02 +08:00
|
|
|
! iteration number in the diagonalization
|
|
|
|
! k-point
|
|
|
|
! eigenvalues of the diagonalization
|
|
|
|
!
|
|
|
|
! Local variables:
|
|
|
|
!
|
2007-01-22 18:54:24 +08:00
|
|
|
REAL(DP) :: cg_iter
|
2005-03-05 01:46:02 +08:00
|
|
|
! number of iteration in CG
|
2007-01-22 18:54:24 +08:00
|
|
|
REAL(DP), EXTERNAL :: erf
|
2005-03-05 01:46:02 +08:00
|
|
|
! error function
|
2007-01-22 18:54:24 +08:00
|
|
|
INTEGER :: ig, ntry, notconv
|
2005-03-05 01:46:02 +08:00
|
|
|
! counter on G vectors
|
|
|
|
! number or repeated call to diagonalization in case of non convergence
|
|
|
|
! number of notconverged elements
|
2005-04-15 23:43:33 +08:00
|
|
|
INTEGER, ALLOCATABLE :: btype(:)
|
|
|
|
! type of band: valence (1) or conduction (0)
|
2007-01-22 18:54:24 +08:00
|
|
|
REAl(DP), ALLOCATABLE :: h_prec(:)
|
|
|
|
! preconditioning matrix (diagonal)
|
2005-03-05 01:46:02 +08:00
|
|
|
|
|
|
|
call start_clock ('hdiag')
|
|
|
|
|
2007-01-22 18:54:24 +08:00
|
|
|
allocate (h_prec( npwx), btype(nbnd))
|
2005-04-15 23:43:33 +08:00
|
|
|
btype(:) = 1
|
2005-03-05 01:46:02 +08:00
|
|
|
!
|
|
|
|
! various initializations
|
|
|
|
!
|
|
|
|
call init_us_2 (npw, igk, xk_, vkb)
|
|
|
|
!
|
|
|
|
! sets the kinetic energy
|
|
|
|
!
|
|
|
|
do ig = 1, npw
|
|
|
|
g2kin (ig) =((xk_ (1) + g (1, igk (ig) ) ) **2 + &
|
|
|
|
(xk_ (2) + g (2, igk (ig) ) ) **2 + &
|
|
|
|
(xk_ (3) + g (3, igk (ig) ) ) **2 ) * tpiba2
|
|
|
|
enddo
|
|
|
|
!
|
|
|
|
! if (qcutz > 0.d0) then
|
|
|
|
! do ig = 1, npw
|
|
|
|
! g2kin (ig) = g2kin (ig) + qcutz * (1.d0 + erf ( (g2kin (ig) &
|
|
|
|
! - ecfixed) / q2sigma) )
|
|
|
|
! enddo
|
|
|
|
! endif
|
|
|
|
|
|
|
|
!
|
|
|
|
! Conjugate-Gradient diagonalization
|
|
|
|
!
|
2007-12-11 23:51:36 +08:00
|
|
|
h_prec=1.0_DP
|
2005-03-05 01:46:02 +08:00
|
|
|
do ig = 1, npw
|
2007-01-22 18:54:24 +08:00
|
|
|
h_prec (ig) = max (1.d0, g2kin (ig) )
|
2005-03-05 01:46:02 +08:00
|
|
|
enddo
|
|
|
|
ntry = 0
|
|
|
|
10 continue
|
2005-04-15 23:43:33 +08:00
|
|
|
if (ntry > 0) then
|
2007-11-23 18:03:07 +08:00
|
|
|
call rotate_wfc &
|
|
|
|
( npwx, npw, nbnd, gstart, nbnd, evc, npol, okvan, evc, et_ )
|
2005-03-05 01:46:02 +08:00
|
|
|
avg_iter = avg_iter + 1.d0
|
|
|
|
endif
|
2007-12-11 23:51:36 +08:00
|
|
|
call ccgdiagg (npwx, npw, nbnd, npol, evc, et_, btype, h_prec, eth_ns, &
|
2005-03-05 01:46:02 +08:00
|
|
|
max_iter, .true., notconv, cg_iter)
|
|
|
|
avg_iter = avg_iter + cg_iter
|
|
|
|
ntry = ntry + 1
|
|
|
|
|
|
|
|
if (ntry.le.5.and.notconv.gt.0) goto 10
|
|
|
|
|
2007-01-22 18:54:24 +08:00
|
|
|
deallocate (btype, h_prec)
|
2005-03-05 01:46:02 +08:00
|
|
|
call stop_clock ('hdiag')
|
|
|
|
|
|
|
|
return
|
|
|
|
end subroutine hdiag
|
|
|
|
|