calculation='bands' introduced (but not finished)

leftover stuff from calculation='raman' removed


git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@2932 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
giannozz 2006-03-20 22:44:35 +00:00
parent 58c4445f85
commit ca10d9a79f
8 changed files with 54 additions and 196 deletions

View File

@ -202,7 +202,7 @@ MODULE control_flags
lsmd, &! if .TRUE. the calc. is string dynamics
lwf, &! if .TRUE. the calc. is with wannier functions
lphonon, &! if .TRUE. the calc. is phonon
lraman, &! if .TRUE. the calc. is raman
lbands, &! if .TRUE. the calc. is band structure
lconstrain, &! if .TRUE. the calc. is constraint
ldamped, &! if .TRUE. the calc. is a damped dynamics
lrescale_t, &! if .TRUE. the ionic temperature is rescaled

View File

@ -94,10 +94,11 @@ MODULE input_parameters
! title = 'title for this simulation'
CHARACTER(LEN=80) :: calculation = 'none'
! calculation = 'scf', 'relax', 'md', 'cp'*, 'vc-relax', 'vc-md',
! 'vc-cp', 'neb', 'smd', 'cp-wf', 'fpmd', 'metadyn'
! Specify the type of the simulation
! 'scf' = electron minimization
! 'scf' = electron minimization/selfconsistency
! 'nscf' = nonselfconsistent calculation of electron states
! 'bands' = as above, band structure calculation only
! 'phonon' = as above, plus data needed for a phonon calculation
! 'relax' = ionic minimization
! 'cp' = Car-Parrinello MD
! 'md' = Car-Parrinello MD
@ -113,7 +114,7 @@ MODULE input_parameters
CHARACTER(LEN=80) :: calculation_allowed(15)
DATA calculation_allowed / 'scf', 'nscf', 'relax', 'md', 'cp', &
'vc-relax', 'vc-md', 'vc-cp', 'phonon', 'raman', 'neb', 'smd', &
'vc-relax', 'vc-md', 'vc-cp', 'phonon', 'bands', 'neb', 'smd', &
'cp-wf', 'fpmd', 'metadyn' /
! Allowed value for calculation parameters
@ -1311,17 +1312,6 @@ MODULE input_parameters
NAMELIST / phonon / modenum, xqq, nq1, nq2, nq3, tr2_ph
!=----------------------------------------------------------------------------=!
! RAMAN Namelist Input Parameters
!=----------------------------------------------------------------------------=!
REAL(DP) :: b_length
! length of the b-vector
LOGICAL :: lcart
! cartesian directions
NAMELIST / raman / b_length, lcart
!=----------------------------------------------------------------------------=!
! WANNIER Namelist Input Parameters
!=----------------------------------------------------------------------------=!

View File

@ -508,26 +508,6 @@ MODULE read_namelists_module
!
!=----------------------------------------------------------------------=!
!
! Variables initialization for Namelist RAMAN
!
!----------------------------------------------------------------------
SUBROUTINE raman_defaults( prog )
!----------------------------------------------------------------------
!
IMPLICIT NONE
!
CHARACTER(LEN=2) :: prog ! ... specify the calling program
!
!
b_length = 0.d0
lcart = .false.
!
RETURN
!
END SUBROUTINE
!
!=----------------------------------------------------------------------=!
!
! Variables initialization for Namelist WANNIER
!
!----------------------------------------------------------------------
@ -980,30 +960,6 @@ MODULE read_namelists_module
!
!=----------------------------------------------------------------------------=!
!
! Broadcast variables values for Namelist RAMAN
!
!=----------------------------------------------------------------------------=!
!
!----------------------------------------------------------------------
SUBROUTINE raman_bcast()
!----------------------------------------------------------------------
!
USE io_global, ONLY: ionode_id
USE mp, ONLY: mp_bcast
!
IMPLICIT NONE
!
!
CALL mp_bcast( b_length, ionode_id )
CALL mp_bcast( lcart, ionode_id )
!
RETURN
!
END SUBROUTINE
!
!
!=----------------------------------------------------------------------------=!
!
! Broadcast variables values for Namelist WANNIER
!
!=----------------------------------------------------------------------------=!
@ -1516,25 +1472,6 @@ MODULE read_namelists_module
!
!=----------------------------------------------------------------------=!
!
! Check input values for Namelist RAMAN
!
!=----------------------------------------------------------------------=!
!
!----------------------------------------------------------------------
SUBROUTINE raman_checkin( prog )
!--------------------------------------------------------------------
!
IMPLICIT NONE
!
CHARACTER(LEN=2) :: prog ! ... specify the calling program
!
!
RETURN
!
END SUBROUTINE
!
!=----------------------------------------------------------------------=!
!
! Check input values for Namelist WANNIER
!
!=----------------------------------------------------------------------=!
@ -1582,10 +1519,7 @@ MODULE read_namelists_module
ion_dynamics = 'none'
cell_dynamics = 'none'
END IF
CASE ('nscf')
! IF( prog == 'CP' ) &
! CALL errore( sub_name,' calculation '//calculation// &
! & ' not implemented ',1)
CASE ('nscf', 'bands')
IF( prog == 'CP' ) occupations = 'bogus'
IF( prog == 'CP' ) electron_dynamics = 'damp'
IF( prog == 'PW' ) startingpot = 'file'
@ -1595,10 +1529,8 @@ MODULE read_namelists_module
& ' not implemented ',1)
IF( prog == 'PW' ) startingpot = 'file'
CASE ('raman')
IF( prog == 'CP' ) &
CALL errore( sub_name,' calculation '//TRIM(calculation)// &
& ' not implemented ',1)
IF( prog == 'PW' ) startingpot = 'file'
CALL errore( sub_name,' calculation '//TRIM(calculation)// &
& ' no longer implemented ',1)
CASE ( 'cp-wf' )
IF( prog == 'CP' ) THEN
electron_dynamics = 'damp'
@ -1694,7 +1626,6 @@ MODULE read_namelists_module
END IF
!
IF ( calculation == 'nscf' .OR. &
calculation == 'raman'.OR. &
calculation == 'phonon' ) THEN
!
startingpot = 'file'
@ -1883,24 +1814,6 @@ MODULE read_namelists_module
CALL phonon_bcast()
CALL phonon_checkin( prog )
!
! ... RAMAN NAMELIST
!
CALL raman_defaults( prog )
ios = 0
IF( ionode ) THEN
IF( TRIM( calculation ) == 'raman' ) THEN
READ( 5, raman, iostat = ios )
END IF
END IF
CALL mp_bcast( ios, ionode_id )
IF( ios /= 0 ) THEN
CALL errore( ' read_namelists ', &
& ' reading namelist raman ', ABS(ios) )
END IF
!
CALL raman_bcast()
CALL raman_checkin( prog )
!
! ... WANNIER NAMELIST
!
CALL wannier_defaults( prog )

View File

@ -42,7 +42,7 @@ SUBROUTINE electrons()
USE scf, ONLY : rho, vr, vltot, vrs, rho_core
USE control_flags, ONLY : mixing_beta, tr2, ethr, niter, nmix, &
iprint, istep, lscf, lmd, conv_elec, &
restart, reduce_io, iverbosity
lbands, restart, reduce_io, iverbosity
USE io_files, ONLY : prefix, iunwfc, iunocc, nwordwfc, iunpath, &
output_drho, iunefield
USE ldaU, ONLY : ns, nsnew, eth, Hubbard_U, &
@ -807,13 +807,13 @@ SUBROUTINE electrons()
!
END DO
!
IF ( lgauss ) THEN
IF ( lgauss .AND. .NOT. lbands ) THEN
!
ef = efermig( et, nbnd, nks, nelec, wk, degauss, ngauss, 0, isk )
!
WRITE( stdout, 9040 ) ef * rytoev
!
ELSE IF ( ltetra ) THEN
ELSE IF ( ltetra .AND. .NOT. lbands ) THEN
!
ef = efermit( et, nbnd, nks, nelec, nspin, ntetra, tetra, 0, isk )
!

View File

@ -99,8 +99,6 @@ SUBROUTINE iosys()
nelec_ => nelec, &
nelup_ => nelup, &
neldw_ => neldw, &
b_length_ => b_length, &
lcart_ => lcart, &
tot_charge_ => tot_charge, &
tot_magnetization_ => tot_magnetization, &
multiplicity_ => multiplicity
@ -138,7 +136,7 @@ SUBROUTINE iosys()
nosym_ => nosym, &
modenum_ => modenum, &
reduce_io, langevin_rescaling, ethr, lscf, lbfgs, &
lmd, lpath, lneb, lsmd, lphonon, ldamped, lraman, &
lmd, lpath, lneb, lsmd, lphonon, ldamped, lbands, &
lrescale_t, lmetadyn, lconstrain, lcoarsegrained, &
restart, twfcollect
!
@ -243,10 +241,6 @@ SUBROUTINE iosys()
!
USE input_parameters, ONLY : phonon, modenum, xqq
!
! ... RAMAN namelist
!
USE input_parameters, ONLY : b_length, lcart
!
! ... "path" specific
!
USE input_parameters, ONLY : pos, full_phs_path_flag
@ -592,7 +586,7 @@ SUBROUTINE iosys()
lsmd = .FALSE.
lmovecell = .FALSE.
lphonon = .FALSE.
lraman = .FALSE.
lbands = .FALSE.
lbfgs = .FALSE.
ldamped = .FALSE.
lforce = tprnfor
@ -609,6 +603,12 @@ SUBROUTINE iosys()
lforce = .FALSE.
nstep = 1
!
CASE( 'bands' )
!
lforce = .FALSE.
lbands = .TRUE.
nstep = 1
!
CASE( 'phonon' )
!
lforce = .FALSE.
@ -742,11 +742,6 @@ SUBROUTINE iosys()
& ': ion_dynamics=' // TRIM( ion_dynamics ) // &
& ' not supported', 1 )
!
CASE( 'raman' )
!
lraman = .TRUE.
nstep = 1
!
CASE( 'neb' )
!
lscf = .TRUE.
@ -1156,9 +1151,6 @@ SUBROUTINE iosys()
modenum_ = modenum
xqq_ = xqq
!
b_length_ = b_length
lcart_ = lcart
!
! ... general "path" variables
!
ds_ = ds

View File

@ -16,7 +16,7 @@ SUBROUTINE punch( what )
USE klist, ONLY : nks, nkstot
USE lsda_mod, ONLY : nspin
USE scf, ONLY : rho
USE control_flags, ONLY : reduce_io, lscf
USE control_flags, ONLY : reduce_io, lscf, lbands
USE wvfct, ONLY : et, wg, nbnd
USE wavefunctions_module, ONLY : evc, evc_nc
USE io_files, ONLY : prefix, iunpun, iunwfc, nwordwfc
@ -58,7 +58,7 @@ SUBROUTINE punch( what )
! ... with few k-points is followed by a non-self-consistent one with added
! ... k-points, whose weight is set to zero.
!
IF ( .NOT. lscf ) CALL sum_band()
IF ( .NOT. lscf .AND. .NOT. lbands ) CALL sum_band()
!
! ... Write: general variables (including dimensions of the arrays),
! ... atomic positions, forces, k-points, eigenvalues

View File

@ -131,8 +131,7 @@ MODULE klist
nelec, &! number of electrons
nelup, &! number of spin-up electrons (if two_fermi_energies=t)
neldw, &! number of spin-dw electrons (if two_fermi_energies=t)
tot_charge, &! total charge
b_length ! length of the b vectors
tot_charge
INTEGER :: &
ngk(npk), &! number of plane waves for each k point
nks, &! number of k points in this pool
@ -143,7 +142,6 @@ MODULE klist
LOGICAL :: &
lgauss, &! if .TRUE.: use gaussian broadening
lxkcry, &! if .TRUE.:k-pnts in cryst. basis accepted in input
lcart, &! if .TRUE.: b vectors in cartesian coordinates
two_fermi_energies ! if .TRUE.: nelup and neldw set ef_up and ef_dw
! separately
!

View File

@ -49,7 +49,7 @@ SUBROUTINE setup()
USE gvect, ONLY : gcutm, ecutwfc, dual, nr1, nr2, nr3
USE gsmooth, ONLY : doublegrid, gcutms
USE klist, ONLY : xk, wk, xqq, nks, nelec, degauss, lgauss, &
lxkcry, nkstot, b_length, lcart, &
lxkcry, nkstot, &
nelup, neldw, two_fermi_energies, &
tot_charge, tot_magnetization, multiplicity
USE lsda_mod, ONLY : lsda, nspin, current_spin, isk, &
@ -63,7 +63,7 @@ SUBROUTINE setup()
USE wvfct, ONLY : nbnd, nbndx, gamma_only
USE control_flags, ONLY : tr2, ethr, alpha0, beta0, lscf, lmd, lpath, &
lphonon, david, isolve, niter, noinv, nosym, &
modenum, lraman
modenum, lbands
USE relax, ONLY : starting_diag_threshold
USE cellmd, ONLY : calc
USE uspp_param, ONLY : psd, betar, nbeta, dion, jjj, lll, tvanp
@ -345,7 +345,7 @@ SUBROUTINE setup()
!
ltest = ( ethr == 0.D0 )
!
IF ( lphonon .or. lraman ) THEN
IF ( lphonon ) THEN
!
! ... in the case of a phonon calculation ethr can not be specified
! ... in the input file
@ -689,56 +689,41 @@ SUBROUTINE setup()
! ... Input k-points are assumed to be given in the IBZ of the Bravais
! ... lattice, with the full point symmetry of the lattice.
! ... If some symmetries are missing in the crystal, "sgama" computes
! ... the missing k-points. If nosym is true (see above) we do not use
! ... any point-group symmetry and leave k-points unchanged.
! ... the missing k-points. If nosym is true (see above) or if calculating
! ... bands, we do not use any point-group symmetry and leave k-points unchanged.
!
input_nks = nks
!
CALL sgama( nrot, nat, s, sname, t_rev, at, bg, tau, ityp, nsym, nr1,&
nr2, nr3, irt, ftau, npk, nks, xk, wk, invsym, minus_q, &
xqq, modenum, noncolin, m_loc )
!
CALL checkallsym( nsym, s, nat, tau, ityp, at, &
bg, nr1, nr2, nr3, irt, ftau )
!
! ... if dynamics is done the system should have no symmetries
! ... (inversion symmetry alone is allowed)
!
IF ( lmd .AND. ( nsym == 2 .AND. .NOT. invsym .OR. nsym > 2 ) &
.AND. .NOT. ( calc == 'mm' .OR. calc == 'nm' ) ) &
CALL infomsg( 'setup', 'Dynamics, you should have no symmetries', -1 )
!
! ... Calculate quantities used in tetrahedra method
!
IF ( ltetra ) THEN
IF ( .NOT. lbands ) THEN
!
ntetra = 6 * nk1 * nk2 * nk3
CALL sgama( nrot, nat, s, sname, t_rev, at, bg, tau, ityp, nsym, nr1,&
nr2, nr3, irt, ftau, npk, nks, xk, wk, invsym, minus_q, &
xqq, modenum, noncolin, m_loc )
!
ALLOCATE( tetra( 4, ntetra ) )
CALL checkallsym( nsym, s, nat, tau, ityp, at, &
bg, nr1, nr2, nr3, irt, ftau )
!
CALL tetrahedra( nsym, s, minus_q, at, bg, npk, k1, k2, k3, &
nk1, nk2, nk3, nks, xk, wk, ntetra, tetra )
! ... if dynamics is done the system should have no symmetries
! ... (inversion symmetry alone is allowed)
!
ELSE
IF ( lmd .AND. ( nsym == 2 .AND. .NOT. invsym .OR. nsym > 2 ) &
.AND. .NOT. ( calc == 'mm' .OR. calc == 'nm' ) ) &
CALL infomsg( 'setup', 'Dynamics, you should have no symmetries', -1 )
!
ntetra = 0
! ... Calculate quantities used in tetrahedra method
!
END IF
!
! ... non scf calculation: do not change the number of k-points
! ... to account for reduced symmetry, unless you need to
! ... ( as in phonon or raman or DOS calculations, or whenever the
! ... Fermi energy has to be calculated )
!
ltest = ( nks /= input_nks ) .AND. &
( .NOT. ( ltetra .OR. lgauss ) ) .AND. &
( .NOT. lscf ) .AND. ( .NOT. ( lphonon .OR. lraman ) )
IF ( ltest ) THEN
!
WRITE( stdout, '(/,5X,"Only input k-points are used ", &
& "(inequivalent points not generated)",/)' )
!
nks = input_nks
IF ( ltetra ) THEN
!
ntetra = 6 * nk1 * nk2 * nk3
!
ALLOCATE( tetra( 4, ntetra ) )
!
CALL tetrahedra( nsym, s, minus_q, at, bg, npk, k1, k2, k3, &
nk1, nk2, nk3, nks, xk, wk, ntetra, tetra )
!
ELSE
!
ntetra = 0
!
END IF
!
END IF
!
@ -746,10 +731,6 @@ SUBROUTINE setup()
!
IF ( lphonon ) CALL set_kplusq( xk, wk, xqq, nks, npk )
!
! ... raman calculation: add k+b to the list of k
!
IF ( lraman ) CALL set_kplusb(ibrav, xk, wk, b_length, nks, npk, lcart)
!
#if defined (EXX)
IF ( dft_is_hybrid() ) CALL exx_grid_init()
#endif
@ -799,22 +780,6 @@ SUBROUTINE setup()
!
ENDIF
!
IF ( lraman ) THEN
!
IF( lcart ) THEN
!
kunit = 7
!
ELSE
!
IF ( ibrav == 1 ) kunit = 7
IF ( ibrav == 2 ) kunit = 9
IF ( ibrav == 3 ) kunit = 13
!
END IF
!
END IF
!
! ... distribute the k-points (and their weights and spin indices)
!
CALL divide_et_impera( xk, wk, isk, lsda, nkstot, nks )