TABLE OF CONTENTS


ABINIT/erlxconv [ Functions ]

[ Top ] [ Functions ]

NAME

  erlxconv

FUNCTION

  FIXME: add description.

COPYRIGHT

  Copyright (C) 2016-2017 ABINIT group (DW)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

  argin(sizein)=description

OUTPUT

  argout(sizeout)=description

SIDE EFFECTS

NOTES

PARENTS

      mover

CHILDREN

      wrtout

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 
 40 subroutine erlxconv(hist,iexit,itime,itime_hist,ntime,tolmxde)
 41     
 42  use defs_basis
 43  use m_errors
 44  use m_profiling_abi
 45 
 46  use m_fstrings,  only : indent
 47  use m_abihist, only : abihist,abihist_findIndex
 48 
 49 !This section has been created automatically by the script Abilint (TD).
 50 !Do not modify the following lines by hand.
 51 #undef ABI_FUNC
 52 #define ABI_FUNC 'erlxconv'
 53  use interfaces_14_hidewrite
 54 !End of the abilint section
 55 
 56  implicit none
 57 
 58 !Arguments ------------------------------------
 59 !scalars
 60  integer,intent(in) :: itime,itime_hist,ntime
 61  integer,intent(inout) :: iexit
 62  real(dp), intent(in) :: tolmxde
 63 !arrays
 64  type(abihist),intent(inout) :: hist
 65 
 66 !Local variables-------------------------------
 67  integer :: ihist,ihist_prev,ihist_prev2
 68  real(dp) :: ediff1,ediff2,maxediff
 69  character(len=500) :: message
 70 ! *************************************************************************
 71 
 72  if (itime_hist<3) then
 73    write(message, '(a,a,a)' ) ch10,&
 74 &   ' erlxconv : minimum 3 Broyd/MD steps to check convergence of energy in relaxations',ch10
 75    call wrtout(std_out,message,'COLL')
 76  else
 77    ihist = hist%ihist
 78    ihist_prev  = abihist_findIndex(hist,-1)
 79    ihist_prev2 = abihist_findIndex(hist,-2)
 80    ediff1 = hist%etot(ihist) - hist%etot(ihist_prev)
 81    ediff2 = hist%etot(ihist) - hist%etot(ihist_prev2)
 82    if ((abs(ediff1)<tolmxde).and.(abs(ediff2)<tolmxde)) then
 83      write(message, '(a,a,i4,a,a,a,a,a,es11.4,a,a)' ) ch10,&
 84 &     ' At Broyd/MD step',itime,', energy is converged : ',ch10,&
 85 &     '  the difference in energy with respect to the two ',ch10,&
 86 &     '  previous steps is < tolmxde=',tolmxde,' ha',ch10
 87      call wrtout(ab_out,message,'COLL')
 88      call wrtout(std_out,message,'COLL')
 89      iexit=1
 90    else
 91      maxediff = max(abs(ediff1),abs(ediff2))
 92      if(iexit==1)then
 93        write(message, '(a,a,a,a,i5,a,a,a,es11.4,a,es11.4,a,a)' ) ch10,&
 94 &       ' erlxconv : WARNING -',ch10,&
 95 &       '  ntime=',ntime,' was not enough Broyd/MD steps to converge energy: ',ch10,&
 96 &       '  max difference in energy =',maxediff,' > tolmxde=',tolmxde,' ha',ch10
 97        call wrtout(std_out,message,'COLL')
 98        call wrtout(ab_out,message,'COLL')
 99 
100        write(std_out,"(8a)")ch10,&
101 &       "--- !RelaxConvergenceWarning",ch10,&
102 &       "message: | ",ch10,TRIM(indent(message)),ch10,&
103 &       "..."
104      else
105        write(message, '(a,a,i4,a,a,a,es11.4,a,es11.4,a,a)' ) ch10,&
106 &       ' erlxconv : at Broyd/MD step',itime,', energy has not converged yet. ',ch10,&
107 &       '  max difference in energy=',maxediff,' > tolmxde=',tolmxde,' ha',ch10
108        call wrtout(std_out,message,'COLL')
109      end if
110    end if
111  end if
112 
113 end subroutine erlxconv