Vdw table now handled as PP are. The table is copied in the ".save" dir

and then re-initialized from the ".save" dir after restart or when a
post processing calls pw_restart. cp_restart modified to respect the new
write_xc call.


git-svn-id: http://qeforge.qe-forge.org/svn/q-e/trunk/espresso@7484 c92efa57-630b-4861-b058-cf58834340f0
This commit is contained in:
riccardo 2011-02-09 11:03:38 +00:00
parent ee90ffe31b
commit 4d10d94165
4 changed files with 108 additions and 6 deletions

View File

@ -373,7 +373,7 @@ MODULE cp_restart
!-------------------------------------------------------------------------------
!
dft_name = get_dft_name()
CALL write_xc( DFT = dft_name, NSP = nsp, LDA_PLUS_U = .FALSE. )
CALL write_xc( DFT = dft_name, NSP = nsp, LDA_PLUS_U = .FALSE., IS_VDW = .FALSE. )
!
!-------------------------------------------------------------------------------
! ... OCCUPATIONS

View File

@ -874,6 +874,37 @@ MODULE xml_io_base
END SUBROUTINE write_ions
!
!------------------------------------------------------------------------
! SUBROUTINE write_vdw_df(vdw_table_name, resource_dir, copy_dirname)
!
! CHARACTER(LEN=*), INTENT(IN) :: vdw_table_name
! CHARACTER(LEN=*), INTENT(IN) :: resource_dir
! CHARACTER(LEN=*), INTENT(IN) :: copy_dirname
! !
! INTEGER :: i, flen
! CHARACTER(LEN=256) :: file_table
!
!
! CALL iotk_write_begin( iunpun, "VDW_KERNEL" )
! CALL iotk_write_dat( iunpun, "VDW_KERNEL_NAME", vdw_table_name )
! CALL iotk_write_end( iunpun, "VDW_KERNEL" )
!
! flen = LEN_TRIM( resource_dir )
! IF ( resource_dir(flen:flen) /= '/' ) THEN
! !
! file_table = resource_dir(1:flen) // '/' // vdw_table_name
! !
! ELSE
! !
! file_table = resource_dir(1:flen) // vdw_table_name
! !
! END IF
! !
! CALL copy_file( TRIM( file_table ), TRIM( copy_dirname ) // "/" // TRIM( vdw_table_name ) )
!
! END SUBROUTINE write_vdw_df
!
!
!------------------------------------------------------------------------
SUBROUTINE write_symmetry( ibrav, symm_type, nrot, nsym, invsym, noinv, &
time_reversal, no_t_rev, &
nr1, nr2, nr3, ftau, s, sname, irt, nat, t_rev )
@ -1163,7 +1194,8 @@ MODULE xml_io_base
!
!------------------------------------------------------------------------
SUBROUTINE write_xc( dft, nsp, lda_plus_u, &
Hubbard_lmax, Hubbard_l, Hubbard_U, Hubbard_alpha )
Hubbard_lmax, Hubbard_l, Hubbard_U, Hubbard_alpha, &
is_vdw, vdw_table_name, pseudo_dir, dirname )
!------------------------------------------------------------------------
!
CHARACTER(LEN=*), INTENT(IN) :: dft
@ -1172,7 +1204,11 @@ MODULE xml_io_base
INTEGER, OPTIONAL, INTENT(IN) :: Hubbard_lmax
INTEGER, OPTIONAL, INTENT(IN) :: Hubbard_l(:)
REAL(DP), OPTIONAL, INTENT(IN) :: Hubbard_U(:), Hubbard_alpha(:)
LOGICAL, INTENT(IN) :: is_vdw
CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: vdw_table_name, pseudo_dir, dirname
!
INTEGER :: i, flen
CHARACTER(LEN=256) :: file_table
!
CALL iotk_write_begin( iunpun, "EXCHANGE_CORRELATION" )
!
@ -1202,6 +1238,40 @@ MODULE xml_io_base
CALL iotk_write_dat( iunpun, "HUBBARD_ALPHA", Hubbard_alpha(1:nsp) )
!
END IF
!
! Vdw kernel table
!
CALL iotk_write_dat( iunpun, "VDW_DF", is_vdw )
IF ( is_vdw ) THEN
IF ( .NOT. PRESENT( vdw_table_name ) .OR. &
.NOT. PRESENT( pseudo_dir ) .OR. &
.NOT. PRESENT( dirname )) &
CALL errore( 'write_exchange_correlation', &
' variable vdw_table_name not present', 1 )
CALL iotk_write_dat( iunpun, "VDW_KERNEL_NAME", vdw_table_name )
!
! Copy the file in .save directory
!
flen = LEN_TRIM( pseudo_dir )
IF ( pseudo_dir(flen:flen) /= '/' ) THEN
!
file_table = pseudo_dir(1:flen) // '/' // vdw_table_name
!
ELSE
!
file_table = pseudo_dir(1:flen) // vdw_table_name
!
END IF
!
CALL copy_file( TRIM( file_table ), TRIM( dirname ) // "/" // TRIM( vdw_table_name ) )
END IF
!
CALL iotk_write_end( iunpun, "EXCHANGE_CORRELATION" )
!

View File

@ -26,6 +26,7 @@ MODULE pw_restart
USE mp_global, ONLY : my_pool_id, intra_image_comm, intra_pool_comm
USE mp, ONLY : mp_bcast, mp_sum, mp_max
USE parser, ONLY : version_compare
!
IMPLICIT NONE
!
@ -101,7 +102,8 @@ MODULE pw_restart
USE noncollin_module, ONLY : angle1, angle2, i_cons, mcons, bfield, &
lambda
USE ions_base, ONLY : amass
USE funct, ONLY : get_dft_name
USE funct, ONLY : get_dft_name, dft_is_vdW
USE kernel_table, ONLY : vdw_table_name
USE scf, ONLY : rho
USE extfield, ONLY : tefield, dipfield, edir, &
emaxpos, eopreg, eamp
@ -128,7 +130,7 @@ MODULE pw_restart
INTEGER :: ike, iks, npw_g, ispin
INTEGER, ALLOCATABLE :: ngk_g(:)
INTEGER, ALLOCATABLE :: igk_l2g(:,:), igk_l2g_kdip(:,:), mill_g(:,:)
LOGICAL :: lwfc
LOGICAL :: lwfc, is_vdw
REAL(DP), ALLOCATABLE :: raux(:)
!
!
@ -376,10 +378,13 @@ MODULE pw_restart
!-------------------------------------------------------------------------------
!
dft_name = get_dft_name()
is_vdw = dft_is_vdW()
!
CALL write_xc( DFT = dft_name, NSP = nsp, LDA_PLUS_U = lda_plus_u, &
HUBBARD_LMAX = Hubbard_lmax, HUBBARD_L = Hubbard_l, &
HUBBARD_U = Hubbard_U, HUBBARD_ALPHA = Hubbard_alpha )
HUBBARD_U = Hubbard_U, HUBBARD_ALPHA = Hubbard_alpha, &
IS_VDW = is_vdw, VDW_TABLE_NAME = vdw_table_name, &
PSEUDO_DIR = pseudo_dir, DIRNAME = dirname)
#ifdef EXX
CALL write_exx( x_gamma_extrapolation, nq1, nq2, nq3, &
exxdiv_treatment, yukawa, ecutvcut, &
@ -2236,6 +2241,7 @@ MODULE pw_restart
USE funct, ONLY : enforce_input_dft
USE ldaU, ONLY : lda_plus_u, Hubbard_lmax, &
Hubbard_l, Hubbard_U, Hubbard_alpha
USE kernel_table, ONLY : vdw_table_name
!
IMPLICIT NONE
!
@ -2244,7 +2250,7 @@ MODULE pw_restart
!
CHARACTER(LEN=20) :: dft_name
INTEGER :: nsp_
LOGICAL :: found, nomsg = .true.
LOGICAL :: found, nomsg = .true., is_vdw
!
ierr = 0
IF ( lxc_read ) RETURN
@ -2283,6 +2289,18 @@ MODULE pw_restart
CALL iotk_scan_dat( iunpun, "HUBBARD_ALPHA", Hubbard_alpha(1:nsp_) )
!
END IF
!
! Vdw DF
!
CALL iotk_scan_dat( iunpun, "VDW_DF", is_vdw, FOUND = found )
IF ( .NOT. found ) is_vdw = .FALSE.
!
IF ( is_vdw ) THEN
!
CALL iotk_scan_dat( iunpun, "VDW_KERNEL_NAME", vdw_table_name )
!
END IF
!
CALL iotk_scan_end( iunpun, "EXCHANGE_CORRELATION" )
!
@ -2292,6 +2310,7 @@ MODULE pw_restart
!
CALL mp_bcast( dft_name, ionode_id, intra_image_comm )
CALL mp_bcast( lda_plus_u, ionode_id, intra_image_comm )
CALL mp_bcast( is_vdw, ionode_id, intra_image_comm )
!
IF ( lda_plus_u ) THEN
!
@ -2301,6 +2320,11 @@ MODULE pw_restart
CALL mp_bcast( Hubbard_alpha, ionode_id, intra_image_comm )
!
END IF
IF ( is_vdw ) THEN
CALL mp_bcast( vdw_table_name, ionode_id, intra_image_comm )
END IF
! discard any further attempt to set a different dft
CALL enforce_input_dft( dft_name, nomsg )
!

View File

@ -54,6 +54,8 @@ SUBROUTINE read_file()
USE io_global, ONLY : stdout
USE dfunct, ONLY : newd
USE control_flags, ONLY : gamma_only
USE funct, ONLY : dft_is_vdW
USE kernel_table, ONLY : initialize_kernel_table
!
IMPLICIT NONE
!
@ -183,6 +185,12 @@ SUBROUTINE read_file()
!
CALL readpp()
!
! ... read the vdw kernel table if needed
!
if (dft_is_vdW()) then
call initialize_kernel_table()
endif
!
okvan = ANY ( upf(:)%tvanp )
okpaw = ANY ( upf(1:nsp)%tpawp )
!