quantum-espresso/CPV/runcp.f90

828 lines
27 KiB
Fortran

!
! Copyright (C) 2002-2005 FPMD-CPV groups
! 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 .
!
! AB INITIO COSTANT PRESSURE MOLECULAR DYNAMICS
!=----------------------------------------------------------------------------=!
MODULE runcp_module
!=----------------------------------------------------------------------------=!
IMPLICIT NONE
PRIVATE
SAVE
PUBLIC :: runcp, runcp_force_pairing
PUBLIC :: runcp_uspp, runcp_ncpp
!=----------------------------------------------------------------------------=!
CONTAINS
!=----------------------------------------------------------------------------=!
! ----------------------------------------------
! BEGIN manual
SUBROUTINE runcp( ttprint, tortho, tsde, cm, c0, cp, cdesc, &
vpot, eigr, fi, ekinc, timerd, timeorto, ht, ei, bec, vnosee )
! This subroutine performs a Car-Parrinello or Steepest-Descent step
! on the electronic variables, computing forces on electrons and,
! when required, the eigenvalues of the Hamiltonian
!
! On output "cp" contains the new plave waves coefficients, while
! "cm" and "c0" are not changed
! ----------------------------------------------
! END manual
! ... declare modules
USE kinds
USE mp_global, ONLY: mpime, nproc
USE mp, ONLY: mp_sum
USE electrons_module, ONLY: pmss, eigs, nb_l
USE cp_electronic_mass, ONLY: emass
USE descriptors_module, ONLY: get_local_dims, owner_of, local_index
USE wave_functions, ONLY : rande, cp_kinetic_energy, gram
USE wave_base, ONLY : frice
USE wave_base, ONLY: hpsi
USE cell_module, ONLY: boxdimensions
USE time_step, ONLY: delt
USE forces, ONLY: dforce
USE orthogonalize, ONLY: ortho
USE wave_types, ONLY: wave_descriptor
USE pseudo_projector, ONLY: projector
USE control_flags, ONLY: tnosee
USE control_flags, ONLY: tdamp
USE wave_constrains, ONLY: update_lambda
USE reciprocal_space_mesh, ONLY: gkmask_l
IMPLICIT NONE
! ... declare subroutine arguments
LOGICAL :: ttprint, tortho, tsde
COMPLEX(DP) :: cm(:,:,:,:), c0(:,:,:,:), cp(:,:,:,:)
TYPE (wave_descriptor), INTENT(IN) :: cdesc
COMPLEX(DP) :: eigr(:,:)
REAL(DP), INTENT(IN) :: fi(:,:,:)
REAL(DP), INTENT(IN) :: bec(:,:)
TYPE (boxdimensions), INTENT(IN) :: ht
REAL (DP) :: vpot(:,:,:,:)
REAL(DP) :: ei(:,:,:)
REAL(DP) :: timerd, timeorto
REAL(DP) :: ekinc(:)
REAL(DP), INTENT(IN) :: vnosee
! ... declare other variables
REAL(DP) :: s1, s2, s3, s4
INTEGER :: ik, nx, nb_lx, ierr, nkl, is
COMPLEX(DP), ALLOCATABLE :: cgam(:,:,:)
REAL(DP), ALLOCATABLE :: gam(:,:,:)
REAL(DP), EXTERNAL :: cclock
! ... end of declarations
! ----------------------------------------------
s1 = cclock()
nb_lx = MAX( nb_l(1), nb_l(2) )
nb_lx = MAX( nb_lx, 1 )
IF( cdesc%gamma ) THEN
ALLOCATE( cgam(1,1,1), gam( nb_lx, SIZE( c0, 2 ), cdesc%nspin ), STAT=ierr)
ELSE
ALLOCATE( cgam( nb_lx, SIZE( c0, 2 ), cdesc%nspin ), gam(1,1,1), STAT=ierr)
END IF
IF( ierr /= 0 ) CALL errore(' runcp ', ' allocating gam, prod ', ierr)
ekinc = 0.0d0
timerd = 0.0d0
timeorto = 0.0d0
! Compute electronic forces and move electrons
CALL runcp_ncpp( cm, c0, cp, cdesc, vpot, eigr, fi, bec, vnosee, &
gam, cgam, lambda = ttprint )
! Compute eigenstate
!
IF( ttprint ) THEN
DO is = 1, cdesc%nspin
nx = cdesc%nbt( is )
nkl = cdesc%nkl
DO ik = 1, nkl
CALL eigs( nx, gam(:,:,is), cgam(:,:,is), tortho, fi(:,ik,is), ei(:,ik,is), cdesc%gamma )
END DO
END DO
END IF
s2 = cclock()
timerd = s2 - s1
! Orthogonalize the new wave functions "cp"
IF( tortho ) THEN
CALL ortho(c0, cp, cdesc, pmss, emass)
ELSE
CALL gram(cp, cdesc)
END IF
s3 = cclock()
timeorto = s3 - s2
! Compute fictitious kinetic energy of the electrons at time t
DO is = 1, cdesc%nspin
ekinc(is) = cp_kinetic_energy( is, cp(:,:,:,is), cm(:,:,:,is), cdesc, pmss, delt)
END DO
DEALLOCATE( cgam, gam, STAT=ierr)
IF( ierr /= 0 ) CALL errore(' runcp ', ' deallocating 1 ', ierr)
RETURN
END SUBROUTINE runcp
!=----------------------------------------------------------------------------------=!
! ----------------------------------------------
! BEGIN manual
SUBROUTINE runcp_ncpp( cm, c0, cp, cdesc, &
vpot, eigr, fi, bec, vnosee, gam, cgam, lambda, fromscra, diis, restart )
! This subroutine performs a Car-Parrinello or Steepest-Descent step
! on the electronic variables, computing forces on electrons and,
! when required, the eigenvalues of the Hamiltonian
!
! On output "cp" contains the new plave waves coefficients, while
! "cm" and "c0" are not changed
! ----------------------------------------------
! END manual
! ... declare modules
USE kinds
USE mp_global, ONLY: mpime, nproc
USE mp, ONLY: mp_sum
USE electrons_module, ONLY: pmss
USE cp_electronic_mass, ONLY: emass
USE wave_base, ONLY: frice, wave_steepest, wave_verlet
USE time_step, ONLY: delt
USE forces, ONLY: dforce
USE wave_types, ONLY: wave_descriptor
USE wave_constrains, ONLY: update_lambda
USE control_flags, ONLY: tnosee, tdamp, tsde
USE pseudo_projector, ONLY: projector
USE reciprocal_space_mesh, ONLY: gkmask_l
IMPLICIT NONE
! ... declare subroutine arguments
COMPLEX(DP) :: cm(:,:,:,:), c0(:,:,:,:), cp(:,:,:,:)
COMPLEX(DP) :: cgam(:,:,:)
REAL(DP) :: gam(:,:,:)
TYPE (wave_descriptor), INTENT(IN) :: cdesc
COMPLEX(DP) :: eigr(:,:)
REAL(DP), INTENT(IN) :: fi(:,:,:)
REAL (DP) :: vpot(:,:,:,:)
REAL (DP), INTENT(IN) :: bec(:,:)
REAL(DP), INTENT(IN) :: vnosee
LOGICAL, OPTIONAL, INTENT(IN) :: lambda, fromscra, diis, restart
! ... declare other variables
REAL(DP) :: svar1, svar2, tmpfac, annee
INTEGER :: i, ig, nx, ngw, nb, ierr, is
INTEGER :: iflag
COMPLEX(DP), ALLOCATABLE :: c2(:), c3(:)
REAL(DP), ALLOCATABLE :: svar3(:)
LOGICAL :: tlam, ttsde
! ... end of declarations
! ----------------------------------------------
IF( PRESENT( lambda ) ) THEN
tlam = lambda
ELSE
tlam = .FALSE.
END IF
iflag = 0
IF( PRESENT( fromscra ) ) THEN
IF( fromscra ) iflag = 1
END IF
IF( PRESENT( restart ) ) THEN
IF( restart ) iflag = 2
END IF
ngw = cdesc%ngwl
ALLOCATE( c2(ngw), c3(ngw), svar3(ngw), STAT = ierr )
IF( ierr /= 0 ) CALL errore(' runcp_ncpp ', ' allocating c2, c3, svar3 ', ierr)
! ... determines friction dynamically according to the Nose' dynamics
!
IF( tnosee ) THEN
annee = vnosee * delt * 0.5d0
ELSE IF ( tdamp ) THEN
annee = frice
ELSE
annee = 0.0d0
END IF
tmpfac = 1.d0 / (1.d0 + annee)
IF( iflag == 0 ) THEN
ttsde = tsde
svar1 = 2.d0 * tmpfac
svar2 = 1.d0 - svar1
svar3( 1:ngw ) = delt * delt / pmss( 1:ngw ) * tmpfac
ELSE IF ( iflag == 1 ) THEN
ttsde = .TRUE.
svar1 = 1.d0
svar2 = 2.d0
svar3( 1:ngw ) = delt * delt / pmss( 1:ngw )
ELSE IF ( iflag == 2 ) THEN
ttsde = .FALSE.
svar1 = 1.d0
svar2 = 0.d0
svar3 = delt * delt / pmss( 1:ngw ) * 0.5d0
END IF
! OPEN(unit=17,form='unformatted',status='old')
! READ(17) vpot ! debug
! CLOSE(17)
DO is = 1, cdesc%nspin
nx = cdesc%nbt( is )
IF( nx > SIZE( fi, 1 ) ) &
CALL errore(' runcp ',' inconsistent occupation numbers ', 1)
nb = nx - MOD(nx, 2)
DO i = 1, nb, 2
CALL dforce( i, is, c0(:,:,1,is), cdesc, fi(:,1,is), c2, c3, vpot(:,:,:,is), eigr, bec )
IF( tlam ) THEN
CALL update_lambda( i, gam( :, :,is), c0(:,:,1,is), cdesc, c2 )
CALL update_lambda( i+1, gam( :, :,is), c0(:,:,1,is), cdesc, c3 )
END IF
IF( iflag == 2 ) THEN
c0(:,i,1,is) = cp(:,i,1,is)
c0(:,i+1,1,is) = cp(:,i+1,1,is)
END IF
IF ( ttsde ) THEN
CALL wave_steepest( cp(:,i,1,is), c0(:,i,1,is), svar3, c2 )
CALL wave_steepest( cp(:,i+1,1,is), c0(:,i+1,1,is), svar3, c3 )
ELSE
cp(:,i,1,is) = cm(:,i,1,is)
cp(:,i+1,1,is) = cm(:,i+1,1,is)
CALL wave_verlet( cp(:,i,1,is), c0(:,i,1,is), svar1, svar2, svar3, c2 )
CALL wave_verlet( cp(:,i+1,1,is), c0(:,i+1,1,is), svar1, svar2, svar3, c3 )
END IF
IF( cdesc%gzero ) cp(1,i,1,is) = DBLE( cp(1,i,1,is) )
IF( cdesc%gzero ) cp(1,i+1,1,is)= DBLE( cp(1,i+1,1,is) )
END DO
IF( MOD(nx,2) /= 0) THEN
nb = nx
CALL dforce( nx, is, c0(:,:,1,is), cdesc, fi(:,1,is), c2, c3, vpot(:,:,:,is), eigr, bec )
IF( tlam ) THEN
CALL update_lambda( nb, gam( :, :,is), c0(:,:,1,is), cdesc, c2 )
END IF
IF( iflag == 2 ) THEN
c0(:,nb,1,is) = cp(:,nb,1,is)
END IF
IF ( ttsde ) THEN
CALL wave_steepest( cp(:,nb,1,is), c0(:,nb,1,is), svar3, c2 )
ELSE
cp(:,nb,1,is) = cm(:,nb,1,is)
CALL wave_verlet( cp(:,nb,1,is), c0(:,nb,1,is), svar1, svar2, svar3, c2 )
END IF
IF( cdesc%gzero ) cp(1,nb,1,is) = DBLE( cp(1,nb,1,is) )
END IF
END DO
DEALLOCATE(svar3, c2, c3, STAT=ierr)
IF( ierr /= 0 ) CALL errore(' runcp_ncpp ', ' deallocating 1 ', ierr)
RETURN
END SUBROUTINE runcp_ncpp
!=----------------------------------------------------------------------------------=!
!cdesc is the desciptor for the wf
!eigr==e^ig*r f is the occupation number
!fnl if the factor non local
SUBROUTINE runcp_force_pairing(ttprint, tortho, tsde, cm, c0, cp, cdesc, &
vpot, eigr, fi, ekinc, timerd, timeorto, ht, ei, bec, vnosee)
! same as runcp, except that electrons are paired forcedly
! i.e. this handles a state dependant Hamiltonian for the paired and unpaired electrons
! unpaired is assumed to exist, to be unique, and located in highest index band
! ----------------------------------------------
! END manual
! ... declare modules
USE kinds
USE mp_global, ONLY: mpime, nproc, group
USE mp, ONLY: mp_sum
USE electrons_module, ONLY: pmss, eigs, nb_l, nupdwn, nspin
USE cp_electronic_mass, ONLY: emass
USE descriptors_module, ONLY: get_local_dims, owner_of, local_index
USE wave_functions, ONLY : rande, cp_kinetic_energy, gram
USE wave_base, ONLY: frice, wave_steepest, wave_verlet
USE wave_base, ONLY: hpsi
USE cell_module, ONLY: boxdimensions
USE time_step, ONLY: delt
USE forces, ONLY: dforce
USE orthogonalize, ONLY: ortho
USE wave_types, ONLY: wave_descriptor
USE control_flags, ONLY: tnosee
USE control_flags, ONLY: tdamp
USE constants, ONLY: au
USE io_global, ONLY: ionode
USE wave_constrains, ONLY: update_lambda
USE reciprocal_space_mesh, ONLY: gkmask_l
IMPLICIT NONE
! ... declare subroutine arguments
LOGICAL :: ttprint, tortho, tsde
COMPLEX(DP) :: cm(:,:,:,:), c0(:,:,:,:), cp(:,:,:,:)
TYPE (wave_descriptor), INTENT(IN) :: cdesc
COMPLEX(DP) :: eigr(:,:)
REAL(DP), INTENT(INOUT) :: fi(:,:,:)
TYPE (boxdimensions), INTENT(IN) :: ht
REAL (DP) :: vpot(:,:,:,:)
REAL(DP) :: ei(:,:,:)
REAL(DP), INTENT(IN) :: bec(:,:)
REAL(DP) :: timerd, timeorto
REAL(DP) :: ekinc(:)
REAL(DP), INTENT(IN) :: vnosee
! ... declare other variables
REAL(DP) :: s3, s4
REAL(DP) :: svar1, svar2, tmpfac, annee
INTEGER :: i, ik,ig, nx, ngw, nb, j, nb_g, nb_lx, ierr, nkl, ibl
INTEGER :: ispin_wfc
REAL(DP), ALLOCATABLE :: occup(:,:), occdown(:,:), occsum(:)
REAL(DP) :: intermed, intermed2
COMPLEX(DP) :: intermed3, intermed4
COMPLEX(DP), ALLOCATABLE :: c2(:)
COMPLEX(DP), ALLOCATABLE :: c3(:)
COMPLEX(DP), ALLOCATABLE :: c4(:)
COMPLEX(DP), ALLOCATABLE :: c5(:)
COMPLEX(DP), ALLOCATABLE :: cgam(:,:)
REAL(DP), ALLOCATABLE :: svar3(:)
REAL(DP), ALLOCATABLE :: gam(:,:)
REAL(DP), ALLOCATABLE :: ei_t(:,:,:)
REAL(DP), EXTERNAL :: cclock
! ... end of declarations
! ----------------------------------------------
IF( nspin == 1 ) &
CALL errore(' runcp_forced_pairing ',' inconsistent nspin ', 1)
nkl = cdesc%nkl
IF( nkl /= SIZE( fi, 2 ) ) &
CALL errore(' runcp_forced_pairing ',' inconsistent number of kpoints ', 1)
ngw = cdesc%ngwl
ALLOCATE(c2(ngw), c3(ngw), c4(ngw), c5(ngw), svar3(ngw), STAT=ierr)
IF( ierr /= 0 ) CALL errore(' runcp_forced_pairing ', ' allocating c2, c3, svar3 ', ierr)
! ... determines friction dynamically according to the Nose' dynamics
IF( tnosee ) THEN
annee = vnosee * delt * 0.5d0
ELSE IF ( tdamp ) THEN
annee = frice
ELSE
annee = 0.0d0
END IF
tmpfac = 1.d0 / (1.d0 + annee)
svar1 = 2.d0 * tmpfac
svar2 = 1.d0 - svar1
svar3(1:ngw) = delt * delt / pmss(1:ngw) * tmpfac
ekinc = 0.0d0
timerd = 0.0d0
timeorto = 0.0d0
nx = cdesc%nbt( 1 )
IF( nx /= SIZE( fi, 1 ) ) &
CALL errore(' runcp_forced_pairing ',' inconsistent occupation numbers ', 1)
IF( nupdwn(1) /= (nupdwn(2) + 1) ) &
CALL errore(' runcp_forced_pairing ',' inconsistent spin numbers ', 1)
nb_g = cdesc%nbt( 1 )
nb_lx = MAX( nb_l(1), nb_l(2) )
nb_lx = MAX( nb_lx, 1 )
IF( cdesc%gamma ) THEN
ALLOCATE(cgam(1,1), gam(nb_lx,nb_g), STAT=ierr)
ELSE
ALLOCATE(cgam(nb_lx,nb_g), gam(1,1), STAT=ierr)
END IF
IF( ierr /= 0 ) CALL errore(' runcp_forced_pairing ', ' allocating gam, prod ', ierr)
ALLOCATE( occup(nx,nkl), occdown(nx,nkl), STAT=ierr )
if ( ierr/=0 ) CALL errore(' runcp_forced_pairing ', 'allocating occup, occdown', ierr)
ALLOCATE (ei_t(nx,nkl,2), STAT=ierr)
IF( ierr /= 0 ) CALL errore(' runcp_forced_pairing ', 'allocating iei_t', ierr)
occup = 0.D0
occdown = 0.D0
occup( 1:nupdwn(1), 1:nkl ) = fi( 1:nupdwn(1), 1:nkl, 1 )
occdown( 1:nupdwn(2), 1:nkl ) = fi( 1:nupdwn(2), 1:nkl, 2 )
! ciclo sui punti K
KAPPA: DO ik = 1, nkl
s4 = cclock()
! inizia a trattare lo stato unpaired
! per lo spin_up unpaired
!
intermed = sum ( c0( :, nupdwn(1), ik, 1 ) * CONJG( c0( :, nupdwn(1), ik, 1 ) ) )
! prodotto delle wf relative all'unpaired el
! lavoro sugli n processori e' per quetso che sommo ...
! vengono messi nella variabile temporanea ei_t(:,:,2)
! ei_t(:,:,1) viene utilizzato in seguito solo per il controllo/conto su <psi|H H|psi><psiunp|psiunp>
! questo e' dovuto al fatto che non posso calcolare gli autovalori con eigs a causa della diversa
! occupazione: lo stato unp dovrebbe gia' essere di suo uno stato di KS
! cmq NON LO POSSO RUOTARE per come e' scritta la rho = sum_{i_1,N} |psi_i|**2 + |psi_unp|**2
!
CALL mp_sum( intermed, group)
ei_t(:,ik,2) = intermed ! <Phiunpaired|Phiunpaired>
! l'autoval dello spin up spaiato la mette in ei; memoria temporanea???
!
! da qui e' per la trattazione degli el. paired, come in LSD dato che utilizzo
! vpot (:,:,:,1) e vpot(:,:,:,2) -nota e' def in spazio R(x,y,z)
! mentre dire c0(:,:,:1) o c0(:,:,:,2) e' la stessa cosa
! accoppiamento che segue e' solo per motivi tecnici
! se il numero di bande e' pari o dispari
! indip dal conto sic o dal particolare problema fisico
! per semplicita' ... di conto considero bande pari
! e poi l'ultima come "dissaccoppiata"
! ripeto questo e' un accoppiamento di bande e non di elettroni
! non faccio alcuna distinzione finora fra gli el
nb = nupdwn(1)
DO i = 1, nb, 2
!
! dforce calcola la forza c2 e c3 sulle bande i e i+1 (sono reali => ne fa due alla volta)
! per il vpot (da potential ed e' il potetnziale di KS) in spin up e in down
!
CALL dforce( i, 1, c0(:,:,1,1), cdesc, fi(:,1,1), c2, c3, vpot(:,:,:,1), eigr, bec )
if( i <= nupdwn(2) ) then
CALL dforce( i, 2, c0(:,:,1,1), cdesc, fi(:,1,2), c4, c5, vpot(:,:,:,2), eigr, bec )
else
c4 = 0.0d0
c5 = 0.0d0
end if
!
! accoppia c2 e c3 da vpot (spin 1) stato i e i+1
! c4 c5 2
! per lo stesso stato con spin diverso
! qui calcolo la forza H|psi> ma e' gia' dato il contributo sia up che dwn
! e quindi qui ho occupazione "fi==2" per gli stati paired; mentre rimane "fi==1" per l'unpaired
c2 = occup(i , ik)*c2 + occdown(i , ik)*c4
c3 = occup(i+1, ik)*c3 + occdown(i+1, ik)*c5
!
! se l'unpaired e' nell'ultima banda "pari"
! allora e' lo stato i+1 e andra' in c3 che per def si trovera' ad avere occdwn=0.d0
! combina in c2 la forza degli spin up/down relativa alla banda i
! conbina in c3 la forza " " " " " " " " " i+1
!
IF( ttprint ) then
!
! c2 e' l'array di comb. lin. dH/dpsi stato i con spin_up e dwn
! c3 i+1
! anche solita divisione sul gamma == matrice lambda dei moltiplicatori di Lagrange
! e faccio il prodotto <psi|dH/dpsi> == <psi|H|psi>
!
CALL update_lambda( i, gam( :, :), c0(:,:,ik,1), cdesc, c2 )
CALL update_lambda( i+1, gam( :, :), c0(:,:,ik,1), cdesc, c3 )
intermed = sum ( c2* CONJG(c2) )
intermed2 = sum ( c3* CONJG(c3) )
intermed3 = sum ( c2* CONJG( c0(:,nupdwn(1),ik,1) ) )
intermed4 = sum ( c3* CONJG( c0(:,nupdwn(1),ik,1) ) )
CALL mp_sum ( intermed, group )
CALL mp_sum ( intermed2, group )
CALL mp_sum ( intermed3, group )
CALL mp_sum ( intermed4, group )
ei_t(i ,ik,1) = intermed * ei_t(i ,ik,2) ! <Phi|H H|Phi>*<Phiunpaired|Phiunpaired>
ei_t(i+1,ik,1) = intermed2 * ei_t(i+1,ik,2)
ei_t(i ,ik,2) = abs (intermed3)
ei_t(i+1,ik,2) = abs (intermed4)
END IF ! ttprint
IF ( tsde ) THEN
CALL wave_steepest( cp(:,i,ik,1), c0(:,i,ik,1), svar3, c2 )
CALL wave_steepest( cp(:,i+1,ik,1), c0(:,i+1,ik,1), svar3, c3 )
ELSE
cp(:,i,ik,1) = cm(:,i,ik,1)
cp(:,i+1,ik,1) = cm(:,i+1,ik,1)
CALL wave_verlet( cp(:,i,ik,1), c0(:,i,ik,1), svar1, svar2, svar3, c2 )
CALL wave_verlet( cp(:,i+1,ik,1), c0(:,i+1,ik,1), svar1, svar2, svar3, c3 )
END IF
IF( cdesc%gzero ) cp(1,i,ik,1) = DBLE( cp(1,i,ik,1) )
IF( cdesc%gzero ) cp(1,i+1,ik,1)= DBLE( cp(1,i+1,ik,1) )
END DO ! bande
IF( MOD(nb,2) /= 0) THEN
!
! devo trattare solo l'tulima banda che conterra' di sicuro l'el unpaired
! in c2 ho quindi la forza relativa all'el unpaired
! per questo conservo in ei_t(:,:,2) il suo autovalore
!
CALL dforce( nb, 1, c0(:,:,1,1), cdesc, fi(:,1,1), c2, c3, vpot(:,:,:,1), eigr, bec )
if( nb <= nupdwn(2) ) then
CALL dforce( nb, 2, c0(:,:,1,1), cdesc, fi(:,1,2), c4, c5, vpot(:,:,:,2), eigr, bec )
else
c4 = 0.0d0
c5 = 0.0d0
end if
IF( ttprint ) THEN
CALL update_lambda( nb, gam( :, :), c0(:,:,ik,1), cdesc, c2 )
if ( nupdwn(1) > nupdwn(2) ) then
intermed = sum ( c2 * CONJG(c2) )
intermed3 = sum ( c2 * CONJG( c0(:, nupdwn(1), ik, 1) ) )
CALL mp_sum ( intermed, group )
CALL mp_sum ( intermed3, group )
ei_t(nb,ik,1) = intermed * ei_t(nb,ik,2) ! <Phi|H H|Phi>*<Phiunpaired|Phiunpaired>
ei_t(nb,ik,2) = abs (intermed3)
endif
END IF
IF ( tsde ) THEN
CALL wave_steepest( cp(:,nb,ik,1), c0(:,nb,ik,1), svar3, c2 )
ELSE
cp(:,nb,ik,1) = cm(:,nb,ik,1)
CALL wave_verlet( cp(:,nb,ik,1), c0(:,nb,ik,1), svar1, svar2, svar3, c2 )
END IF
IF( cdesc%gzero ) cp(1,nb,ik,1) = DBLE( cp(1,nb,ik,1) )
END IF
IF( ttprint ) THEN
IF( ionode .AND. ( nupdwn(2) > 0 ) ) THEN
WRITE(6,1006)
WRITE(6,1003) ik, 1
WRITE(6,1004)
WRITE(6,1007) ( ei_t( i, ik, 2 ) * au, i = 1, nupdwn(2) )
END IF
DO i = 1, nupdwn(2)
ei_t(i, ik, 1) = ei_t(i, ik, 2) * ei_t(i,ik,2) / ei_t(i, ik, 1)
END DO
IF( ionode ) THEN
WRITE(6,1002) ik, 1
WRITE(6,1004)
WRITE(6,1007) ( ei_t( i, ik, 1), i = 1, nupdwn(1) )
WRITE(6,1005) ei_t( nb, ik, 2)
END IF
1002 FORMAT(/,3X,'presence in unpaired space (%), kp = ',I3, ' , spin = ',I2,/)
1003 FORMAT(/,3X,'energie cross-terms <Phunpaired| H|Phipaired>, kp = ',I3, ' , spin = ',I2,/)
1004 FORMAT(/,3X,'componente ei_t(i,,1)==<Phi|H H|Phi>*<Phiunpaired|Phiunpaired> su spin up')
1005 FORMAT(/,3X,'eigenvalue (ei_t) unpaired electron: ei_unp = ',F8.4,/)
1006 FORMAT(/,3X,'eigenvalues (ei_t) for states UP == DWN')
1007 FORMAT(/,3X,10F8.4)
ALLOCATE( occsum( SIZE( fi, 1 ) ) )
occsum(:) = occup(:,ik) + occdown(:,ik)
! CALCOLO GLI AUTOVAL di KS per le wf paired
! impongo il vincolo del force_pairing che spin up e dwn siano uguali
! infine metto gli autovalori in ei
if( cdesc%gamma ) then
CALL eigs(nupdwn(2), gam, cgam, tortho, occsum, ei(:,ik,1), cdesc%gamma)
else
CALL eigs(nupdwn(2), gam, cgam, tortho, occsum, ei(:,ik,1), cdesc%gamma)
endif
DEALLOCATE( occsum )
DO i = 1, nupdwn(2)
ei( i, ik, 2 ) = ei( i , ik, 1)
END DO
ei(nupdwn(1), ik, 1) = ei_t(nupdwn(1), ik, 2)
ei(nupdwn(1), ik, 2) = 0.D0
ei_t(nupdwn(1), ik, 2) = 0.D0
ENDIF
END DO KAPPA
s3 = cclock()
timerd = timerd + s3 - s4
IF( tortho ) THEN
CALL ortho( 1, c0(:,:,:,1), cp(:,:,:,1), cdesc, pmss, emass)
ELSE
CALL gram(1, cp(:,:,:,1), cdesc)
END IF
s4 = cclock()
timeorto = timeorto + s4 - s3
! Compute fictitious kinetic energy of the electrons at time t
ekinc(1) = cp_kinetic_energy( 1, cp(:,:,:,1), cm(:,:,:,1), cdesc, pmss, delt)
ekinc(2) = cp_kinetic_energy( 2, cp(:,:,:,1), cm(:,:,:,1), cdesc, pmss, delt)
DEALLOCATE( ei_t, svar3, c2, c3, c4, c5, cgam, gam, occup, occdown, STAT=ierr)
IF( ierr /= 0 ) CALL errore(' runcp_force_pairing ', ' deallocating ', ierr)
RETURN
END SUBROUTINE runcp_force_pairing
!=----------------------------------------------------------------------------=!
SUBROUTINE runcp_uspp( nfi, fccc, ccc, ema0bg, dt2bye, rhos, bec, c0, cm, &
fromscra, restart )
!
use wave_base, only: wave_steepest, wave_verlet
use wave_base, only: frice
use control_flags, only: tnosee, tbuff, lwf, tsde
!use uspp, only : nhsa=> nkb, betae => vkb, rhovan => becsum, deeq
use uspp, only : deeq, betae => vkb
use reciprocal_vectors, only : gstart
use electrons_base, only : n=>nbsp
use wannier_subroutines, only: ef_potential
use efield_module, only: dforce_efield, tefield
use gvecw, only: ngw
!
IMPLICIT NONE
integer, intent(in) :: nfi
real(8) :: fccc, ccc
real(8) :: ema0bg(:), dt2bye
real(8) :: rhos(:,:)
real(8) :: bec(:,:)
complex(8) :: c0(:,:), cm(:,:)
logical, optional, intent(in) :: fromscra
logical, optional, intent(in) :: restart
!
real(8) :: verl1, verl2, verl3
real(8), allocatable:: emadt2(:)
real(8), allocatable:: emaver(:)
complex(8), allocatable:: c2(:), c3(:)
integer :: i
integer :: iflag
logical :: ttsde
iflag = 0
IF( PRESENT( fromscra ) ) THEN
IF( fromscra ) iflag = 1
END IF
IF( PRESENT( restart ) ) THEN
IF( restart ) iflag = 2
END IF
!
!==== set friction ====
!
IF( iflag == 0 ) THEN
if( tnosee ) then
verl1 = 2.0d0 * fccc
verl2 = 1.0d0 - verl1
verl3 = 1.0d0 * fccc
else
verl1=2./(1.+frice)
verl2=1.-verl1
verl3=1./(1.+frice)
end if
ELSE IF( iflag == 1 .OR. iflag == 2 ) THEN
verl1 = 1.0d0
verl2 = 0.0d0
END IF
allocate(c2(ngw))
allocate(c3(ngw))
ALLOCATE( emadt2( ngw ) )
ALLOCATE( emaver( ngw ) )
IF( iflag == 0 ) THEN
emadt2 = dt2bye * ema0bg
emaver = emadt2 * verl3
ttsde = tsde
ELSE IF( iflag == 1 ) THEN
ccc = 0.5d0 * dt2bye
if(tsde) ccc = dt2bye
emadt2 = ccc * ema0bg
emaver = emadt2
ttsde = .TRUE.
ELSE IF( iflag == 2 ) THEN
emadt2 = dt2bye * ema0bg
emaver = emadt2 * 0.5d0
ttsde = .FALSE.
END IF
if( lwf ) then
call ef_potential( nfi, rhos, bec, deeq, betae, c0, cm, emadt2, emaver, verl1, verl2, c2, c3 )
else
do i=1,n,2
call dforce(bec,betae,i,c0(1,i),c0(1,i+1),c2,c3,rhos)
if( tefield ) then
CALL dforce_efield ( bec, i, c0, c2, c3, rhos)
end if
IF( iflag == 2 ) THEN
cm(:,i) = c0(:,i)
cm(:,i+1) = c0(:,i+1)
END IF
if( ttsde ) then
CALL wave_steepest( cm(:, i ), c0(:, i ), emadt2, c2 )
CALL wave_steepest( cm(:, i+1), c0(:, i+1), emadt2, c3 )
else
CALL wave_verlet( cm(:, i ), c0(:, i ), verl1, verl2, emaver, c2 )
CALL wave_verlet( cm(:, i+1), c0(:, i+1), verl1, verl2, emaver, c3 )
endif
if ( gstart == 2 ) then
cm(1, i)=CMPLX(DBLE(cm(1, i)),0.d0)
cm(1,i+1)=CMPLX(DBLE(cm(1,i+1)),0.d0)
end if
end do
end if
IF( iflag == 0 ) THEN
ccc = fccc * dt2bye
END IF
DEALLOCATE( emadt2 )
DEALLOCATE( emaver )
deallocate(c2)
deallocate(c3)
!
!==== end of loop which updates electronic degrees of freedom
!
! buffer for wavefunctions is unit 21
!
if(tbuff) rewind 21
END SUBROUTINE runcp_uspp
!=----------------------------------------------------------------------------=!
END MODULE runcp_module
!=----------------------------------------------------------------------------=!