TABLE OF CONTENTS


ABINIT/symchk [ Functions ]

[ Top ] [ Functions ]

NAME

 symchk

FUNCTION

 Symmetry checker for atomic coordinates.
 Checks for translated atomic coordinate tratom(3) to agree
 with some coordinate xred(3,iatom) where atomic types agree too.
 All coordinates are "reduced", i.e. given in terms of primitive
 reciprocal translations.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR)
 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

 natom=number of atoms in unit cell
 tratom(3)=reduced coordinates for a single atom which presumably
   result from the application of a symmetry operation to an atomic
   coordinate
 trtypat=type of atom (integer) translated to tratom
 typat(natom)=types of all atoms in unit cell (integer)
 xred(3,natom)=reduced coordinates for all atoms in unit cell

OUTPUT

 difmin(3)=minimum difference between apparently equivalent atoms
   (give value separately for each coordinate)--note that value
   may be NEGATIVE so take abs later if needed
 eatom=atom label of atom which is SAME as tratom to within a primitive
   cell translation ("equivalent atom")
 transl(3)=primitive cell translation to make iatom same as tratom (integers)

PARENTS

      m_polynomial_coeff,symatm

CHILDREN

SOURCE

 44 #if defined HAVE_CONFIG_H
 45 #include "config.h"
 46 #endif
 47 
 48 #include "abi_common.h"
 49 
 50 
 51 subroutine symchk(difmin,eatom,natom,tratom,transl,trtypat,typat,xred)
 52 
 53  use defs_basis
 54  use m_profiling_abi
 55  use m_errors
 56 
 57 !This section has been created automatically by the script Abilint (TD).
 58 !Do not modify the following lines by hand.
 59 #undef ABI_FUNC
 60 #define ABI_FUNC 'symchk'
 61 !End of the abilint section
 62 
 63  implicit none
 64 
 65 !Arguments ------------------------------------
 66 !scalars
 67  integer,intent(in) :: natom,trtypat
 68  integer,intent(out) :: eatom
 69 !arrays
 70  integer,intent(in) :: typat(natom)
 71  integer,intent(out) :: transl(3)
 72  real(dp),intent(in) :: tratom(3),xred(3,natom)
 73  real(dp),intent(out) :: difmin(3)
 74 
 75 !Local variables-------------------------------
 76 !scalars
 77  integer :: iatom,jatom,trans1,trans2,trans3
 78  real(dp) :: test,test1,test2,test3,testmn
 79 
 80 ! *************************************************************************
 81 
 82 !Start testmn out at large value
 83  testmn=1000000.d0
 84 
 85 !Loop through atoms--
 86 !when types agree, check for agreement after primitive translation
 87  jatom=1
 88  do iatom=1,natom
 89    if (trtypat/=typat(iatom)) cycle
 90 
 91 !  Check all three components
 92    test1=tratom(1)-xred(1,iatom)
 93    test2=tratom(2)-xred(2,iatom)
 94    test3=tratom(3)-xred(3,iatom)
 95 !  Find nearest integer part of difference
 96    trans1=nint(test1)
 97    trans2=nint(test2)
 98    trans3=nint(test3)
 99 !  Check whether, after translation, they agree
100    test1=test1-dble(trans1)
101    test2=test2-dble(trans2)
102    test3=test3-dble(trans3)
103 
104    test=abs(test1)+abs(test2)+abs(test3)
105    if (test<tol10) then
106 !    Note that abs() is not taken here
107      difmin(1)=test1
108      difmin(2)=test2
109      difmin(3)=test3
110      jatom=iatom
111      transl(1)=trans1
112      transl(2)=trans2
113      transl(3)=trans3
114 !    Break out of loop when agreement is within tolerance
115      exit
116    else
117 !    Keep track of smallest difference if greater than tol10
118      if (test<testmn) then
119        testmn=test
120 !      Note that abs() is not taken here
121        difmin(1)=test1
122        difmin(2)=test2
123        difmin(3)=test3
124        jatom=iatom
125        transl(1)=trans1
126        transl(2)=trans2
127        transl(3)=trans3
128      end if
129    end if
130 
131 !  End loop over iatom. Note a "cycle" and an "exit" inside the loop
132  end do
133 
134  eatom=jatom
135 
136 end subroutine symchk