TABLE OF CONTENTS


ABINIT/symatm [ Functions ]

[ Top ] [ Functions ]

NAME

 symatm

FUNCTION

 This routine has been improved using ideas of p. 649 of notes,
 implementing suggestion of Andrew Horsfield: replace search for
 equivalent atoms using direct primitive cell translations by
 use of dot product relation which must produce an integer.
 Relation: $[inv(S(i))*(x(a)-tnons(i)) - x(inv(S)(i,a))] = integer$
 where $S(i) =$ symmetry matrix in real space, tnons=nonsymmorphic translation
 (may be 0 0 0), and $x(inv(S)(i,a))$ is sought atom into which $x(a)$ gets
 rotated by $inv(S)$.  Integer gives primitive translation coordinates to get
 back to original unit cell.
 Equivalent to $S*t(b)+tnons-x(a)=another$ $integer$ for $x(b)=x(inv(S))$.
 For each symmetry operation, find the number of the position to
 which each atom is sent in the unit cell by the INVERSE of the
 symmetry operation inv(symrel); i.e. this is the atom which, when acted
 upon by the given symmetry element isym, gets transformed into atom iatom.
 This routine uses the fact that inv(symrel)=trans(symrec),
 the inverse of the symmetry operation expressed in the basis of real
 space primitive translations equals the transpose of the same symmetry
 operation expressed in the basis of reciprocal space primitive transl.
 $xred(nu,indsym(4,isym,ia))=symrec(mu,nu,isym)*(xred(mu,ia)-tnons(mu,isym))
 - transl(mu)$ where $transl$ is also a set of integers and
 where translation transl places coordinates within unit cell (note sign).
 Note that symrec is the set of arrays which are actually input here.
 These arrays have integer elements.
 tnons is the nonsymmorphic translation or else is zero.
 If nsym=1 (i.e. only the identity symmetry is present) then
 indsym merely takes each atom into itself.
 The array of integer translations "transl" gets included within array
 "indsym" as seen below.

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 cell.
 nsym=number of space group symmetries.
 symrec(3,3,nsym)=symmetries expressed in terms of their action on
                  reciprocal space primitive translations (integer).
 tnons(3,nsym)=nonsymmorphic translations for each symmetry (would
               be 0 0 0 each for a symmorphic space group)
 typat(natom)=integer identifying type of atom.
 xred(3,natom)=reduced coordinates of atoms in terms of real space
               primitive translations
 tolsym=tolerance for the symmetries

OUTPUT

 indsym(4,nsym,natom)=indirect indexing array described above: for each
                      isym,iatom, fourth element is label of atom into
                      which iatom is sent by INVERSE of symmetry operation
                      isym; first three elements are the primitive
                      translations which must be subtracted after the
                      transformation to get back to the original unit cell.

PARENTS

      get_npert_rbz,ingeo,initberry,initorbmag,m_ab7_symmetry,m_crystal,m_ddb
      m_polynomial_coeff,m_tdep_sym,setsym,thmeig

CHILDREN

      symchk,wrtout

SOURCE

 72 #if defined HAVE_CONFIG_H
 73 #include "config.h"
 74 #endif
 75 
 76 #include "abi_common.h"
 77 
 78 subroutine symatm(indsym,natom,nsym,symrec,tnons,tolsym,typat,xred)
 79 
 80  use defs_basis
 81  use m_errors
 82  use m_profiling_abi
 83 
 84 !This section has been created automatically by the script Abilint (TD).
 85 !Do not modify the following lines by hand.
 86 #undef ABI_FUNC
 87 #define ABI_FUNC 'symatm'
 88  use interfaces_14_hidewrite
 89  use interfaces_41_geometry, except_this_one => symatm
 90 !End of the abilint section
 91 
 92  implicit none
 93 
 94 !Arguments ------------------------------------
 95 !scalars
 96  integer,intent(in) :: natom,nsym
 97  real(dp), intent(in) :: tolsym
 98 !arrays
 99  integer,intent(in) :: symrec(3,3,nsym),typat(natom)
100  integer,intent(out) :: indsym(4,nsym,natom)
101  real(dp),intent(in) :: tnons(3,nsym),xred(3,natom)
102 
103 !Local variables-------------------------------
104 !scalars
105  integer :: eatom,errout,iatom,ii,isym,mu
106  real(dp) :: difmax,err
107  character(len=500) :: message
108 !arrays
109  integer :: transl(3)
110  real(dp) :: difmin(3),tratom(3)
111 
112 ! *************************************************************************
113 
114  err=zero
115  errout=0
116 
117  do isym=1,nsym
118    do iatom=1,natom
119 
120      do mu=1,3 ! Apply inverse transformation to original coordinates. Note transpose of symrec.
121        tratom(mu) = dble(symrec(1,mu,isym))*(xred(1,iatom)-tnons(1,isym))&
122 &       +dble(symrec(2,mu,isym))*(xred(2,iatom)-tnons(2,isym))&
123 &       +dble(symrec(3,mu,isym))*(xred(3,iatom)-tnons(3,isym))
124      end do
125 !    
126 !    Find symmetrically equivalent atom
127      call symchk(difmin,eatom,natom,tratom,transl,typat(iatom),typat,xred)
128 !    
129 !    Put information into array indsym: translations and label
130      indsym(1,isym,iatom)=transl(1)
131      indsym(2,isym,iatom)=transl(2)
132      indsym(3,isym,iatom)=transl(3)
133      indsym(4,isym,iatom)=eatom
134 !    
135 !    Keep track of maximum difference between transformed coordinates and
136 !    nearest "target" coordinate
137      difmax=max(abs(difmin(1)),abs(difmin(2)),abs(difmin(3)))
138      err=max(err,difmax)
139 
140      if (difmax>tolsym) then ! Print warnings if differences exceed tolerance
141        write(message, '(3a,i3,a,i6,a,i3,a,a,3es12.4,3a)' )&
142 &       'Trouble finding symmetrically equivalent atoms',ch10,&
143 &       'Applying inv of symm number',isym,' to atom number',iatom,'  of typat',typat(iatom),ch10,&
144 &       'gives tratom=',tratom(1:3),'.',ch10,&
145 &       'This is further away from every atom in crystal than the allowed tolerance.'
146        MSG_WARNING(message)
147 
148        write(message, '(a,3i3,a,a,3i3,a,a,3i3)' ) &
149 &       '  The inverse symmetry matrix is',symrec(1,1:3,isym),ch10,&
150 &       '                                ',symrec(2,1:3,isym),ch10,&
151 &       '                                ',symrec(3,1:3,isym)
152        call wrtout(std_out,message,'COLL')
153        write(message, '(a,3f13.7)' )'  and the nonsymmorphic transl. tnons =',(tnons(mu,isym),mu=1,3)
154 
155        call wrtout(std_out,message,'COLL')
156        write(message, '(a,1p,3e11.3,a,a,i5)' ) &
157 &       '  The nearest coordinate differs by',difmin(1:3),ch10,&
158 &       '  for indsym(nearest atom)=',indsym(4,isym,iatom)
159        call wrtout(std_out,message,'COLL')
160 !      
161 !      Use errout to reduce volume of error diagnostic output
162        if (errout==0) then
163          write(message,'(6a)') ch10,&
164 &         '  This indicates that when symatm attempts to find atoms  symmetrically',ch10, &
165 &         '  related to a given atom, the nearest candidate is further away than some',ch10,&
166 &         '  tolerance.  Should check atomic coordinates and symmetry group input data.'
167          call wrtout(std_out,message,'COLL')
168          errout=1
169        end if
170 
171      end if !difmax>tol
172    end do !iatom
173  end do !isym
174 
175  if (.FALSE.) then
176    do iatom=1,natom
177      write(message, '(a,i0,a)' )' symatm: atom number ',iatom,' is reached starting at atom'
178      call wrtout(std_out,message,'COLL')
179      do ii=1,(nsym-1)/24+1
180        if(natom<100)then
181          write(message, '(1x,24i3)' ) (indsym(4,isym,iatom),isym=1+(ii-1)*24,min(nsym,ii*24))
182        else 
183          write(message, '(1x,24i6)' ) (indsym(4,isym,iatom),isym=1+(ii-1)*24,min(nsym,ii*24))
184        end if
185        call wrtout(std_out,message,'COLL')
186      end do
187    end do
188  end if
189 
190  if (err>tolsym) then
191    write(message, '(1x,a,1p,e14.5,a,e12.4)' )'symatm: maximum (delta t)=',err,' is larger than tol=',tolsym
192    call wrtout(std_out,message,'COLL')
193  end if
194 
195 !Stop execution if error is really big
196  if (err>0.01d0) then
197    write(message,'(5a)')&
198 &   'Largest error (above) is so large (0.01) that either input atomic coordinates (xred)',ch10,&
199 &   'are wrong or space group symmetry data is wrong.',ch10,&
200 &   'Action : correct your input file.'
201    MSG_ERROR(message)
202  end if
203 
204 end subroutine symatm