TABLE OF CONTENTS
ABINIT/symptgroup [ Functions ]
NAME
symptgroup
FUNCTION
Derive the name of the point group (+holohedry), from symrel. Warning: might have to change the holohedry hR to hP, if hexagonal axes
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=actual number of symmetries symrel(3,3,nsym)=nsym symmetry operations in real space in terms of primitive translations
OUTPUT
iholohedry=holohedry number ptgroup=symmetry point group
PARENTS
symanal,symbrav
CHILDREN
symdet
SOURCE
33 #if defined HAVE_CONFIG_H 34 #include "config.h" 35 #endif 36 37 #include "abi_common.h" 38 39 subroutine symptgroup(iholohedry,nsym,ptgroup,symrel) 40 41 use defs_basis 42 use m_profiling_abi 43 use m_errors 44 45 !This section has been created automatically by the script Abilint (TD). 46 !Do not modify the following lines by hand. 47 #undef ABI_FUNC 48 #define ABI_FUNC 'symptgroup' 49 use interfaces_41_geometry, except_this_one => symptgroup 50 !End of the abilint section 51 52 implicit none 53 54 !Arguments ------------------------------------ 55 !scalars 56 integer,intent(in) :: nsym 57 integer,intent(out) :: iholohedry 58 character(len=5),intent(out) :: ptgroup 59 !arrays 60 integer,intent(in) :: symrel(3,3,nsym) 61 62 !Local variables------------------------------- 63 !scalars 64 integer :: inversion,iorder,isym 65 character(len=500) :: message 66 !arrays 67 integer :: identity(3,3),matrix(3,3),n_axes(-6:6),trial(3,3) 68 integer,allocatable :: determinant(:),order(:),root_invers(:) 69 character(len=2),allocatable :: ptsym(:) 70 71 !************************************************************************** 72 73 !DEBUG 74 !write(std_out,*)' symptgroup : enter' 75 !do isym=1,nsym 76 !write(std_out,'(i3,2x,9i3)' )isym,symrel(:,:,isym) 77 !end do 78 !ENDDEBUG 79 80 identity(:,:)=0 81 identity(1,1)=1 ; identity(2,2)=1 ; identity(3,3)=1 82 n_axes(:)=0 83 84 ABI_ALLOCATE(determinant,(nsym)) 85 ABI_ALLOCATE(order,(nsym)) 86 ABI_ALLOCATE(ptsym,(nsym)) 87 ABI_ALLOCATE(root_invers,(nsym)) 88 89 !Get the determinant 90 call symdet(determinant,nsym,symrel) 91 92 !Get the order of each the symmetry operation, as well as the maximal order 93 !Also, examine whether each symmetry operation is the inversion, or a root 94 !of the inversion (like -3) 95 !Finally, decide which kind of point symmetry operation it is 96 do isym=1,nsym 97 98 trial(:,:)=identity(:,:) 99 matrix(:,:)=symrel(:,:,isym) 100 order(isym)=0 101 root_invers(isym)=0 102 do iorder=1,6 103 trial=matmul(matrix,trial) 104 if(sum((trial-identity)**2)==0)then 105 order(isym)=iorder 106 exit 107 end if 108 if(sum((trial+identity)**2)==0)then 109 root_invers(isym)=iorder 110 if(iorder==1)inversion=isym 111 end if 112 end do 113 if(order(isym)==0)then 114 write(message, '(a,i0,a)' )' The symmetry operation number',isym,' is not a root of unity' 115 MSG_BUG(message) 116 end if 117 118 ! determinant, order and root_invers are enough to determine the 119 ! kind of symmetry operation 120 ptsym(isym)='no' 121 select case(order(isym)) 122 case(1) 123 ptsym(isym)=' 1' ; n_axes(1)=n_axes(1)+1 124 case(2) 125 if(determinant(isym)== 1)then 126 ptsym(isym)=' 2' ; n_axes(2)=n_axes(2)+1 127 else if(determinant(isym)==-1 .and. root_invers(isym)==1)then 128 ptsym(isym)='-1' ; n_axes(-1)=n_axes(-1)+1 129 else if(determinant(isym)==-1 .and. root_invers(isym)==0)then 130 ptsym(isym)='-2' ; n_axes(-2)=n_axes(-2)+1 131 end if 132 case(3) 133 ptsym(isym)=' 3' ; n_axes(3)=n_axes(3)+1 134 case(4) 135 if(determinant(isym)== 1)then 136 ptsym(isym)=' 4' ; n_axes(4)=n_axes(4)+1 137 else if(determinant(isym)==-1)then 138 ptsym(isym)='-4' ; n_axes(-4)=n_axes(-4)+1 139 end if 140 case(6) 141 if(determinant(isym)== 1)then 142 ptsym(isym)=' 6' ; n_axes(6)=n_axes(6)+1 143 else if(determinant(isym)==-1 .and. root_invers(isym)==3)then 144 ptsym(isym)='-3' ; n_axes(-3)=n_axes(-3)+1 145 else if(determinant(isym)==-1 .and. root_invers(isym)==0)then 146 ptsym(isym)='-6' ; n_axes(-6)=n_axes(-6)+1 147 end if 148 end select 149 150 if(ptsym(isym)=='no')then 151 write(message,'(a,i4,a,a,a,i4,a,a,i4,a,a,i4)' )& 152 & 'The symmetry operation number',isym,' could not be identified',ch10,& 153 & 'order(isym) =',order(isym),ch10,& 154 & 'determinant(isym)=',determinant(isym),ch10,& 155 & 'root_invers(isym)=',root_invers(isym) 156 MSG_BUG(message) 157 end if 158 159 end do 160 161 iholohedry=0 162 if (sum((n_axes-(/0,0,0,0,0,0, 0 ,1,0,0,0,0,0/))**2)==0)then 163 ptgroup=' 1' ; iholohedry=1 164 else if(sum((n_axes-(/0,0,0,0,0,1, 0 ,1,0,0,0,0,0/))**2)==0)then 165 ptgroup=' -1' ; iholohedry=1 166 167 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,1,0,0,0,0/))**2)==0)then 168 ptgroup=' 2' ; iholohedry=2 169 else if(sum((n_axes-(/0,0,0,0,1,0, 0 ,1,0,0,0,0,0/))**2)==0)then 170 ptgroup=' -2' ; iholohedry=2 171 else if(sum((n_axes-(/0,0,0,0,1,1, 0 ,1,1,0,0,0,0/))**2)==0)then 172 ptgroup=' 2/m' ; iholohedry=2 173 174 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,3,0,0,0,0/))**2)==0)then 175 ptgroup=' 222' ; iholohedry=3 176 else if(sum((n_axes-(/0,0,0,0,2,0, 0 ,1,1,0,0,0,0/))**2)==0)then 177 ptgroup=' mm2' ; iholohedry=3 178 else if(sum((n_axes-(/0,0,0,0,3,1, 0 ,1,3,0,0,0,0/))**2)==0)then 179 ptgroup=' mmm' ; iholohedry=3 180 181 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,1,0,2,0,0/))**2)==0)then 182 ptgroup=' 4' ; iholohedry=4 183 else if(sum((n_axes-(/0,0,2,0,0,0, 0 ,1,1,0,0,0,0/))**2)==0)then 184 ptgroup=' -4' ; iholohedry=4 185 else if(sum((n_axes-(/0,0,2,0,1,1, 0 ,1,1,0,2,0,0/))**2)==0)then 186 ptgroup=' 4/m' ; iholohedry=4 187 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,5,0,2,0,0/))**2)==0)then 188 ptgroup=' 422' ; iholohedry=4 189 else if(sum((n_axes-(/0,0,0,0,4,0, 0 ,1,1,0,2,0,0/))**2)==0)then 190 ptgroup=' 4mm' ; iholohedry=4 191 else if(sum((n_axes-(/0,0,2,0,2,0, 0 ,1,3,0,0,0,0/))**2)==0)then 192 ptgroup=' -42m' ; iholohedry=4 193 else if(sum((n_axes-(/0,0,2,0,5,1, 0 ,1,5,0,2,0,0/))**2)==0)then 194 ptgroup='4/mmm' ; iholohedry=4 195 196 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,0,2,0,0,0/))**2)==0)then 197 ptgroup=' 3' ; iholohedry=5 198 else if(sum((n_axes-(/0,0,0,2,0,1, 0 ,1,0,2,0,0,0/))**2)==0)then 199 ptgroup=' -3' ; iholohedry=5 200 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,3,2,0,0,0/))**2)==0)then 201 ptgroup=' 32' ; iholohedry=5 202 else if(sum((n_axes-(/0,0,0,0,3,0, 0 ,1,0,2,0,0,0/))**2)==0)then 203 ptgroup=' 3m' ; iholohedry=5 204 else if(sum((n_axes-(/0,0,0,2,3,1, 0 ,1,3,2,0,0,0/))**2)==0)then 205 ptgroup=' -3m' ; iholohedry=5 206 207 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,1,2,0,0,2/))**2)==0)then 208 ptgroup=' 6' ; iholohedry=6 209 else if(sum((n_axes-(/2,0,0,0,1,0, 0 ,1,0,2,0,0,0/))**2)==0)then 210 ptgroup=' -6' ; iholohedry=6 211 else if(sum((n_axes-(/2,0,0,2,1,1, 0 ,1,1,2,0,0,2/))**2)==0)then 212 ptgroup=' 6/m' ; iholohedry=6 213 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,7,2,0,0,2/))**2)==0)then 214 ptgroup=' 622' ; iholohedry=6 215 else if(sum((n_axes-(/0,0,0,0,6,0, 0 ,1,1,2,0,0,2/))**2)==0)then 216 ptgroup=' 6mm' ; iholohedry=6 217 else if(sum((n_axes-(/2,0,0,0,4,0, 0 ,1,3,2,0,0,0/))**2)==0)then 218 ptgroup=' -62m' ; iholohedry=6 219 else if(sum((n_axes-(/2,0,0,2,7,1, 0 ,1,7,2,0,0,2/))**2)==0)then 220 ptgroup='6/mmm' ; iholohedry=6 221 222 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,3,8,0,0,0/))**2)==0)then 223 ptgroup=' 23' ; iholohedry=7 224 else if(sum((n_axes-(/0,0,0,8,3,1, 0 ,1,3,8,0,0,0/))**2)==0)then 225 ptgroup=' m-3' ; iholohedry=7 226 else if(sum((n_axes-(/0,0,0,0,0,0, 0 ,1,9,8,6,0,0/))**2)==0)then 227 ptgroup=' 432' ; iholohedry=7 228 else if(sum((n_axes-(/0,0,6,0,6,0, 0 ,1,3,8,0,0,0/))**2)==0)then 229 ptgroup=' -43m' ; iholohedry=7 230 else if(sum((n_axes-(/0,0,6,8,9,1, 0 ,1,9,8,6,0,0/))**2)==0)then 231 ptgroup=' m-3m' ; iholohedry=7 232 233 end if 234 235 if(iholohedry==0)then 236 MSG_ERROR_CLASS('Could not find the point group', "TolSymError") 237 end if 238 239 !DEBUG 240 !do isym=1,nsym 241 !write(std_out,'(a,3i5)' )& 242 !& ' symptgroup : isym,determinant,order=',isym,determinant(isym),order(isym) 243 !end do 244 245 !write(std_out,'(a,13i3)' )' symptgroup : n_axes(-6:6)=',n_axes(-6:6) 246 !write(std_out,*)' iholohedry, ptgroup=',iholohedry,',',ptgroup 247 !ENDDEBUG 248 249 ABI_DEALLOCATE(determinant) 250 ABI_DEALLOCATE(order) 251 ABI_DEALLOCATE(ptsym) 252 ABI_DEALLOCATE(root_invers) 253 254 end subroutine symptgroup