2003-01-20 05:58:50 +08:00
|
|
|
!
|
2007-10-24 03:47:26 +08:00
|
|
|
! Copyright (C) 2001-2007 Quantum-Espresso group
|
2003-01-20 05:58:50 +08:00
|
|
|
! 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 .
|
|
|
|
!
|
2004-06-26 01:25:37 +08:00
|
|
|
#include "f_defs.h"
|
2004-06-12 21:44:18 +08:00
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2004-11-04 21:41:23 +08:00
|
|
|
SUBROUTINE orthoatwfc
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! This routine is meant to orthogonalize all the atomic wfcs. This is
|
|
|
|
! useful when we want to compute the occupation of the atomic orbitals
|
|
|
|
! in order to make lda+U calculations
|
|
|
|
!
|
2004-06-12 21:44:18 +08:00
|
|
|
USE kinds, ONLY : DP
|
2003-11-04 18:53:05 +08:00
|
|
|
USE io_global, ONLY : stdout
|
2006-11-05 11:02:28 +08:00
|
|
|
USE io_files, ONLY : iunat, iunsat, nwordatwfc, iunigk
|
2004-06-12 21:44:18 +08:00
|
|
|
USE ions_base, ONLY : nat
|
|
|
|
USE basis, ONLY : natomwfc
|
2007-01-23 00:38:47 +08:00
|
|
|
USE klist, ONLY : nks, xk, ngk
|
2006-11-29 01:25:00 +08:00
|
|
|
USE ldaU, ONLY : swfcatom, U_projection
|
2004-06-12 21:44:18 +08:00
|
|
|
USE wvfct, ONLY : npwx, npw, igk, gamma_only
|
|
|
|
USE uspp, ONLY : nkb, vkb
|
2005-01-06 00:43:26 +08:00
|
|
|
USE becmod, ONLY : becp, rbecp, becp_nc
|
|
|
|
USE noncollin_module, ONLY : noncolin, npol
|
2004-02-12 20:32:06 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
|
|
|
!
|
2004-11-04 21:41:23 +08:00
|
|
|
INTEGER :: ik, ibnd, info, i, j, k, na, nb, nt, isym, n, ntemp, m, &
|
2005-01-06 00:43:26 +08:00
|
|
|
l, lm, ltot, ntot, ipol
|
2003-01-20 05:58:50 +08:00
|
|
|
! the k point under consideration
|
|
|
|
! counter on bands
|
2005-08-28 22:09:42 +08:00
|
|
|
REAL(DP) :: t0, scnds
|
2003-01-20 05:58:50 +08:00
|
|
|
! cpu time spent
|
2004-11-04 21:41:23 +08:00
|
|
|
LOGICAL :: orthogonalize_wfc
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2005-08-28 22:09:42 +08:00
|
|
|
COMPLEX(DP) :: temp, t (5)
|
|
|
|
COMPLEX(DP) , ALLOCATABLE :: wfcatom (:,:), work (:,:), overlap (:,:)
|
|
|
|
REAL(DP) , ALLOCATABLE :: e (:)
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2003-02-08 00:04:36 +08:00
|
|
|
t0 = scnds ()
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2005-01-06 00:43:26 +08:00
|
|
|
IF (noncolin) THEN
|
2006-11-29 01:25:00 +08:00
|
|
|
ALLOCATE (wfcatom( npwx*npol, natomwfc))
|
2005-01-06 00:43:26 +08:00
|
|
|
ELSE
|
|
|
|
ALLOCATE (wfcatom( npwx, natomwfc))
|
|
|
|
END IF
|
2004-11-04 21:41:23 +08:00
|
|
|
ALLOCATE (overlap( natomwfc , natomwfc))
|
|
|
|
ALLOCATE (work ( natomwfc , natomwfc))
|
|
|
|
ALLOCATE (e ( natomwfc))
|
2004-04-03 00:05:17 +08:00
|
|
|
|
2004-11-04 21:41:23 +08:00
|
|
|
IF (U_projection=="file") THEN
|
2004-04-03 00:05:17 +08:00
|
|
|
WRITE( stdout,*) 'LDA+U Projector read from file '
|
2004-11-04 21:41:23 +08:00
|
|
|
RETURN
|
|
|
|
END IF
|
2004-04-03 00:05:17 +08:00
|
|
|
|
2004-11-04 21:41:23 +08:00
|
|
|
IF (U_projection=="atomic") THEN
|
|
|
|
orthogonalize_wfc = .FALSE.
|
2004-04-03 00:05:17 +08:00
|
|
|
WRITE( stdout,*) 'Atomic wfc used for LDA+U Projector are NOT orthogonalized'
|
2004-11-04 21:41:23 +08:00
|
|
|
ELSE IF (U_projection=="ortho-atomic") THEN
|
|
|
|
orthogonalize_wfc = .TRUE.
|
2004-04-03 00:05:17 +08:00
|
|
|
WRITE( stdout,*) 'Atomic wfc used for LDA+U Projector are orthogonalized'
|
2004-11-04 21:41:23 +08:00
|
|
|
IF (gamma_only) THEN
|
2004-04-03 00:05:17 +08:00
|
|
|
WRITE( stdout,*) 'Gamma-only calculation for this case not implemented'
|
2004-11-04 21:41:23 +08:00
|
|
|
STOP
|
|
|
|
END IF
|
2007-07-24 00:47:12 +08:00
|
|
|
ELSE IF (U_projection=="norm-atomic") THEN
|
|
|
|
orthogonalize_wfc = .TRUE.
|
|
|
|
WRITE( stdout,*) 'Atomic wfc used for LDA+U Projector are normalized but NOT orthogonalized'
|
|
|
|
IF (gamma_only) THEN
|
|
|
|
WRITE( stdout,*) 'Gamma-only calculation for this case not implemented'
|
|
|
|
STOP
|
|
|
|
END IF
|
2004-11-04 21:41:23 +08:00
|
|
|
ELSE
|
|
|
|
WRITE( stdout,*) "U_projection_type =", U_projection
|
|
|
|
CALL errore ("orthoatwfc"," this U_projection_type is not valid",1)
|
|
|
|
END IF
|
2004-04-03 00:05:17 +08:00
|
|
|
|
2004-05-12 05:08:21 +08:00
|
|
|
! Allocate the array becp = <beta|wfcatom>
|
2004-11-04 21:41:23 +08:00
|
|
|
IF ( gamma_only ) THEN
|
|
|
|
ALLOCATE (rbecp (nkb,natomwfc))
|
|
|
|
ELSE
|
2005-01-06 00:43:26 +08:00
|
|
|
IF (noncolin) THEN
|
|
|
|
ALLOCATE ( becp_nc (nkb, npol, natomwfc))
|
|
|
|
ELSE
|
|
|
|
ALLOCATE ( becp (nkb,natomwfc))
|
|
|
|
END IF
|
2004-11-04 21:41:23 +08:00
|
|
|
END IF
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2004-11-04 21:41:23 +08:00
|
|
|
IF (nks > 1) REWIND (iunigk)
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2004-11-04 21:41:23 +08:00
|
|
|
DO ik = 1, nks
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2007-01-23 00:38:47 +08:00
|
|
|
npw = ngk (ik)
|
|
|
|
IF (nks > 1) READ (iunigk) igk
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2003-01-20 05:58:50 +08:00
|
|
|
overlap(:,:) = (0.d0,0.d0)
|
|
|
|
work(:,:) = (0.d0,0.d0)
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2007-11-05 21:39:28 +08:00
|
|
|
CALL atomic_wfc (ik, wfcatom)
|
2006-11-05 11:02:28 +08:00
|
|
|
!
|
|
|
|
! write atomic wfc on unit iunat
|
|
|
|
!
|
2006-11-29 01:25:00 +08:00
|
|
|
CALL davcio (wfcatom, nwordatwfc, iunat, ik, 1)
|
2006-11-05 11:02:28 +08:00
|
|
|
|
2004-11-04 21:41:23 +08:00
|
|
|
CALL init_us_2 (npw, igk, xk (1, ik), vkb)
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2004-11-04 21:41:23 +08:00
|
|
|
IF ( gamma_only ) THEN
|
|
|
|
CALL pw_gemm ('Y', nkb, natomwfc, npw, vkb, npwx, &
|
2004-02-12 20:32:06 +08:00
|
|
|
wfcatom, npwx, rbecp, nkb)
|
2004-11-04 21:41:23 +08:00
|
|
|
ELSE
|
2005-01-06 00:43:26 +08:00
|
|
|
IF (noncolin) THEN
|
2006-11-29 01:25:00 +08:00
|
|
|
CALL ccalbec_nc(nkb,npwx,npw,npol,natomwfc,becp_nc,vkb,wfcatom)
|
2005-01-06 00:43:26 +08:00
|
|
|
ELSE
|
|
|
|
CALL ccalbec (nkb, npwx, npw, natomwfc, becp, vkb, wfcatom)
|
|
|
|
END IF
|
2004-11-04 21:41:23 +08:00
|
|
|
ENDIF
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2005-01-06 00:43:26 +08:00
|
|
|
IF (noncolin) THEN
|
2006-11-29 01:25:00 +08:00
|
|
|
CALL s_psi_nc (npwx, npw, natomwfc, wfcatom, swfcatom)
|
2005-01-06 00:43:26 +08:00
|
|
|
ELSE
|
|
|
|
CALL s_psi (npwx, npw, natomwfc, wfcatom, swfcatom)
|
|
|
|
ENDIF
|
2003-01-20 05:58:50 +08:00
|
|
|
|
2006-04-27 16:40:16 +08:00
|
|
|
IF (orthogonalize_wfc) THEN
|
2006-03-23 02:58:38 +08:00
|
|
|
!
|
|
|
|
! calculate overlap matrix
|
|
|
|
!
|
|
|
|
IF (noncolin) THEN
|
|
|
|
CALL ZGEMM ('c', 'n', natomwfc, natomwfc, npwx*npol, (1.d0, 0.d0), &
|
2006-11-29 01:25:00 +08:00
|
|
|
wfcatom, npwx, swfcatom, npwx, (0.d0,0.d0), overlap, natomwfc)
|
2006-03-23 02:58:38 +08:00
|
|
|
ELSE
|
2006-11-29 01:25:00 +08:00
|
|
|
CALL ZGEMM ('c', 'n', natomwfc, natomwfc, npw, (1.d0, 0.d0), &
|
|
|
|
wfcatom, npwx, swfcatom, npwx, (0.d0, 0.d0), overlap, natomwfc)
|
2006-03-23 02:58:38 +08:00
|
|
|
END IF
|
2003-02-21 22:57:00 +08:00
|
|
|
#ifdef __PARA
|
2006-03-23 02:58:38 +08:00
|
|
|
CALL reduce (2 * natomwfc * natomwfc, overlap)
|
2003-01-20 05:58:50 +08:00
|
|
|
#endif
|
2007-07-24 00:47:12 +08:00
|
|
|
IF (U_projection=="norm-atomic") THEN
|
|
|
|
DO i = 1, natomwfc
|
|
|
|
DO j = i+1, natomwfc
|
|
|
|
overlap(i,j) = cmplx(0.d0,0.d0)
|
|
|
|
overlap(j,i) = cmplx(0.d0,0.d0)
|
|
|
|
ENDDO
|
|
|
|
ENDDO
|
|
|
|
END IF
|
2006-03-23 02:58:38 +08:00
|
|
|
!
|
|
|
|
! find O^-.5
|
|
|
|
!
|
|
|
|
CALL cdiagh (natomwfc, overlap, natomwfc, e, work)
|
|
|
|
DO i = 1, natomwfc
|
|
|
|
e (i) = 1.d0 / dsqrt (e (i) )
|
|
|
|
ENDDO
|
|
|
|
DO i = 1, natomwfc
|
|
|
|
DO j = i, natomwfc
|
|
|
|
temp = (0.d0, 0.d0)
|
|
|
|
DO k = 1, natomwfc
|
|
|
|
temp = temp + e (k) * work (j, k) * CONJG (work (i, k) )
|
|
|
|
ENDDO
|
|
|
|
overlap (i, j) = temp
|
|
|
|
IF (j.NE.i) overlap (j, i) = CONJG (temp)
|
2004-11-04 21:41:23 +08:00
|
|
|
ENDDO
|
2006-03-23 02:58:38 +08:00
|
|
|
ENDDO
|
|
|
|
!
|
|
|
|
! trasform atomic orbitals O^-.5 psi
|
|
|
|
!
|
|
|
|
DO i = 1, npw
|
|
|
|
work(:,1) = (0.d0,0.d0)
|
|
|
|
IF (noncolin) THEN
|
|
|
|
DO ipol=1,npol
|
2006-11-29 01:25:00 +08:00
|
|
|
j = i + (ipol-1)*npwx
|
2006-03-23 02:58:38 +08:00
|
|
|
CALL ZGEMV ('n',natomwfc,natomwfc,(1.d0,0.d0),overlap, &
|
2006-11-29 01:25:00 +08:00
|
|
|
natomwfc,swfcatom(j,1),npwx*npol, &
|
2006-03-23 02:58:38 +08:00
|
|
|
(0.d0,0.d0),work,1)
|
2006-11-29 01:25:00 +08:00
|
|
|
CALL ZCOPY (natomwfc,work,1,swfcatom(j,1),npwx*npol)
|
2006-03-23 02:58:38 +08:00
|
|
|
END DO
|
|
|
|
ELSE
|
|
|
|
CALL ZGEMV ('n', natomwfc, natomwfc, (1.d0, 0.d0) , overlap, &
|
|
|
|
natomwfc, swfcatom (i, 1) , npwx, (0.d0, 0.d0) , work, 1)
|
|
|
|
CALL ZCOPY (natomwfc, work, 1, swfcatom (i, 1), npwx)
|
|
|
|
END IF
|
|
|
|
ENDDO
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2006-04-27 16:40:16 +08:00
|
|
|
END IF ! orthogonalize_wfc
|
|
|
|
|
2006-11-05 11:02:28 +08:00
|
|
|
!
|
2006-11-29 01:25:00 +08:00
|
|
|
! write S * atomic wfc to unit iunsat
|
2006-11-05 11:02:28 +08:00
|
|
|
!
|
2006-11-29 01:25:00 +08:00
|
|
|
CALL davcio (swfcatom, nwordatwfc, iunsat, ik, 1)
|
2004-02-12 20:32:06 +08:00
|
|
|
|
2004-11-04 21:41:23 +08:00
|
|
|
ENDDO
|
|
|
|
DEALLOCATE (overlap)
|
|
|
|
DEALLOCATE (work)
|
|
|
|
DEALLOCATE (e)
|
2006-11-29 01:25:00 +08:00
|
|
|
DEALLOCATE (wfcatom)
|
2004-11-04 21:41:23 +08:00
|
|
|
IF ( gamma_only ) THEN
|
|
|
|
DEALLOCATE (rbecp)
|
|
|
|
ELSE
|
2005-01-06 00:43:26 +08:00
|
|
|
IF (noncolin) THEN
|
|
|
|
DEALLOCATE (becp_nc)
|
|
|
|
ELSE
|
|
|
|
DEALLOCATE (becp)
|
|
|
|
END IF
|
2004-11-04 21:41:23 +08:00
|
|
|
END IF
|
2004-02-12 20:32:06 +08:00
|
|
|
!
|
2004-11-04 21:41:23 +08:00
|
|
|
RETURN
|
2004-02-12 20:32:06 +08:00
|
|
|
|
|
|
|
END SUBROUTINE orthoatwfc
|