update of release notes and small fixes to xml

This commit is contained in:
Pietro Delugas 2019-02-28 15:21:19 +01:00
parent c731fe86b7
commit bd1fedf22f
6 changed files with 35 additions and 7 deletions

View File

@ -1,4 +1,13 @@
New in development version:
New in version 6.4:
* Experimental and specific for gamma_only case: specifing nscdm the SCDM localization is performed
only for iterations multiples of nscdm, in the intermediate iterations the localized orbitals are
derived with parallel transport from the last SCDM localization.
* Added experimental version of SCDM localization with many K_POINTS. The calculation using SCDM
is set as in the gamma-only case just specifing localization_thr to any value greater than 0 in
the system namelist.
* It is now possible to limit the number of xml step elements printed out for relaxation or
molecular dynamics simulation setting the envinroment variable MAX_XML_STEPS, useful in case
of very long trajectories to avoid issues due to too large file size.
* EPW works with ultrasoft pseudopotentials (F. Giustino, S. Poncé, R. Margine)
* New code hp.x to compute Hubbard parameters using density-functional
perturbation theory (experimental stage) (I. Timrov, N. Marzari, and M. Cococcioni,

View File

@ -121,7 +121,9 @@ MODULE input_parameters
INTEGER :: iprint = 10
! number of steps/scf iterations between successive writings
! of relevant physical quantities to standard output
INTEGER :: max_xml_steps = 0
! max number of steps between successive appending of an xml step
! in the xml data file, default 0 means all steps are printed.
INTEGER :: isave = 100
! number of steps between successive savings of
! information needed to restart the run (see "ndr", "ndw")

View File

@ -1317,7 +1317,7 @@ CONTAINS
step_counter = step_counter+1
!
step_obj%tagname="step"
step_obj%n_step = step_counter
step_obj%n_step = i_step
!
CALL qes_init( scf_conv_obj,"scf_conv", scf_has_converged, n_scf_steps, scf_error )
!

View File

@ -60,6 +60,7 @@ MODULE read_namelists_module
IMPLICIT NONE
!
CHARACTER(LEN=2) :: prog ! ... specify the calling program
CHARACTER(LEN=20) :: temp_string
!
!
IF ( prog == 'PW' ) THEN
@ -105,6 +106,10 @@ MODULE read_namelists_module
pseudo_dir = TRIM( pseudo_dir ) // '/espresso/pseudo/'
END IF
!
! ... max number of md steps added to the xml file. Needs to be limited for very long
! md simulations
CALL get_environment_variable('MAX_XML_STEPS', temp_string)
IF ( TRIM(temp_string) .NE. ' ') READ(temp_string, *) max_xml_steps
refg = 0.05_DP
max_seconds = 1.E+7_DP
ekin_conv_thr = 1.E-6_DP

View File

@ -27,6 +27,7 @@ USE force_mod, ONLY: force, sigma
USE control_flags,ONLY: nstep, n_scf_steps, scf_error, conv_elec
USE fcp_variables,ONLY: fcp_mu, lfcpopt, lfcpdyn
USE extfield, ONLY: gate, etotgatefield, tefield, etotefield
USE input_parameters, ONLY: max_xml_steps
!-----------------------------------------------------------------------------
! END_GLOBAL_VARIABLES
!-----------------------------------------------------------------------------
@ -53,6 +54,17 @@ REAL(DP),TARGET :: potstat_contr_tgt, fcp_force_tgt, fcp_tot_charge_
REAL(DP),POINTER :: potstat_contr_ptr, fcp_force_ptr, fcp_tot_charge_ptr,&
demet_ptr, degauss_ptr, gatefield_en_ptr, efield_corr_ptr
!
INTEGER :: stride = 1, max_xml_steps_
IF ( max_xml_steps > 0 ) THEN
stride = nstep/max_xml_steps
max_xml_steps_ = max_xml_steps+2
ELSE
max_xml_steps_ = nstep
END IF
IF (.NOT. ( i_step == 1 .OR. MOD(i_step-1, stride) == 0 .OR. i_step == nstep)) RETURN
NULLIFY(potstat_contr_ptr, fcp_force_ptr, fcp_tot_charge_ptr, demet_ptr, degauss_ptr, &
gatefield_en_ptr, efield_corr_ptr)
!
@ -65,8 +77,8 @@ END IF
IF ( lfcpopt .OR. lfcpdyn ) THEN
potstat_contr_tgt = ef * tot_charge / e2
potstat_contr_ptr => potstat_contr_tgt
!FIXME ( again shouldn't we use Hartree units for this ? )
fcp_force_tgt = fcp_mu - ef
!
fcp_force_tgt = (fcp_mu - ef)/e2
fcp_force_ptr => fcp_force_tgt
!
fcp_tot_charge_tgt = tot_charge
@ -81,7 +93,7 @@ IF (tefield) THEN
efield_corr_tgt = etotefield/e2
efield_corr_ptr => efield_corr_tgt
END IF
CALL qexsd_step_addstep ( i_step, nstep, nsp, atm, ityp, nat, alat*tau, alat, alat*at(:,1), &
CALL qexsd_step_addstep ( i_step, max_xml_steps_, nsp, atm, ityp, nat, alat*tau, alat, alat*at(:,1), &
alat*at(:,2), alat*at(:,3), etot/e2, eband/e2, ehart/e2, vtxc/e2, etxc/e2, &
ewld/e2, degauss_ptr, demet_ptr, force/e2, sigma/e2, conv_elec, n_scf_steps, scf_error, &
FCP_FORCE = fcp_force_ptr , FCP_TOT_CHARGE = fcp_tot_charge_ptr,&

View File

@ -246,7 +246,7 @@ MODULE pw_restart_new
n_opt_steps = istep
END IF
scf_has_converged = conv_elec
n_scf_steps = n_scf_steps
n_scf_steps_ = n_scf_steps
CASE ("nscf", "bands" )
n_opt_steps = 0
scf_has_converged = .FALSE.