TABLE OF CONTENTS


ABINIT/xfh_update [ Functions ]

[ Top ] [ Functions ]

NAME

 xfh_update

FUNCTION

 Update the contents of the history xfhist taking values
 from xred, acell, rprim, fred_corrected and strten

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR, JCC, SE)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors,
 see ~abinit/doc/developers/contributors.txt .

INPUTS

OUTPUT

SIDE EFFECTS

PARENTS

      mover

CHILDREN

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 
 40 subroutine xfh_update(ab_xfh,acell,fred_corrected,natom,rprim,strten,xred)
 41 
 42  use m_profiling_abi
 43  use defs_basis
 44  use m_abimover
 45 
 46 !This section has been created automatically by the script Abilint (TD).
 47 !Do not modify the following lines by hand.
 48 #undef ABI_FUNC
 49 #define ABI_FUNC 'xfh_update'
 50 !End of the abilint section
 51 
 52 implicit none
 53 
 54 !Arguments ------------------------------------
 55 !scalars
 56 
 57 type(ab_xfh_type),intent(inout) :: ab_xfh
 58 integer,intent(in) :: natom
 59 
 60 !arrays
 61 real(dp),intent(in) :: acell(3)
 62 real(dp),intent(in) :: xred(3,natom)
 63 real(dp),intent(in) :: rprim(3,3)
 64 real(dp),intent(in) :: fred_corrected(3,natom)
 65 real(dp),intent(in) :: strten(6)
 66 
 67 !Local variables-------------------------------
 68 !scalars
 69 !integer :: kk
 70 
 71 !*********************************************************************
 72 
 73 !DEBUG
 74 !write (ab_out,*) '---WROTE TO XFHIST---'
 75 
 76 !write (ab_out,*) 'XRED'
 77 !do kk=1,natom
 78 !write (ab_out,*) xred(:,kk)
 79 !end do
 80 !write (ab_out,*) 'FRED'
 81 !do kk=1,natom
 82 !write (ab_out,*) fred_corrected(:,kk)
 83 !end do
 84 !write(ab_out,*) 'RPRIM'
 85 !do kk=1,3
 86 !write(ab_out,*) rprim(:,kk)
 87 !end do
 88 !write(ab_out,*) 'ACELL'
 89 !write(ab_out,*) acell(:)
 90 !DEBUG
 91 
 92  ab_xfh%nxfh=ab_xfh%nxfh+1
 93 
 94  ab_xfh%xfhist(:,1:natom,1,ab_xfh%nxfh)=xred(:,:)
 95  ab_xfh%xfhist(:,natom+1,1,ab_xfh%nxfh)=acell(:)
 96  ab_xfh%xfhist(:,natom+2:natom+4,1,ab_xfh%nxfh)=rprim(:,:)
 97  ab_xfh%xfhist(:,1:natom,2,ab_xfh%nxfh)=fred_corrected(:,:)
 98  ab_xfh%xfhist(:,natom+2,2,ab_xfh%nxfh)=strten(1:3)
 99  ab_xfh%xfhist(:,natom+3,2,ab_xfh%nxfh)=strten(4:6)
100 
101 end subroutine xfh_update