TABLE OF CONTENTS
ABINIT/remove_inversion [ Functions ]
NAME
remove_inversion
FUNCTION
Remove the inversion symmetry from a symmetry set as well all the improper rotations (if present)
COPYRIGHT
Copyright (C) 2008-2018 ABINIT group (MG) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
INPUTS
nsym=initial number of symmetries symrel(3,3,nsym)=Initial set of symmetry operarations in real space tnons(3,nsym)=Initial fractional translations
OUTPUT
nsym_out=Number of symmetries in the set without improper rotation symrel_out(:,:) [pointer] = output symmetries without improper rotations tnons_out(:) [pointer] = fractional translations associated to symrel_out pinv=-1 if the inversion has been removed, 1 otherwise
NOTES
Note the use of pointers, memory is allocated inside the procedure and passed back to the caller. Thus memory deallocation is relegated to the caller. To be on the safe side the pointers should be nullified before entering.
PARENTS
m_crystal,m_io_kss,outkss
CHILDREN
set2unit,symdet,wrtout
SOURCE
40 #if defined HAVE_CONFIG_H 41 #include "config.h" 42 #endif 43 44 #include "abi_common.h" 45 46 subroutine remove_inversion(nsym,symrel,tnons,nsym_out,symrel_out,tnons_out,pinv) 47 48 use defs_basis 49 use m_profiling_abi 50 use m_errors 51 52 use m_numeric_tools, only : isinteger, set2unit 53 54 !This section has been created automatically by the script Abilint (TD). 55 !Do not modify the following lines by hand. 56 #undef ABI_FUNC 57 #define ABI_FUNC 'remove_inversion' 58 use interfaces_14_hidewrite 59 use interfaces_41_geometry, except_this_one => remove_inversion 60 !End of the abilint section 61 62 implicit none 63 64 !Arguments ------------------------------------ 65 !scalars 66 integer,intent(in) :: nsym 67 integer,intent(out) :: nsym_out,pinv 68 !arrays 69 integer,intent(in) :: symrel(3,3,nsym) 70 integer,pointer :: symrel_out(:,:,:) 71 real(dp),intent(in) :: tnons(3,nsym) 72 real(dp),pointer :: tnons_out(:,:) 73 74 !Local variables------------------------------- 75 !scalars 76 integer :: is,is2,is_discarded,is_inv,is_retained,nsym2 77 logical :: found 78 character(len=500) :: msg 79 !arrays 80 integer :: determinant(nsym),inversion(3,3),symrel2(3,3,nsym) 81 real(dp) :: dtnons(3),tnons2(3,nsym) 82 83 ! ********************************************************************* 84 85 msg = 'Removing inversion related symmetrie from initial set' 86 MSG_WARNING(msg) 87 ! 88 !=== Find the occurence of the inversion symmetry === 89 call set2unit(inversion) ; inversion=-inversion 90 91 is_inv=0 ; found=.FALSE. 92 do while (is_inv<nsym .and. .not.found) 93 is_inv=is_inv+1 ; found=ALL(symrel(:,:,is_inv)==inversion) 94 end do 95 if (found) then 96 write(msg,'(a,i3)')' The inversion is symmetry operation no. ',is_inv 97 else 98 write(msg,'(a)')' The inversion was not found in the symmetries list.' 99 end if 100 call wrtout(std_out,msg,'COLL') 101 ! 102 !=== Find the symmetries that are related through the inversion symmetry === 103 call symdet(determinant,nsym,symrel) 104 nsym2=0 105 do is=1,nsym-1 106 do is2=is+1,nsym 107 108 dtnons(:)=tnons(:,is2)-tnons(:,is)-tnons(:,is_inv) 109 found=ALL(symrel(:,:,is)==-symrel(:,:,is2)).and.isinteger(dtnons,tol8) 110 111 if (found) then 112 nsym2=nsym2+1 113 ! * Retain symmetries with positive determinant 114 if (ALL(tnons(:,is2)<tol8).and.ALL(tnons(:,is)<tol8)) then 115 is_retained=is2 ; is_discarded=is 116 if (determinant(is)==1) then 117 is_retained=is ; is_discarded=is2 118 end if 119 else if (ALL(tnons(:,is2)<tol8)) then 120 is_retained=is2 ; is_discarded=is 121 else 122 is_retained=is ; is_discarded=is2 123 end if 124 125 symrel2(:,:,nsym2)=symrel(:,:,is_retained) 126 tnons2 (:,nsym2)=tnons (:,is_retained) 127 write(msg,'(a,i3,a,i3,3a,i3,a)')& 128 & ' Symmetry operations no. ',is,' and no. ',is2,& 129 & ' are related through the inversion.',ch10,& 130 & ' Symmetry operation no. ',is_discarded,' will be suppressed.' 131 call wrtout(std_out,msg,'COLL') 132 end if ! found 133 134 end do !is2 135 end do !is 136 137 if (nsym2/=(nsym/2).or.nsym==1) then 138 write(msg,'(a)')' Program uses the original set of symmetries ' 139 call wrtout(std_out,msg,'COLL') 140 nsym_out=nsym 141 ABI_ALLOCATE(symrel_out,(3,3,nsym)) 142 ABI_ALLOCATE(tnons_out,(3,nsym)) 143 symrel_out(:,:,:)=symrel(:,:,1:nsym) 144 tnons_out(:,:)=tnons(:,1:nsym) 145 pinv=1 146 else 147 write(msg,'(a)')' Inversion related operations have been suppressed from symmetries list.' 148 call wrtout(std_out,msg,'COLL') 149 nsym_out=nsym2 150 ABI_ALLOCATE(symrel_out,(3,3,nsym2)) 151 ABI_ALLOCATE(tnons_out,(3,nsym2)) 152 symrel_out(:,:,:)=symrel2(:,:,1:nsym2) 153 tnons_out(:,:)=tnons(:,1:nsym2) 154 pinv=-1 155 end if 156 157 end subroutine remove_inversion