2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-01-13 22:49:17 +08:00
|
|
|
! Copyright (C) 2001-2004 PWSCF 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 .
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
#include "machine.h"
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
!----------------------------------------------------------------------------
|
2004-01-13 22:49:17 +08:00
|
|
|
SUBROUTINE update_pot()
|
2003-12-10 22:57:07 +08:00
|
|
|
!----------------------------------------------------------------------------
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-01-13 22:49:17 +08:00
|
|
|
! ... update potential, use the integer variable order to decide the way
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... order = 0 copy the old potential (nothing is done)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... order = 1 subtract old atomic charge density and sum the new
|
|
|
|
! ... if dynamics is done the routine extrapolates also
|
|
|
|
! ... the difference between the the scf charge and the
|
|
|
|
! ... atomic one,
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... order = 2 extrapolate the wavefunctions:
|
2004-03-29 16:42:37 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... |psi(t+dt)> = 2*|psi(t)> - |psi(t-dt)>
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... order = 3 extrapolate the wavefunctions with the second-order
|
|
|
|
! ... formula:
|
2004-03-29 16:42:37 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... |psi(t+dt)> = |psi(t) +
|
2004-03-29 16:42:37 +08:00
|
|
|
! ... + alpha0*( |psi(t)> - |psi(t-dt)> )
|
|
|
|
! ... + beta0* ( |psi(t-dt)> - |psi(t-2*dt)> )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-29 16:42:37 +08:00
|
|
|
! ... where alpha0 and beta0 are calculated in
|
|
|
|
! ... "find_alpha_and_beta()" so that |tau'-tau(t+dt)| is
|
|
|
|
! ... minimum;
|
|
|
|
! ... tau' and tau(t+dt) are respectively the atomic positions
|
|
|
|
! ... at time t+dt and the extrapolated one:
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-29 16:42:37 +08:00
|
|
|
! ... tau(t+dt) = tau(t) + alpha0*( tau(t) - tau(t-dt) )
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... + beta0*( tau(t-dt) -tau(t-2*dt) )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
!
|
|
|
|
USE control_flags, ONLY : order, history
|
|
|
|
USE io_files, ONLY : prefix, tmp_dir
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... local variables
|
|
|
|
!
|
|
|
|
INTEGER :: rho_order, wfc_order
|
|
|
|
LOGICAL :: exists
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
CALL start_clock( 'update_pot' )
|
|
|
|
!
|
2004-02-26 19:50:36 +08:00
|
|
|
IF ( order == 0 ) THEN
|
|
|
|
!
|
|
|
|
CALL stop_clock( 'update_pot' )
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
2004-03-15 23:25:20 +08:00
|
|
|
END IF
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... determines the maximum effective order of the extrapolation on the
|
|
|
|
! ... basis of the files that are really available
|
|
|
|
!
|
|
|
|
rho_order = MIN( 1, history )
|
|
|
|
!
|
|
|
|
INQUIRE( FILE = TRIM( tmp_dir ) // &
|
|
|
|
& TRIM( prefix ) // '.oldrho', EXIST = exists )
|
|
|
|
!
|
|
|
|
IF ( exists ) THEN
|
|
|
|
!
|
|
|
|
rho_order = MIN( 2, history )
|
|
|
|
!
|
|
|
|
INQUIRE( FILE = TRIM( tmp_dir ) // &
|
|
|
|
& TRIM( prefix ) // '.oldrho2', EXIST = exists )
|
|
|
|
!
|
|
|
|
IF ( exists ) THEN
|
|
|
|
!
|
|
|
|
rho_order = MIN( 3, history )
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
wfc_order = MIN( 1, history, order )
|
|
|
|
!
|
|
|
|
INQUIRE( FILE = TRIM( tmp_dir ) // &
|
|
|
|
& TRIM( prefix ) // '.oldwfc', EXIST = exists )
|
|
|
|
!
|
|
|
|
IF ( exists ) THEN
|
|
|
|
!
|
|
|
|
wfc_order = MIN( 2, history, order )
|
|
|
|
!
|
|
|
|
INQUIRE( FILE = TRIM( tmp_dir ) // &
|
|
|
|
& TRIM( prefix ) // '.oldwfc2', EXIST = exists )
|
|
|
|
!
|
|
|
|
IF ( exists ) THEN
|
|
|
|
!
|
|
|
|
wfc_order = MIN( 3, history, order )
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
END IF
|
|
|
|
!
|
|
|
|
CALL extrapolate_charge( rho_order )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
IF ( order >= 2 ) CALL extrapolate_wfcs( wfc_order )
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL stop_clock( 'update_pot' )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE update_pot
|
|
|
|
!
|
|
|
|
!
|
|
|
|
!----------------------------------------------------------------------------
|
2004-03-24 17:36:50 +08:00
|
|
|
SUBROUTINE extrapolate_charge( rho_order )
|
2003-12-10 22:57:07 +08:00
|
|
|
!----------------------------------------------------------------------------
|
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
USE io_global, ONLY : stdout
|
|
|
|
USE kinds, ONLY : DP
|
2004-04-02 18:30:15 +08:00
|
|
|
USE cell_base, ONLY : omega, bg, alat
|
2004-03-24 17:36:50 +08:00
|
|
|
USE basis, ONLY : nat, tau, ntyp, ityp
|
|
|
|
USE gvect, ONLY : nrxx, ngm, g, gg, gstart, nr1, nr2, nr3, nl, &
|
|
|
|
eigts1, eigts2, eigts3, nrx1, nrx2, nrx3
|
|
|
|
USE lsda_mod, ONLY : lsda, nspin
|
|
|
|
USE scf, ONLY : rho, rho_core, vr
|
|
|
|
USE control_flags, ONLY : alpha0, beta0, imix
|
|
|
|
USE ener, ONLY : ehart, etxc, vtxc
|
|
|
|
USE cellmd, ONLY : lmovecell, omega_old
|
|
|
|
USE vlocal, ONLY : strf
|
|
|
|
USE io_files, ONLY : prefix
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
INTEGER, INTENT(IN) :: rho_order
|
|
|
|
!
|
|
|
|
! ... local variables
|
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
REAL(KIND=DP), ALLOCATABLE :: work(:), work1(:)
|
2004-01-13 22:49:17 +08:00
|
|
|
! work is the difference between charge density and atomic charge
|
|
|
|
! at time t
|
|
|
|
! work1 is the same thing at time t-dt
|
2003-12-10 22:57:07 +08:00
|
|
|
REAL(KIND=DP) :: charge
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
IF ( rho_order == 0 ) RETURN
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
ALLOCATE( work(nrxx) )
|
2004-03-24 17:36:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
work(:) = 0.D0
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... if order = 1 update the potential subtracting to the charge density
|
|
|
|
! ... the "old" atomic charge and summing the new one
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
WRITE( stdout,'(/5X,"NEW-OLD atomic charge density approx. for the potential")' )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
! ... in the lsda case the magnetization will follow rigidly the density
|
|
|
|
! ... keeping fixed the value of zeta = mag / rho_tot.
|
|
|
|
! ... zeta is set here and put in rho(* ??? while rho(*,1) will contain the
|
|
|
|
! ... total valence charge
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
IF ( lsda ) CALL rho2zeta( rho, rho_core, nrxx, nspin, 1 )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... subtract the old atomic charge density
|
|
|
|
!
|
|
|
|
CALL atomic_rho( work, 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
rho(:,1) = rho(:,1) - work(:)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
IF ( lmovecell ) rho(:,1) = rho(:,1) * omega_old
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... extrapolate the difference between the atomic charge a
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... the self-consistent one
|
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
IF ( rho_order == 1 ) THEN
|
|
|
|
!
|
|
|
|
CALL io_pot( + 1, TRIM( prefix )//'.oldrho', rho, 1 )
|
|
|
|
!
|
|
|
|
ELSE IF ( rho_order == 2 ) THEN
|
|
|
|
!
|
|
|
|
! ... oldrho -> work
|
|
|
|
!
|
|
|
|
CALL io_pot( - 1, TRIM( prefix )//'.oldrho', work, 1 )
|
|
|
|
!
|
|
|
|
! ... rho -> oldrho
|
|
|
|
! ... work -> oldrho2
|
|
|
|
!
|
|
|
|
CALL io_pot( + 1, TRIM( prefix )//'.oldrho', rho, 1 )
|
|
|
|
CALL io_pot( + 1, TRIM( prefix )//'.oldrho2', work, 1 )
|
|
|
|
!
|
|
|
|
! ... alpha0 has been calculated in move_ions
|
|
|
|
!
|
|
|
|
rho(:,1) = rho(:,1) + alpha0 * ( rho(:,1) - work(:) )
|
|
|
|
!
|
|
|
|
ELSE IF ( rho_order == 3 ) THEN
|
|
|
|
!
|
|
|
|
ALLOCATE( work1(nrxx) )
|
|
|
|
!
|
|
|
|
work1(:) = 0.D0
|
|
|
|
!
|
|
|
|
! ... oldrho2 -> work1
|
|
|
|
! ... oldrho -> work
|
|
|
|
!
|
|
|
|
CALL io_pot( - 1, TRIM( prefix )//'.oldrho2', work1, 1 )
|
|
|
|
CALL io_pot( - 1, TRIM( prefix )//'.oldrho', work, 1 )
|
|
|
|
!
|
|
|
|
! ... rho -> oldrho
|
|
|
|
! ... work -> oldrho2
|
|
|
|
!
|
|
|
|
CALL io_pot( + 1, TRIM( prefix )//'.oldrho', rho, 1 )
|
|
|
|
CALL io_pot( + 1, TRIM( prefix )//'.oldrho2', work, 1 )
|
|
|
|
!
|
|
|
|
! ... alpha0 and beta0 have been calculated in move_ions
|
|
|
|
!
|
|
|
|
rho(:,1) = rho(:,1) + alpha0 * ( rho(:,1) - work(:) ) + &
|
|
|
|
beta0 * ( work(:) - work1(:) )
|
|
|
|
!
|
|
|
|
DEALLOCATE( work1 )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
IF ( lmovecell ) rho(:,1) = rho(:,1) / omega
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... calculate structure factors for the new positions
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
IF ( lmovecell ) CALL scale_h()
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
CALL struc_fact( nat, tau, ntyp, ityp, ngm, g, bg, nr1, nr2, nr3, &
|
|
|
|
strf, eigts1, eigts2, eigts3 )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... add atomic charges in the new positions
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL atomic_rho( work, 1 )
|
2004-03-24 17:36:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
rho(:,1) = rho(:,1) + work(:)
|
2004-03-24 17:36:50 +08:00
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
CALL set_rhoc()
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... reset up and down charge densities in the LSDA case
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( lsda ) CALL rho2zeta( rho, rho_core, nrxx, nspin, -1 )
|
|
|
|
!
|
|
|
|
CALL v_of_rho( rho, rho_core, nr1, nr2, nr3, nrx1, nrx2, nrx3, &
|
|
|
|
nrxx, nl, ngm, gstart, nspin, g, gg, alat, omega, &
|
|
|
|
ehart, etxc, vtxc, charge, vr )
|
|
|
|
!
|
|
|
|
! ... write potential (and rho) on file
|
|
|
|
!
|
|
|
|
IF ( imix >= 0 ) CALL io_pot( + 1, TRIM( prefix )//'.rho', rho, nspin )
|
|
|
|
!
|
|
|
|
CALL io_pot( + 1, TRIM( prefix )//'.pot', vr, nspin )
|
|
|
|
!
|
|
|
|
DEALLOCATE( work )
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE extrapolate_charge
|
|
|
|
!
|
|
|
|
!
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
2004-03-24 17:36:50 +08:00
|
|
|
SUBROUTINE extrapolate_wfcs( wfc_order )
|
2003-01-20 05:58:50 +08:00
|
|
|
!-----------------------------------------------------------------------
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... This routine extrapolate the wfc's after a "parallel alignment"
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... of the basis of the t-dt and t time steps, according to a recipe
|
|
|
|
! ... by Mead, Rev. Mod. Phys., vol 64, pag. 51 (1992), eqs. 3.20-3.29
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-29 16:42:37 +08:00
|
|
|
#define ONE (1.D0,0.D0)
|
2003-12-10 22:57:07 +08:00
|
|
|
#define ZERO (0.D0,0.D0)
|
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
USE io_global, ONLY : stdout
|
2004-02-10 22:54:54 +08:00
|
|
|
USE kinds, ONLY : DP
|
2004-01-20 20:26:22 +08:00
|
|
|
USE klist, ONLY : nks
|
2004-03-24 17:36:50 +08:00
|
|
|
USE control_flags, ONLY : isolve, alpha0, beta0, order
|
2004-01-20 20:26:22 +08:00
|
|
|
USE wvfct, ONLY : nbnd, npw, npwx, igk
|
|
|
|
USE io_files, ONLY : nwordwfc, iunigk, iunwfc, iunoldwfc, &
|
2004-03-24 17:36:50 +08:00
|
|
|
iunoldwfc2, prefix
|
2004-01-20 20:26:22 +08:00
|
|
|
USE wavefunctions_module, ONLY : evc
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
IMPLICIT NONE
|
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
INTEGER, INTENT(IN) :: wfc_order
|
|
|
|
!
|
2004-01-13 22:49:17 +08:00
|
|
|
! ... local variables
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
INTEGER :: j, i, ik, zero_ew, lwork, info
|
2004-01-13 22:49:17 +08:00
|
|
|
! do-loop variables
|
|
|
|
! counter on k-points
|
2004-02-10 22:54:54 +08:00
|
|
|
! number of zero 'eigenvalues' of the s_m matrix
|
|
|
|
! used by singular value decomposition (ZGESVD)
|
|
|
|
! flag returned by ZGESVD
|
2004-03-24 17:36:50 +08:00
|
|
|
COMPLEX(KIND=DP), ALLOCATABLE :: s_m(:,:), sp_m(:,:), &
|
|
|
|
u_m(:,:), w_m(:,:), work(:)
|
2004-01-13 22:49:17 +08:00
|
|
|
! the overlap matrix s (eq. 3.24)
|
|
|
|
! its dagger
|
2004-02-10 22:54:54 +08:00
|
|
|
! left unitary matrix in the SVD of sp_m
|
|
|
|
! right unitary matrix in the SVD of sp_m
|
|
|
|
! workspace for ZGESVD
|
2003-12-10 22:57:07 +08:00
|
|
|
COMPLEX(KIND=DP), ALLOCATABLE :: evcold(:,:)
|
2004-01-13 22:49:17 +08:00
|
|
|
! wavefunctions at previous iteration
|
2004-02-10 22:54:54 +08:00
|
|
|
REAL(KIND=DP), ALLOCATABLE :: ew(:), rwork(:)
|
|
|
|
! the eigenvalues of s_m
|
|
|
|
! workspace for ZGESVD
|
2004-03-24 17:36:50 +08:00
|
|
|
LOGICAL :: exst
|
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
IF ( wfc_order == 0 ) THEN
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
RETURN
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
ELSE IF ( wfc_order == 1 ) THEN
|
|
|
|
!
|
|
|
|
CALL diropn( iunoldwfc, TRIM( prefix ) // '.oldwfc', nwordwfc, exst )
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
DO ik = 1, nks
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
! ... "now" -> "old"
|
|
|
|
!
|
|
|
|
CALL davcio( evc, nwordwfc, iunwfc, ik, - 1 )
|
|
|
|
CALL davcio( evc, nwordwfc, iunoldwfc, ik, + 1 )
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
END DO
|
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
CLOSE( UNIT = iunoldwfc, STATUS = 'KEEP' )
|
|
|
|
!
|
|
|
|
ELSE IF ( wfc_order == 2 ) THEN
|
|
|
|
!
|
|
|
|
CALL diropn( iunoldwfc, TRIM( prefix ) // '.oldwfc', nwordwfc, exst )
|
2004-03-29 16:42:37 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
IF ( order > 2 ) &
|
|
|
|
CALL diropn( iunoldwfc2, TRIM( prefix ) // '.oldwfc2', nwordwfc, exst )
|
2004-02-10 22:54:54 +08:00
|
|
|
!
|
|
|
|
ALLOCATE( evcold(npwx,nbnd) )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
WRITE( UNIT = stdout, &
|
|
|
|
FMT = '(5X,"Extrapolating wave-functions (first order) ...")' )
|
|
|
|
!
|
|
|
|
lwork = 5 * nbnd
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
ALLOCATE( s_m(nbnd,nbnd), sp_m(nbnd,nbnd), u_m(nbnd,nbnd), &
|
|
|
|
w_m(nbnd,nbnd), work(lwork), ew(nbnd), rwork(lwork) )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
IF ( nks > 1 ) REWIND( iunigk )
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
|
|
|
zero_ew = 0
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
DO ik = 1, nks
|
|
|
|
!
|
|
|
|
IF ( nks > 1 ) READ( iunigk ) npw, igk
|
|
|
|
!
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc, ik, - 1 )
|
2004-03-24 17:36:50 +08:00
|
|
|
CALL davcio( evc, nwordwfc, iunwfc, ik, - 1 )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
|
|
|
! ... construct s_m = <evcold|evc>
|
|
|
|
!
|
|
|
|
CALL ZGEMM( 'C', 'N', nbnd, nbnd, npw, ONE, evcold, npwx, evc, &
|
|
|
|
npwx, ZERO, s_m, nbnd )
|
2004-03-24 17:36:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
CALL reduce( 2 * nbnd * nbnd, s_m )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... construct sp_m
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
DO i = 1, nbnd
|
|
|
|
!
|
|
|
|
sp_m(:,i) = CONJG( s_m (i,:) )
|
|
|
|
!
|
|
|
|
END DO
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... the unitary matrix [sp_m*s_m]^(-1/2)*sp_m (eq. 3.29)
|
|
|
|
! ... by means the singular value decomposition (SVD) of
|
|
|
|
! ... sp_m = u_m * diag(ew) * w_m
|
|
|
|
! ... becomes u_m * w_m
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
CALL ZGESVD( 'A', 'A', nbnd, nbnd, sp_m, nbnd, ew, u_m, nbnd, &
|
2004-03-29 16:42:37 +08:00
|
|
|
w_m, nbnd, work, lwork, rwork, info )
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
|
|
|
! ... check on eigenvalues
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
DO i = 1, nbnd
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
|
|
|
IF ( ew(i) < 0.1D0 ) zero_ew = zero_ew + 1
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
END DO
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... use sp_m to store u_m * w_m
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
CALL ZGEMM( 'N', 'N', nbnd, nbnd, nbnd, ONE, u_m, nbnd, w_m, &
|
|
|
|
nbnd, ZERO, sp_m, nbnd )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... now use evcold as workspace to calculate "aligned" wavefcts:
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... evcold_i = sum_j evc_j*sp_m_ji (eq.3.21)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
CALL ZGEMM( 'N', 'N', npw, nbnd, nbnd, ONE, evc, npwx, sp_m, &
|
2003-12-10 22:57:07 +08:00
|
|
|
nbnd, ZERO, evcold, npwx )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... save on file the aligned wavefcts
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
CALL davcio( evcold, nwordwfc, iunwfc, ik, + 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... re-read from file the wavefcts at (t-dt)
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
CALL davcio( evc, nwordwfc, iunoldwfc, ik, - 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
! ... extrapolate the wfc's,
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
evc = 2.D0 * evcold - evc
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... move the files: "old" -> "old1" and "now" -> "old"
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
IF ( order > 2 ) THEN
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc, ik, - 1 )
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc2, ik, + 1 )
|
2004-02-26 19:50:36 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
END IF
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
CALL davcio( evcold, nwordwfc, iunwfc, ik, - 1 )
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc, ik, + 1 )
|
|
|
|
!
|
|
|
|
! ... save evc on file iunwfc
|
|
|
|
!
|
|
|
|
CALL davcio( evc, nwordwfc, iunwfc, ik, 1 )
|
|
|
|
!
|
|
|
|
END DO
|
|
|
|
!
|
|
|
|
IF ( zero_ew > 0 ) &
|
|
|
|
WRITE( stdout, '(/,5X,"Message from extrapolate_wfcs: ",/, &
|
|
|
|
& 5X,"the matrix <psi(t-dt)|psi(t)> has ",I2, &
|
|
|
|
& " zero eigenvalues")' ) zero_ew
|
|
|
|
!
|
|
|
|
DEALLOCATE( s_m, sp_m, u_m, w_m, work, ew, rwork )
|
|
|
|
!
|
|
|
|
DEALLOCATE( evcold )
|
|
|
|
!
|
|
|
|
CLOSE( UNIT = iunoldwfc, STATUS = 'KEEP' )
|
|
|
|
IF ( order > 2 ) &
|
|
|
|
CLOSE( UNIT = iunoldwfc2, STATUS = 'KEEP' )
|
|
|
|
!
|
|
|
|
ELSE
|
|
|
|
!
|
|
|
|
! ... case : wfc_order = 3
|
|
|
|
!
|
|
|
|
CALL diropn( iunoldwfc, TRIM( prefix ) // '.oldwfc', nwordwfc, exst )
|
|
|
|
CALL diropn( iunoldwfc2, TRIM( prefix ) // '.oldwfc2', nwordwfc, exst )
|
|
|
|
!
|
|
|
|
ALLOCATE( evcold(npwx,nbnd) )
|
|
|
|
!
|
|
|
|
WRITE( UNIT = stdout, &
|
|
|
|
FMT = '(5X,"Extrapolating wave-functions (second order) ...")' )
|
|
|
|
!
|
|
|
|
lwork = 5 * nbnd
|
|
|
|
!
|
|
|
|
ALLOCATE( s_m(nbnd,nbnd), sp_m(nbnd,nbnd), u_m(nbnd,nbnd), &
|
|
|
|
w_m(nbnd,nbnd), work(lwork), ew(nbnd), rwork(lwork) )
|
|
|
|
!
|
|
|
|
IF ( nks > 1 ) REWIND( iunigk )
|
|
|
|
!
|
|
|
|
zero_ew = 0
|
|
|
|
!
|
|
|
|
DO ik = 1, nks
|
|
|
|
!
|
|
|
|
IF ( nks > 1 ) READ( iunigk ) npw, igk
|
|
|
|
!
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc, ik, - 1 )
|
|
|
|
CALL davcio( evc, nwordwfc, iunwfc, ik, - 1 )
|
|
|
|
!
|
|
|
|
! ... construct s_m = <evcold|evc>
|
|
|
|
!
|
|
|
|
CALL ZGEMM( 'C', 'N', nbnd, nbnd, npw, ONE, evcold, npwx, evc, &
|
|
|
|
npwx, ZERO, s_m, nbnd )
|
|
|
|
!
|
|
|
|
CALL reduce( 2 * nbnd * nbnd, s_m )
|
|
|
|
!
|
|
|
|
! ... construct sp_m
|
|
|
|
!
|
|
|
|
DO i = 1, nbnd
|
|
|
|
!
|
|
|
|
sp_m(:,i) = CONJG( s_m (i,:) )
|
|
|
|
!
|
|
|
|
END DO
|
|
|
|
!
|
|
|
|
! ... the unitary matrix [sp_m*s_m]^(-1/2)*sp_m (eq. 3.29)
|
|
|
|
! ... by means the singular value decomposition (SVD) of
|
|
|
|
! ... sp_m = u_m * diag(ew) * w_m
|
|
|
|
! ... becomes u_m * w_m
|
|
|
|
!
|
|
|
|
CALL ZGESVD( 'A', 'A', nbnd, nbnd, sp_m, nbnd, ew, u_m, nbnd, &
|
|
|
|
w_m, nbnd, work, lwork, rwork, info )
|
|
|
|
!
|
|
|
|
! ... check on eigenvalues
|
|
|
|
!
|
|
|
|
DO i = 1, nbnd
|
|
|
|
!
|
|
|
|
IF ( ew(i) < 0.1D0 ) zero_ew = zero_ew + 1
|
|
|
|
!
|
|
|
|
END DO
|
|
|
|
!
|
|
|
|
! ... use sp_m to store u_m * w_m
|
|
|
|
!
|
|
|
|
CALL ZGEMM( 'N', 'N', nbnd, nbnd, nbnd, ONE, u_m, nbnd, w_m, &
|
|
|
|
nbnd, ZERO, sp_m, nbnd )
|
|
|
|
!
|
|
|
|
! ... now use evcold as workspace to calculate "aligned" wavefcts:
|
|
|
|
!
|
|
|
|
! ... evcold_i = sum_j evc_j*sp_m_ji (eq.3.21)
|
|
|
|
!
|
|
|
|
CALL ZGEMM( 'N', 'N', npw, nbnd, nbnd, ONE, evc, npwx, sp_m, &
|
|
|
|
nbnd, ZERO, evcold, npwx )
|
|
|
|
!
|
|
|
|
! ... save on file the aligned wavefcts
|
|
|
|
!
|
|
|
|
CALL davcio( evcold, nwordwfc, iunwfc, ik, + 1 )
|
|
|
|
!
|
|
|
|
! ... re-read from file the wavefcts at (t-dt)
|
|
|
|
!
|
|
|
|
CALL davcio( evc, nwordwfc, iunoldwfc, ik, - 1 )
|
|
|
|
!
|
|
|
|
! ... extrapolate the wfc's,
|
|
|
|
! ... if wfc_order == 3 use the second order extrapolation formula
|
|
|
|
! ... alpha0 and beta0 are calculated in "move_ions"
|
|
|
|
!
|
|
|
|
evc = ( 1 + alpha0 ) * evcold + ( beta0 - alpha0 ) * evc
|
|
|
|
!
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc2, ik, - 1 )
|
|
|
|
!
|
|
|
|
evc = evc - beta0 * evcold
|
|
|
|
!
|
|
|
|
! ... move the files: "old" -> "old1" and "now" -> "old"
|
|
|
|
!
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc, ik, - 1 )
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc2, ik, + 1 )
|
|
|
|
CALL davcio( evcold, nwordwfc, iunwfc, ik, - 1 )
|
|
|
|
CALL davcio( evcold, nwordwfc, iunoldwfc, ik, + 1 )
|
2003-01-20 05:58:50 +08:00
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
! ... save evc on file iunwfc
|
|
|
|
!
|
|
|
|
CALL davcio( evc, nwordwfc, iunwfc, ik, 1 )
|
|
|
|
!
|
|
|
|
END DO
|
|
|
|
!
|
2004-01-20 20:26:22 +08:00
|
|
|
IF ( zero_ew > 0 ) &
|
2004-02-10 22:54:54 +08:00
|
|
|
WRITE( stdout, '(/,5X,"Message from extrapolate_wfcs: ",/, &
|
2004-01-20 20:26:22 +08:00
|
|
|
& 5X,"the matrix <psi(t-dt)|psi(t)> has ",I2, &
|
2004-02-10 22:54:54 +08:00
|
|
|
& " zero eigenvalues")' ) zero_ew
|
2004-01-20 20:26:22 +08:00
|
|
|
!
|
2004-02-10 22:54:54 +08:00
|
|
|
DEALLOCATE( s_m, sp_m, u_m, w_m, work, ew, rwork )
|
|
|
|
!
|
|
|
|
DEALLOCATE( evcold )
|
2003-12-10 22:57:07 +08:00
|
|
|
!
|
2004-03-24 17:36:50 +08:00
|
|
|
CLOSE( UNIT = iunoldwfc, STATUS = 'KEEP' )
|
|
|
|
CLOSE( UNIT = iunoldwfc2, STATUS = 'KEEP' )
|
|
|
|
!
|
2003-12-10 22:57:07 +08:00
|
|
|
END IF
|
|
|
|
!
|
|
|
|
RETURN
|
|
|
|
!
|
|
|
|
END SUBROUTINE extrapolate_wfcs
|