TABLE OF CONTENTS


ABINIT/symrelrot [ Functions ]

[ Top ] [ Functions ]

NAME

 symrelrot

FUNCTION

 Transform the symmetry matrices symrel expressed in the coordinate system rprimd,
 to symmetry matrices symrel expressed in the new coordinate system rprimd_new

COPYRIGHT

 Copyright (C) 2000-2018 ABINIT group (XG)
 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

 nsym=number of symmetries
 rprimd(3,3)=dimensional primitive translations for real space (bohr)
 rprimd_new(3,3)=new dimensional primitive translations for real space (bohr)

OUTPUT

  (see side effects)

SIDE EFFECTS

 Input/Output
 symrel(3,3,nsym)=symmetry operations in real space in terms
 of primitive translations rprimd at input and rprimd_new at output

PARENTS

      ingeo,m_esymm,symbrav,symlatt,symspgr

CHILDREN

      matr3inv

SOURCE

 38 #if defined HAVE_CONFIG_H
 39 #include "config.h"
 40 #endif
 41 
 42 #include "abi_common.h"
 43 
 44 
 45 subroutine symrelrot(nsym,rprimd,rprimd_new,symrel,tolsym)
 46 
 47  use defs_basis
 48  use m_errors
 49  use m_profiling_abi
 50 
 51 !This section has been created automatically by the script Abilint (TD).
 52 !Do not modify the following lines by hand.
 53 #undef ABI_FUNC
 54 #define ABI_FUNC 'symrelrot'
 55  use interfaces_32_util
 56 !End of the abilint section
 57 
 58  implicit none
 59 
 60 !Arguments ------------------------------------
 61 !scalars
 62  integer,intent(in) :: nsym
 63  real(dp),intent(in) :: tolsym
 64 !arrays
 65  integer,intent(inout) :: symrel(3,3,nsym)
 66  real(dp),intent(in) :: rprimd(3,3),rprimd_new(3,3)
 67 
 68 !Local variables-------------------------------
 69 !scalars
 70  integer :: ii,isym,jj
 71  real(dp) :: val
 72  character(len=500) :: message
 73 !arrays
 74  real(dp) :: coord(3,3),coordinvt(3,3),matr1(3,3),matr2(3,3),rprimd_invt(3,3)
 75 
 76 !**************************************************************************
 77 
 78 !Compute the coordinates of rprimd_new in the system defined by rprimd(:,:)
 79  call matr3inv(rprimd,rprimd_invt)
 80  do ii=1,3
 81    coord(:,ii)=rprimd_new(1,ii)*rprimd_invt(1,:)+ &
 82 &   rprimd_new(2,ii)*rprimd_invt(2,:)+ &
 83 &   rprimd_new(3,ii)*rprimd_invt(3,:)
 84  end do
 85 
 86 !Transform symmetry matrices in the system defined by rprimd_new
 87  call matr3inv(coord,coordinvt)
 88  do isym=1,nsym
 89    do ii=1,3
 90      matr1(:,ii)=symrel(:,1,isym)*coord(1,ii)+&
 91 &     symrel(:,2,isym)*coord(2,ii)+&
 92 &     symrel(:,3,isym)*coord(3,ii)
 93    end do
 94    do ii=1,3
 95      matr2(:,ii)=coordinvt(1,:)*matr1(1,ii)+&
 96 &     coordinvt(2,:)*matr1(2,ii)+&
 97 &     coordinvt(3,:)*matr1(3,ii)
 98    end do
 99 
100 !  Check that the new symmetry matrices are made of integers, and store them
101    do ii=1,3
102      do jj=1,3
103        val=matr2(ii,jj)
104 !      Need to allow for twice tolsym, in case of centered Bravais lattices (but do it for all lattices ...)
105        if(abs(val-nint(val))>two*tolsym)then
106          write(message,'(2a,a,i3,a,a,3es14.6,a,a,3es14.6,a,a,3es14.6)')&
107 &         'One of the components of symrel is non-integer,',ch10,&
108 &         '  for isym=',isym,ch10,&
109 &         '  symrel=',matr2(:,1),ch10,&
110 &         '         ',matr2(:,2),ch10,&
111 &         '         ',matr2(:,3)
112          MSG_ERROR_CLASS(message, "TolSymError")
113        end if
114        symrel(ii,jj,isym)=nint(val)
115      end do
116    end do
117 !  End loop isym
118  end do
119 
120 end subroutine symrelrot