2003-09-24 18:03:55 +08:00
|
|
|
!
|
|
|
|
! Copyright (C) 2003 A. Smogunov
|
|
|
|
! 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-11-01 17:26:40 +08:00
|
|
|
! Generalized to spinor wavefunctions and spin-orbit Oct. 2004 (ADC).
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
#include "f_defs.h"
|
2004-11-01 17:26:40 +08:00
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
SUBROUTINE local
|
2003-09-24 18:03:55 +08:00
|
|
|
!
|
|
|
|
! This subroutine computes 2D eigenfunctions and eigenvalues for
|
|
|
|
! the local potential in each slab and performs 2D reduction of
|
|
|
|
! the plane wave basis set.
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
|
|
|
|
USE io_global, ONLY : stdout, ionode
|
2004-11-01 17:26:40 +08:00
|
|
|
USE pwcom
|
|
|
|
USE noncollin_module, ONLY : npol
|
|
|
|
USE io_files
|
|
|
|
USE cond
|
2004-11-04 21:08:25 +08:00
|
|
|
USE mp_global, ONLY : nproc, me_pool, root_pool
|
2004-03-15 23:25:20 +08:00
|
|
|
USE parallel_include
|
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
IMPLICIT NONE
|
2004-03-15 23:25:20 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
INTEGER :: i, il, j, jl, ixy, ig, jg, ipol, igper, k, &
|
2003-09-24 18:03:55 +08:00
|
|
|
ios, index, number, nprob, nteam, nteamnow, &
|
2004-11-01 17:26:40 +08:00
|
|
|
status, info, kin, kfin, is, js
|
2004-11-04 21:08:25 +08:00
|
|
|
INTEGER, ALLOCATABLE :: fftxy(:,:)
|
|
|
|
REAL(kind=DP), PARAMETER :: eps=1.d-6
|
|
|
|
REAL(kind=DP), ALLOCATABLE :: el(:), gp(:)
|
|
|
|
COMPLEX(kind=DP), ALLOCATABLE :: amat(:,:), amat1(:,:), ymat(:,:), &
|
2003-09-24 18:03:55 +08:00
|
|
|
psibase(:,:), psiprob(:,:)
|
2004-11-04 21:08:25 +08:00
|
|
|
COMPLEX(kind=DP),PARAMETER :: one=(1.d0,0.d0), zero=(0.d0,0.d0)
|
|
|
|
COMPLEX(kind=DP) :: aij, xfact, ZDOTC
|
|
|
|
LOGICAL :: exst
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! To divide the slabs between CPU
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL start_clock('local')
|
|
|
|
CALL slabcpu(nrz, nrzp, nkofz, bdl1, bdl2, bdr1, bdr2, z)
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! If all the information is already contained in the file it reads it.
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
IF (lread_loc) THEN
|
|
|
|
CALL seqopn(4,fil_loc,'unformatted',exst)
|
|
|
|
IF(.NOT.exst) CALL errore ('local','fil_loc not found',1)
|
|
|
|
READ(4) n2d
|
2003-09-24 18:03:55 +08:00
|
|
|
! Allocate variables depending on n2d
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL allocate_cond_2
|
|
|
|
READ(4) ((newbg(ig,il), ig=1, ngper*npol), il=1, n2d)
|
2003-11-04 19:30:32 +08:00
|
|
|
! WRITE( stdout,*) 'ngper, n2d = ', ngper, n2d
|
2004-11-04 21:08:25 +08:00
|
|
|
READ(4) (((psiper(ig,il,k),ig=1,n2d),il=1,n2d), &
|
2003-09-24 18:03:55 +08:00
|
|
|
k=1,nrzp)
|
2004-11-04 21:08:25 +08:00
|
|
|
READ(4) ((zkr(il,k),il=1,n2d),k=1,nrzp)
|
2003-09-24 18:03:55 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
CLOSE(unit=4)
|
|
|
|
RETURN
|
|
|
|
ENDIF
|
2003-09-24 18:03:55 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
ALLOCATE( gp( 2 ) )
|
|
|
|
ALLOCATE( el( ngper * npol ) )
|
|
|
|
ALLOCATE( amat( ngper * npol, ngper * npol ) )
|
|
|
|
ALLOCATE( psibase( ngper * npol, ngper * npol ) )
|
|
|
|
ALLOCATE( psiprob( ngper * npol, ngper * npol ) )
|
|
|
|
ALLOCATE( fftxy(-(nrx-1)/2:nrx/2,-(nry-1)/2:nry/2) )
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! To form fftxy correspondence
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
DO i=1, nrx
|
2003-09-24 18:03:55 +08:00
|
|
|
il=i-1
|
2004-11-04 21:08:25 +08:00
|
|
|
IF (il.GT.nrx/2) il=il-nrx
|
|
|
|
DO j=1, nry
|
2003-09-24 18:03:55 +08:00
|
|
|
jl=j-1
|
2004-11-04 21:08:25 +08:00
|
|
|
IF (jl.GT.nry/2) jl=jl-nry
|
2003-09-24 18:03:55 +08:00
|
|
|
fftxy(il,jl)=i+(j-1)*nrx
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
|
|
|
ENDDO
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! to find kin and kfin
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
DO k=1, nrz
|
|
|
|
IF (z(k).LE.bdl1+eps) kin=k
|
|
|
|
IF (z(k).LE.bdr2-eps) kfin=k
|
|
|
|
ENDDO
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! Starting k and number of CPU
|
|
|
|
!
|
|
|
|
nteam=1
|
2004-11-04 21:08:25 +08:00
|
|
|
|
|
|
|
kin = kin + me_pool
|
|
|
|
nteam = nproc
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! set and solve the eigenvalue equation for each slab
|
|
|
|
!
|
|
|
|
n2d=0
|
|
|
|
nprob=0
|
|
|
|
psibase=(0.d0,0.d0)
|
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
DO WHILE(kin.LE.kfin)
|
2004-11-01 17:26:40 +08:00
|
|
|
amat=(0.d0,0.d0)
|
2004-11-04 21:08:25 +08:00
|
|
|
DO ig=1, ngper
|
|
|
|
DO jg=ig, ngper
|
|
|
|
DO ipol=1, 2
|
2003-09-24 18:03:55 +08:00
|
|
|
gp(ipol)=gper(ipol,ig)-gper(ipol,jg)
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
2003-09-24 18:03:55 +08:00
|
|
|
index=number(gp, at, fftxy, nrx, nry)
|
2004-11-04 21:08:25 +08:00
|
|
|
IF (index.GT.0) THEN
|
|
|
|
DO is=1,npol
|
|
|
|
DO js=1,npol
|
2004-11-01 17:26:40 +08:00
|
|
|
amat(ig+(is-1)*ngper,jg+(js-1)*ngper)=vppot(kin,index,is,js)
|
|
|
|
amat(jg+(js-1)*ngper,ig+(is-1)*ngper)= &
|
|
|
|
DCONJG(amat(ig+(is-1)*ngper,jg+(js-1)*ngper))
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
|
|
|
ENDDO
|
|
|
|
ENDIF
|
|
|
|
ENDDO
|
|
|
|
DO is=1,npol
|
2004-11-01 17:26:40 +08:00
|
|
|
amat(ig+(is-1)*ngper,ig+(is-1)*ngper)= &
|
|
|
|
amat(ig+(is-1)*ngper,ig+(is-1)*ngper)+ &
|
|
|
|
(gper(1,ig)**2 + gper(2,ig)**2)*tpiba2
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
|
|
|
ENDDO
|
|
|
|
CALL hev_ab(ngper*npol, amat, ngper*npol, el, psiprob, &
|
2003-09-24 18:03:55 +08:00
|
|
|
-1.d1, eryd+ewind, nprob)
|
2004-11-01 17:26:40 +08:00
|
|
|
|
|
|
|
! do is=1,ngper*npol
|
|
|
|
! write(6,'("------------------------------",i5,f15.7)') is, el(is)
|
|
|
|
! do ig=1,ngper*npol
|
|
|
|
! if (ig.eq.ngper+1) write(6,'("----")')
|
|
|
|
! write(6,'(i5,2f15.7)') ig, psiprob(ig,is)
|
|
|
|
! enddo
|
|
|
|
! enddo
|
|
|
|
! stop
|
|
|
|
|
2003-09-24 18:03:55 +08:00
|
|
|
#ifdef __PARA
|
2004-11-04 21:08:25 +08:00
|
|
|
IF ( me_pool == root_pool ) THEN
|
|
|
|
CALL mpi_send(nprob,1,MPI_INTEGER,0,17, &
|
2003-09-24 18:03:55 +08:00
|
|
|
MPI_COMM_WORLD,info )
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL errore ('n2d reduction','info<>0 in send',info)
|
|
|
|
CALL mpi_send(psiprob,2*ngper*npol*ngper*npol,MPI_REAL8,0,18, &
|
2003-09-24 18:03:55 +08:00
|
|
|
MPI_COMM_WORLD,info )
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL errore ('n2d reduction','info<>0 in send',info)
|
|
|
|
ELSE
|
|
|
|
CALL gramsh(ngper*npol,nprob,1,nprob, &
|
2003-09-24 18:03:55 +08:00
|
|
|
psibase,psiprob,n2d,epsproj)
|
|
|
|
|
|
|
|
nteamnow=kfin-kin+1
|
2004-11-04 21:08:25 +08:00
|
|
|
IF(nteamnow.GT.nteam) nteamnow=nteam
|
2003-09-24 18:03:55 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
DO ig=1, nteamnow-1
|
|
|
|
CALL mpi_recv(nprob,1,MPI_INTEGER, &
|
2003-09-24 18:03:55 +08:00
|
|
|
ig,17,MPI_COMM_WORLD,status,info )
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL errore ('n2d reduction','info<>0 in recv',info)
|
|
|
|
CALL mpi_recv(psiprob,2*ngper*npol*ngper*npol,MPI_REAL8, &
|
2003-09-24 18:03:55 +08:00
|
|
|
ig,18,MPI_COMM_WORLD,status,info )
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL errore ('n2d reduction','info<>0 in recv',info)
|
|
|
|
CALL gramsh(ngper*npol,nprob,1,nprob, &
|
2003-09-24 18:03:55 +08:00
|
|
|
psibase,psiprob,n2d,epsproj)
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
|
|
|
ENDIF
|
2003-09-24 18:03:55 +08:00
|
|
|
#else
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL gramsh(ngper*npol,nprob,1,nprob,psibase,psiprob,n2d,epsproj)
|
2003-09-24 18:03:55 +08:00
|
|
|
#endif
|
|
|
|
kin=kin+nteam
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
#ifdef __PARA
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL mpi_barrier( )
|
|
|
|
CALL mpi_bcast(n2d,1,MPI_INTEGER,0,MPI_COMM_WORLD,info)
|
|
|
|
CALL errore ('reduction','mpi_bcast 1',info)
|
|
|
|
CALL mpi_bcast(psibase,2*ngper*npol*ngper*npol,MPI_REAL8,0, &
|
2003-09-24 18:03:55 +08:00
|
|
|
MPI_COMM_WORLD,info)
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL errore ('reduction','mpi_bcast 1',info)
|
2003-09-24 18:03:55 +08:00
|
|
|
#endif
|
|
|
|
|
|
|
|
!
|
|
|
|
! Allocate variables depending on n2d
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL allocate_cond_2
|
|
|
|
IF (npol.EQ.2) THEN
|
2004-11-01 17:26:40 +08:00
|
|
|
WRITE( stdout,*) 'ngper, ngper*npol, n2d = ', ngper, ngper*npol, n2d
|
2004-11-04 21:08:25 +08:00
|
|
|
ELSE
|
2004-11-01 17:26:40 +08:00
|
|
|
WRITE( stdout,*) 'ngper, n2d = ', ngper, n2d
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDIF
|
2003-09-24 18:03:55 +08:00
|
|
|
!
|
|
|
|
! Construct components of basis vector set on G_per
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL DCOPY(2*ngper*npol*n2d,psibase,1,newbg,1)
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
!
|
|
|
|
! set and solve the eigenvalue equation for each slab
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
ALLOCATE( amat1( n2d, n2d ) )
|
|
|
|
ALLOCATE( ymat( ngper*npol, n2d ) )
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
! for reduced basis set H'_{ab}=e*^i_aH_{ij}e^j_b
|
2004-11-04 21:08:25 +08:00
|
|
|
DO k=1, nrz
|
|
|
|
IF(nkofz(k).NE.0) THEN
|
2003-09-24 18:03:55 +08:00
|
|
|
ymat=(0.d0,0.d0)
|
|
|
|
!
|
|
|
|
! First compute y_{ib}=H_{ij}e_{jb}
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
DO ig=1, ngper
|
|
|
|
DO jg=1, ngper
|
|
|
|
DO ipol=1, 2
|
2004-11-01 17:26:40 +08:00
|
|
|
gp(ipol) = gper(ipol,ig) - gper(ipol,jg)
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
2004-11-01 17:26:40 +08:00
|
|
|
index=number(gp, at, fftxy, nrx, nry)
|
2004-11-04 21:08:25 +08:00
|
|
|
DO is=1,npol
|
|
|
|
DO js=1,npol
|
|
|
|
IF (index.GT.0) THEN
|
2004-11-01 17:26:40 +08:00
|
|
|
aij=vppot(k,index,is,js)
|
2004-11-04 21:08:25 +08:00
|
|
|
ELSE
|
2004-11-01 17:26:40 +08:00
|
|
|
aij=(0.d0,0.d0)
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDIF
|
|
|
|
IF ((ig.EQ.jg).AND.(is.EQ.js)) &
|
2004-11-01 17:26:40 +08:00
|
|
|
aij=aij+(gper(1,ig)**2+ &
|
|
|
|
gper(2,ig)**2)*tpiba2
|
|
|
|
amat(ig+(is-1)*ngper,jg+(js-1)*ngper)= aij
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDDO
|
|
|
|
ENDDO
|
|
|
|
ENDDO
|
|
|
|
ENDDO
|
|
|
|
CALL ZGEMM('n','n',ngper*npol,n2d,ngper*npol,one,amat,ngper*npol, &
|
2004-11-01 17:26:40 +08:00
|
|
|
newbg,ngper*npol,zero,ymat,ngper*npol)
|
2003-09-24 18:03:55 +08:00
|
|
|
!
|
|
|
|
! and construct H'_{ab}=<e_a|y_b>
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
DO il=1, n2d
|
|
|
|
DO jl=il, n2d
|
2004-11-01 17:26:40 +08:00
|
|
|
amat1(il,jl)=ZDOTC(ngper*npol,newbg(1,il),1,ymat(1,jl),1)
|
2004-11-04 21:08:25 +08:00
|
|
|
amat1(jl,il)=CONJG(amat1(il,jl))
|
|
|
|
ENDDO
|
|
|
|
ENDDO
|
2003-09-24 18:03:55 +08:00
|
|
|
!
|
|
|
|
! Solving the eigenvalue problem and construction zk
|
|
|
|
!
|
|
|
|
info=-1
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL hev_ab(n2d, amat1, n2d, zkr(1,nkofz(k)), &
|
2003-09-24 18:03:55 +08:00
|
|
|
psiper(1,1,nkofz(k)), 0.d0, 0.d0, info)
|
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDIF
|
|
|
|
ENDDO
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
#ifdef __PARA
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL mpi_barrier()
|
2003-09-24 18:03:55 +08:00
|
|
|
#endif
|
|
|
|
|
|
|
|
!
|
|
|
|
! saving the 2D data on the file if lwrite_loc=.t.
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
IF (lwrite_loc) THEN
|
|
|
|
IF(fil_loc.EQ.' ') CALL errore ('local','fil_loc no name',1)
|
|
|
|
CALL seqopn(4,fil_loc,'unformatted',exst)
|
|
|
|
WRITE(4) n2d
|
|
|
|
WRITE(4) ((newbg(ig,il), ig=1, ngper*npol), il=1, n2d)
|
2003-09-24 18:03:55 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
WRITE(4) (((psiper(ig,il,k),ig=1,n2d),il=1,n2d), &
|
2003-09-24 18:03:55 +08:00
|
|
|
k=1,nrzp)
|
2004-11-04 21:08:25 +08:00
|
|
|
WRITE(4) ((zkr(il,k),il=1,n2d),k=1,nrzp)
|
|
|
|
CLOSE(unit=4)
|
|
|
|
ENDIF
|
2003-09-24 18:03:55 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
DEALLOCATE(amat)
|
|
|
|
DEALLOCATE(amat1)
|
|
|
|
DEALLOCATE(ymat)
|
|
|
|
DEALLOCATE(gp)
|
|
|
|
DEALLOCATE(psibase)
|
|
|
|
DEALLOCATE(psiprob)
|
|
|
|
DEALLOCATE(el)
|
|
|
|
DEALLOCATE(fftxy)
|
2003-09-24 18:03:55 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
CALL stop_clock('local')
|
|
|
|
RETURN
|
|
|
|
END SUBROUTINE local
|
2003-09-24 18:03:55 +08:00
|
|
|
!-----------------------------------
|
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
FUNCTION number(gp, at, fftxy, nrx, nry)
|
2003-09-24 18:03:55 +08:00
|
|
|
!
|
|
|
|
! This function receives as input the coordinates of 2D g vector
|
|
|
|
! and write on output its fft position.
|
|
|
|
!
|
2004-11-04 21:08:25 +08:00
|
|
|
IMPLICIT NONE
|
|
|
|
INTEGER :: nrx, nry, fftxy(-(nrx-1)/2:nrx/2, -(nry-1)/2:nry/2), &
|
2003-09-24 18:03:55 +08:00
|
|
|
number, n1, n2
|
2004-11-04 21:08:25 +08:00
|
|
|
REAL(kind=KIND(0.d0)) :: gp(2), at(3,3), x1, x2
|
|
|
|
REAL(kind=KIND(0.d0)), PARAMETER :: eps=1.d-4
|
2003-09-24 18:03:55 +08:00
|
|
|
|
|
|
|
x1=gp(1)*at(1,1)+gp(2)*at(2,1)
|
|
|
|
x2=gp(1)*at(1,2)+gp(2)*at(2,2)
|
2004-11-04 21:08:25 +08:00
|
|
|
n1=NINT(x1)
|
|
|
|
n2=NINT(x2)
|
|
|
|
IF (n1.LE.nrx/2.AND.n1.GE.-(nrx-1)/2.AND. &
|
|
|
|
n2.LE.nry/2.AND.n2.GE.-(nry-1)/2) THEN
|
2003-09-24 18:03:55 +08:00
|
|
|
number=fftxy(n1,n2)
|
2004-11-04 21:08:25 +08:00
|
|
|
ELSE
|
2003-09-24 18:03:55 +08:00
|
|
|
!
|
|
|
|
! The g vector is not inside the 2D mesh
|
|
|
|
!
|
|
|
|
number=-1
|
2004-11-04 21:08:25 +08:00
|
|
|
ENDIF
|
2003-09-24 18:03:55 +08:00
|
|
|
|
2004-11-04 21:08:25 +08:00
|
|
|
RETURN
|
|
|
|
END FUNCTION number
|
2003-09-24 18:03:55 +08:00
|
|
|
|