Default value for nr1b, nr2b, nr3b is set to 0. The nr*b MUST be provided

for US PP, are completely ignored (box grid is not initialized) for NC PP.
This prevents nasty errors in US PP calculations if nr*b are forgotten


git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@4045 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
giannozz 2007-07-18 10:23:06 +00:00
parent ec79afc5ba
commit 7b8c1ad33e
13 changed files with 81 additions and 58 deletions

View File

@ -1107,6 +1107,13 @@ END FUNCTION
REAL(DP) :: x(3), xmod
INTEGER :: nr(3), nrb(3), xint, is, ia, i, isa
!
IF ( nr1b < 1) CALL errore &
('initbox', 'incorrect value for box grid dimensions', 1)
IF ( nr2b < 1) CALL errore &
('initbox', 'incorrect value for box grid dimensions', 2)
IF ( nr3b < 1) CALL errore &
('initbox', 'incorrect value for box grid dimensions', 3)
nr (1)=nr1
nr (2)=nr2
nr (3)=nr3

View File

@ -24,7 +24,7 @@ SUBROUTINE cprmain( tau_out, fion_out, etot_out )
USE core, ONLY : nlcc_any, rhoc
USE uspp_param, ONLY : nhm, nh
USE cvan, ONLY : nvb, ish
USE uspp, ONLY : nkb, vkb, becsum, deeq
USE uspp, ONLY : nkb, vkb, becsum, deeq, okvan
USE energies, ONLY : eht, epseu, exc, etot, eself, enl, &
ekin, atot, entropy, egrand, enthal, &
ekincm, print_energies
@ -256,7 +256,7 @@ SUBROUTINE cprmain( tau_out, fion_out, etot_out )
!
END IF
!
IF ( tfor .OR. thdyn .OR. tfirst ) THEN
IF ( okvan .AND. (tfor .OR. thdyn .OR. tfirst) ) THEN
!
CALL initbox( tau0, taub, irb, ainv, a1, a2, a3 )
!
@ -718,8 +718,10 @@ SUBROUTINE cprmain( tau_out, fion_out, etot_out )
! ... in this case optimize c0 and lambda for smooth
! ... restart with CP
!
CALL initbox( tau0, taub, irb, ainv, a1, a2, a3 )
CALL phbox( taub, eigrb, ainvb )
IF ( okvan) THEN
CALL initbox( tau0, taub, irb, ainv, a1, a2, a3 )
CALL phbox( taub, eigrb, ainvb )
END IF
CALL r_to_s( tau0, taus, na, nsp, ainv )
CALL phfacs( ei1, ei2, ei3, eigr, mill_l, taus, nr1, nr2, nr3, nat )
CALL strucf( sfac, ei1, ei2, ei3, mill_l, ngs )

View File

@ -49,7 +49,7 @@ CONTAINS
USE energies, ONLY : entropy, eself, enl, ekin, enthal, etot, ekincm
USE energies, ONLY : dft_energy_type, debug_energies
USE dener, ONLY : denl, denl6, dekin6, detot
USE uspp, ONLY : vkb, becsum, deeq, nkb
USE uspp, ONLY : vkb, becsum, deeq, nkb, okvan
USE io_global, ONLY : stdout, ionode
USE core, ONLY : nlcc_any
USE gvecw, ONLY : ngw
@ -149,11 +149,12 @@ CONTAINS
!
CALL strucf( sfac, ei1, ei2, ei3, mill_l, ngs )
!
CALL initbox ( tau0, taub, irb, ainv, a1, a2, a3 )
!
CALL phbox( taub, eigrb, ainvb )
IF ( okvan ) THEN
CALL initbox ( tau0, taub, irb, ainv, a1, a2, a3 )
CALL phbox( taub, eigrb, ainvb )
END IF
!
! random initialization
! wfc initialization with random numbers
!
CALL wave_rand_init( cm, nbsp, 1 )
!

View File

@ -159,18 +159,21 @@
!
! generation of little box g-vectors
!
! sets the small box parameters
IF ( nr1b > 0 .AND. nr2b > 0 .AND. nr3b > 0 ) THEN
rat1 = DBLE( nr1b ) / DBLE( nr1 )
rat2 = DBLE( nr2b ) / DBLE( nr2 )
rat3 = DBLE( nr3b ) / DBLE( nr3 )
CALL small_box_set( alat, omega, a1, a2, a3, rat1, rat2, rat3 )
! sets the small box parameters
! now set gcutb
!
gcutb = ecut / tpibab / tpibab
!
CALL ggenb ( b1b, b2b, b3b, nr1b, nr2b, nr3b, nr1bx, nr2bx, nr3bx, gcutb )
rat1 = DBLE( nr1b ) / DBLE( nr1 )
rat2 = DBLE( nr2b ) / DBLE( nr2 )
rat3 = DBLE( nr3b ) / DBLE( nr3 )
CALL small_box_set( alat, omega, a1, a2, a3, rat1, rat2, rat3 )
! now set gcutb
gcutb = ecut / tpibab / tpibab
!
CALL ggenb ( b1b, b2b, b3b, nr1b, nr2b, nr3b, nr1bx, nr2bx, nr3bx, gcutb )
END IF
! ... printout g vector distribution summary
!

View File

@ -66,7 +66,7 @@
tprnsfac, tcarpar, &
tdipole, &
tnosee, tnosep, force_pairing, tconvthrs, convergence_criteria, tionstep, nstepe, &
ekin_conv_thr, ekin_maxiter, conv_elec, lneb, tnoseh, tuspp, etot_conv_thr, tdamp
ekin_conv_thr, ekin_maxiter, conv_elec, lneb, tnoseh, etot_conv_thr, tdamp
USE atoms_type_module, ONLY: atoms_type
USE cell_base, ONLY: press, wmass, boxdimensions, updatecell, cell_force, cell_move, gethinv
USE polarization, ONLY: ddipole
@ -107,7 +107,7 @@
vnhh, xnhh0, xnhhm, xnhhp, qnh, temph
USE cell_base, ONLY: cell_gamma
USE grid_subroutines, ONLY: realspace_grids_init, realspace_grids_para
USE uspp, ONLY: vkb, nkb
USE uspp, ONLY: vkb, nkb, okvan
!
USE reciprocal_vectors, ONLY: &
g, & ! G-vectors square modulus

View File

@ -193,13 +193,14 @@ END FUNCTION calculate_dx
USE mp, ONLY: mp_bcast, mp_sum
USE io_global, ONLY: stdout, ionode, ionode_id
USE uspp, ONLY : okvan
USE uspp_param, ONLY : zp, tvanp, oldvan
USE atom, ONLY: numeric, nlcc, oc, lchi, nchi
USE cvan, ONLY: nvb
use ions_base, only: zv, nsp
use read_upf_module, only: read_pseudo_upf
use read_uspp_module, only: readvan, readrrkj
use control_flags, only: program_name, tuspp
use control_flags, only: program_name
use funct, only: get_iexch, get_icorr, get_igcx, get_igcc, set_dft_from_name, dft_is_hybrid
USE upf_to_internal, ONLY: set_pseudo_upf
@ -315,10 +316,7 @@ END FUNCTION calculate_dx
ELSE
call set_pseudo_upf( is, upf( is ) )
!
IF( upf(is)%tvanp ) THEN
tuspp = .TRUE.
ELSE
tuspp = .FALSE.
IF( .NOT. upf(is)%tvanp ) THEN
CALL upf2ncpp( upf(is), ap(is) )
END IF
!
@ -391,7 +389,7 @@ END FUNCTION calculate_dx
! count u-s vanderbilt species
!
if (tvanp(is)) nvb=nvb+1
!
END IF
if ( xc_type /= 'none' ) then
@ -429,6 +427,8 @@ END FUNCTION calculate_dx
CALL check_types_order()
END IF
okvan = ( nvb > 0 )
RETURN
END SUBROUTINE readpp

View File

@ -502,6 +502,7 @@ SUBROUTINE from_restart_x( &
USE efield_module, ONLY : efield_berry_setup, tefield, &
efield_berry_setup2, tefield2
USE small_box, ONLY : ainvb
USE uspp, ONLY : okvan
!
IMPLICIT NONE
!
@ -647,9 +648,10 @@ SUBROUTINE from_restart_x( &
! ... computes form factors and initializes nl-pseudop. according
! ... to starting cell (from ndr or again standard input)
!
CALL initbox( tau0, taub, irb, ainv, a1, a2, a3 )
!
CALL phbox( taub, eigrb, ainvb )
IF ( okvan) THEN
CALL initbox( tau0, taub, irb, ainv, a1, a2, a3 )
CALL phbox( taub, eigrb, ainvb )
END IF
!
CALL phfacs( ei1, ei2, ei3, eigr, mill_l, taus, nr1, nr2, nr3, nat )
!

View File

@ -22,7 +22,7 @@ SUBROUTINE smdmain( tau, fion_out, etot_out, nat_out )
USE atom, ONLY : nlcc
USE core, ONLY : nlcc_any, rhoc
USE uspp_param, ONLY : nhm
USE uspp, ONLY : nkb, vkb, rhovan => becsum, deeq
USE uspp, ONLY : nkb, vkb, rhovan => becsum, deeq, okvan
USE cvan, ONLY : nvb
USE energies, ONLY : eht, epseu, exc, etot, eself, enl, ekin, &
esr, print_energies
@ -474,9 +474,10 @@ SUBROUTINE smdmain( tau, fion_out, etot_out, nat_out )
CALL phfacs( ei1, ei2, ei3, eigr, mill_l, rep(sm_k)%taus, nr1, nr2, nr3, nat )
!
CALL initbox ( rep(sm_k)%tau0, taub, irb, ainv, a1, a2, a3 )
!
CALL phbox( taub, eigrb, ainvb )
IF ( okvan ) THEN
CALL initbox ( rep(sm_k)%tau0, taub, irb, ainv, a1, a2, a3 )
CALL phbox( taub, eigrb, ainvb )
END IF
!
IF(trane) THEN
IF(sm_k == 1) THEN

View File

@ -45,7 +45,6 @@ MODULE control_flags
PUBLIC :: tksw, trhor, thdyn, iprsta, trhow
PUBLIC :: twfcollect, printwfc
PUBLIC :: lkpoint_dir
PUBLIC :: tuspp
PUBLIC :: program_name
!
! ... declare execution control variables
@ -85,7 +84,6 @@ MODULE control_flags
LOGICAL :: tscreen = .FALSE. ! Use screened coulomb potentials for cluster calculations
LOGICAL :: twfcollect = .FALSE. ! Collect wave function in the restart file at the end of run.
LOGICAL :: lkpoint_dir = .TRUE. ! save each k point in a different directory
LOGICAL :: tuspp = .FALSE. ! Ultra-soft pseudopotential are being used
INTEGER :: printwfc = -1 ! Print wave functions, temporarely used only by ensemble-dft
LOGICAL :: force_pairing = .FALSE. ! Force pairing
LOGICAL :: tchi2 = .FALSE. ! Compute Chi^2

View File

@ -161,15 +161,21 @@
CALL errore(' realspace_grids_init ', ' smooth grid larger than dense grid?',1)
END IF
IF( nr1b < 1 ) nr1b = 1
IF( nr2b < 1 ) nr2b = 1
IF( nr3b < 1 ) nr3b = 1
! no default values for grid box: if nr*b=0, ignore
IF( nr1b > 0 .AND. nr2b > 0 .AND. nr3b > 0 ) THEN
nr1b = good_fft_order( nr1b ) ! small box is not parallelized
nr2b = good_fft_order( nr2b )
nr3b = good_fft_order( nr3b )
nr1bx = good_fft_dimension( nr1b )
ELSE
nr1bx = nr1b
END IF
nr1b = good_fft_order( nr1b ) ! small box is not parallelized
nr2b = good_fft_order( nr2b ) ! small box is not parallelized
nr3b = good_fft_order( nr3b ) ! small box is not parallelized
nr1bx = good_fft_dimension( nr1b )
nr2bx = nr2b
nr3bx = nr3b
nnrbx = nr1bx * nr2bx * nr3bx
@ -253,13 +259,16 @@
WRITE( stdout,*) ' Number of x-y planes for each processors: '
WRITE( stdout, fmt = '( 3X, "nr3sl = ", 10I5 )' ) ( dffts%npp( i ), i = 1, nproc_image )
WRITE( stdout,*)
WRITE( stdout,*) ' Small Box Real Mesh'
WRITE( stdout,*) ' -------------------'
WRITE( stdout,1000) nr1b, nr2b, nr3b, nr1bl, nr2bl, nr3bl, 1, 1, 1
WRITE( stdout,1010) nr1bx, nr2bx, nr3bx
WRITE( stdout,1020) nnrbx
IF ( nr1b > 0 .AND. nr2b > 0 .AND. nr3b > 0 ) THEN
WRITE( stdout,*)
WRITE( stdout,*) ' Small Box Real Mesh'
WRITE( stdout,*) ' -------------------'
WRITE( stdout,1000) nr1b, nr2b, nr3b, nr1bl, nr2bl, nr3bl, 1, 1, 1
WRITE( stdout,1010) nr1bx, nr2bx, nr3bx
WRITE( stdout,1020) nnrbx
END IF
END IF
1000 FORMAT(3X, &

View File

@ -3588,7 +3588,7 @@ SUBROUTINE para_dcholdc( n, a, lda, comm )
!----------------------------------------------------------------------------
!
! ... trivial parallelization (using a parallel version of DGEMV) of
! ... the Cholesky deconposition (equivalent to DPOTF2)
! ... the Cholesky decomposition (equivalent to DPOTF2)
!
USE kinds, ONLY : DP
!
@ -3637,7 +3637,7 @@ SUBROUTINE para_zcholdc( n, a, lda, comm )
!----------------------------------------------------------------------------
!
! ... trivial parallelization (using a parallel version of ZGEMV) of
! ... the Cholesky deconposition (equivalent to ZPOTF2)
! ... the Cholesky decomposition (equivalent to ZPOTF2)
!
USE kinds, ONLY : DP
!
@ -3691,7 +3691,7 @@ END SUBROUTINE para_zcholdc
SUBROUTINE para_dtrtri( n, a, lda, comm )
!----------------------------------------------------------------------------
!
! ... parallel inversion of a lower trinagular matrix done distributing
! ... parallel inversion of a lower triangular matrix done distributing
! ... by columns ( the number of columns assigned to each processor are
! ... chosen to optimize the load balance )
!
@ -3795,7 +3795,7 @@ END SUBROUTINE para_dtrtri
SUBROUTINE para_ztrtri( n, a, lda, comm )
!----------------------------------------------------------------------------
!
! ... parallel inversion of a lower trinagular matrix done distributing
! ... parallel inversion of a lower triangular matrix done distributing
! ... by columns ( the number of columns assigned to each processor are
! ... chosen to optimize the load balance in the limit of large matrices )
!

View File

@ -149,9 +149,9 @@ MODULE read_namelists_module
nr1s = 0
nr2s = 0
nr3s = 0
nr1b = 3
nr2b = 3
nr3b = 3
nr1b = 0
nr2b = 0
nr3b = 0
occupations = 'fixed'
smearing = 'gaussian'
degauss = 0.0_DP

View File

@ -80,7 +80,7 @@ MODULE uspp
nhtolm(:,:) ! correspondence n <-> combined lm index for (l,m)
!
LOGICAL :: &
okvan ! if .TRUE. at least one pseudo is Vanderbilt
okvan = .FALSE. ! if .TRUE. at least one pseudo is Vanderbilt
!
COMPLEX(DP), ALLOCATABLE, TARGET :: &
vkb(:,:) ! all beta functions in reciprocal space