TABLE OF CONTENTS


ABINIT/chkorthsy [ Functions ]

[ Top ] [ Functions ]

FUNCTION

 Check the orthogonality of the symmetry operations
 (lengths and absolute values of scalar products should be preserved)

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

 gprimd(3,3)=dimensional primitive transl. for reciprocal space (bohr**-1)
 rmet=Real space metric.
 nsym=actual number of symmetries
 rprimd(3,3)=dimensional primitive translations for real space (bohr)
 symrel(3,3,1:nsym)=symmetry operations in real space in terms of primitive translations

SIDE EFFECTS

 iexit= if 0 at input, will do the check, and stop if there is a problem, return 0 if no problem
        if 1 at input, will always input, return 0 if no problem, -1 if there is a problem,
                       also, suppresses printing of problem

PARENTS

      chkinp,ingeo,symmetrize_rprimd

CHILDREN

SOURCE

 36 #if defined HAVE_CONFIG_H
 37 #include "config.h"
 38 #endif
 39 
 40 #include "abi_common.h"
 41 
 42 
 43 subroutine chkorthsy(gprimd,iexit,nsym,rmet,rprimd,symrel)
 44 
 45  use defs_basis
 46  use m_errors
 47  use m_profiling_abi
 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 'chkorthsy'
 53 !End of the abilint section
 54 
 55  implicit none
 56 
 57 !Arguments ------------------------------------
 58 !scalars
 59  integer,intent(in) :: nsym
 60  integer,intent(inout) :: iexit
 61 !arrays
 62  integer,intent(in) :: symrel(3,3,nsym)
 63  real(dp),intent(in) :: gprimd(3,3),rmet(3,3),rprimd(3,3)
 64 
 65 !Local variables-------------------------------
 66 !scalars
 67  integer :: ii,isym,jj
 68  real(dp),parameter :: tol=2.0d-8
 69  real(dp) :: residual,rmet2
 70  character(len=500) :: message
 71 !arrays
 72  real(dp) :: prods(3,3),rmet_sym(3,3),rprimd_sym(3,3)
 73 
 74 ! *************************************************************************
 75 
 76 !  DEBUG
 77 !  write(std_out,*)' chkorthsy : enter '
 78 !  write(std_out,*)rprimd(1:3,1)
 79 !  write(std_out,*)rprimd(1:3,2)
 80 !  write(std_out,*)rprimd(1:3,3)
 81 !  ENDDEBUG
 82 
 83  rmet2=zero
 84  do ii=1,3
 85    do jj=1,3
 86      rmet2=rmet2+rmet(ii,jj)*2
 87    end do
 88  end do
 89 
 90 !Loop over all symmetry operations
 91  do isym=1,nsym
 92 
 93 !  Compute symmetric of primitive vectors under point symmetry operations
 94    do ii=1,3
 95      rprimd_sym(:,ii)=symrel(1,ii,isym)*rprimd(:,1)+&
 96 &     symrel(2,ii,isym)*rprimd(:,2)+&
 97 &     symrel(3,ii,isym)*rprimd(:,3)
 98    end do
 99 
100 !  DEBUG
101    !write(std_out,*)' chkorthsy : isym=',isym
102 !  write(std_out,*)rprimd(1:3,1)
103 !  write(std_out,*)rprimd(1:3,2)
104 !  write(std_out,*)rprimd(1:3,3)
105 !  write(std_out,*)rprimd_sym(1:3,1)
106 !  write(std_out,*)rprimd_sym(1:3,2)
107 !  write(std_out,*)rprimd_sym(1:3,3)
108 !  ENDDEBUG
109 
110 !  If the new lattice is the same as the original one,
111 !  the lengths and angles are preserved
112    do ii=1,3
113      rmet_sym(ii,:)=rprimd_sym(1,ii)*rprimd_sym(1,:)+&
114 &     rprimd_sym(2,ii)*rprimd_sym(2,:)+&
115 &     rprimd_sym(3,ii)*rprimd_sym(3,:)
116    end do
117 
118 !  DEBUG
119 !  write(std_out,*)' rmet :'
120 !  write(std_out,*)rmet(1:3,1)
121 !  write(std_out,*)rmet(1:3,2)
122 !  write(std_out,*)rmet(1:3,3)
123 !  write(std_out,*)rmet_sym(1:3,1)
124 !  write(std_out,*)rmet_sym(1:3,2)
125 !  write(std_out,*)rmet_sym(1:3,3)
126 !  ENDDEBUG
127 
128    residual=zero
129    do ii=1,3
130      do jj=1,3
131        residual=residual+(rmet_sym(ii,jj)-rmet(ii,jj))**2
132      end do
133    end do
134 
135    if(sqrt(residual)>tol*sqrt(rmet2))then
136      if(iexit==0)then
137        write(message, '(a,i5,a,a,a,a,a,es12.4,a,a,a,a,a,a,a)' )&
138 &       'The symmetry operation number',isym,' does not preserve',ch10,&
139 &       'vector lengths and angles.',ch10,&
140 &       'The value of the residual is',residual,'.',ch10,&
141 &       'Action : modify rprim, acell and/or symrel so that',ch10,&
142 &       'vector lengths and angles are preserved.',ch10,&
143 &       'Beware, the tolerance on symmetry operations is very small.'
144        MSG_ERROR(message)
145      else
146        iexit=-1
147      end if
148    end if
149 
150 !  Also, the scalar product of rprimd_sym and gprimd must give integer numbers
151    do ii=1,3
152      prods(ii,:)=rprimd_sym(1,ii)*gprimd(1,:)+ &
153 &     rprimd_sym(2,ii)*gprimd(2,:)+ &
154 &     rprimd_sym(3,ii)*gprimd(3,:)
155    end do
156 
157    do ii=1,3
158      do jj=1,3
159        residual=prods(ii,jj)-anint(prods(ii,jj))
160        if(abs(residual)>tol)then
161          if(iexit==0)then
162            write(message, '(a,i0,a,a,a,a,a,a,a)' )&
163 &           'The symmetry operation number',isym,' generates',ch10,&
164 &           'a different lattice.',ch10,&
165 &           'Action : modify rprim, acell and/or symrel so that',ch10,&
166 &           'the lattice is preserved.'
167            MSG_ERROR(message)
168          else
169            iexit=-1
170          end if
171        end if
172      end do
173    end do
174 
175    if(iexit==-1) exit
176 
177  end do ! isym
178 
179  if(iexit==1)iexit=0
180 
181 end subroutine chkorthsy