TABLE OF CONTENTS
ABINIT/irreducible_set_pert [ Functions ]
NAME
irreducible_set_pert
FUNCTION
Determines a set of perturbations that form a basis in that, using symmetry, they can be used to generate all other perturbations that are asked to be calculated (target).
COPYRIGHT
Copyright (C) 1999-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
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. mpert =maximum number of iper natom= number of atoms nsym=number of space group symmetries rfdir(3)=direction for the perturbations rfpert(mpert)=information on the perturbations symrec(3,3,nsym)=3x3 matrices of the group symmetries (reciprocal space) symrel(3,3,nsym)=3x3 matrices of the group symmetries (real space) symq(4,2,nsym)= (integer) three first numbers define the G vector; fourth number is 0 if the q-vector is not preserved, is 1 otherwise second index is one without time-reversal symmetry, two with time-reversal symmetry
OUTPUT
pertsy(3,mpert)= the target perturbation is described by the two last indices (idir, and ipert), the value is 0, 1 or -1, see notes.
NOTES
Output will be in the pertsy array, 0 for non-target perturbations 1 for basis perturbations -1 for perturbations that can be found from basis perturbations
PARENTS
get_npert_rbz,m_dvdb,respfn
CHILDREN
SOURCE
53 #if defined HAVE_CONFIG_H 54 #include "config.h" 55 #endif 56 57 #include "abi_common.h" 58 59 60 subroutine irreducible_set_pert(indsym,mpert,natom,nsym,pertsy,rfdir,rfpert,symq,symrec,symrel) 61 62 use defs_basis 63 use m_profiling_abi 64 65 !This section has been created automatically by the script Abilint (TD). 66 !Do not modify the following lines by hand. 67 #undef ABI_FUNC 68 #define ABI_FUNC 'irreducible_set_pert' 69 !End of the abilint section 70 71 implicit none 72 73 !Arguments ------------------------------- 74 !scalars 75 integer,intent(in) :: mpert,natom,nsym 76 !arrays 77 integer,intent(in) :: indsym(4,nsym,natom),rfdir(3),rfpert(mpert) 78 integer,intent(in) :: symq(4,2,nsym),symrec(3,3,nsym),symrel(3,3,nsym) 79 integer,intent(out) :: pertsy(3,mpert) 80 81 !Local variables ------------------------- 82 !scalars 83 integer :: found,idir1,idisy1,ii,ipert1,ipesy1,isign,isym,itirev,jj 84 !arrays 85 integer :: sym1(3,3) 86 87 ! ********************************************************************* 88 89 !Zero pertsy 90 pertsy(:,:)=0 91 92 do ipert1=1,mpert 93 do idir1=1,3 94 if(rfpert(ipert1)==1.and.rfdir(idir1)==1)then 95 ! write(std_out,*)' for candidate idir =',idir1,' ipert = ',ipert1 96 97 ! Loop on all symmetries, including time-reversal 98 do isym=1,nsym 99 do itirev=1,2 100 isign=3-2*itirev 101 102 if(symq(4,itirev,isym)/=0)then 103 104 found=1 105 106 ! Here select the symmetric of ipert1 107 if(ipert1<=natom)then 108 ipesy1=indsym(4,isym,ipert1) 109 do ii=1,3 110 do jj=1,3 111 sym1(ii,jj)=symrec(ii,jj,isym) 112 end do 113 end do 114 else if(ipert1==(natom+2) .or. ipert1==(natom+6))then 115 ipesy1=ipert1 116 do ii=1,3 117 do jj=1,3 118 sym1(ii,jj)=symrel(ii,jj,isym) 119 end do 120 end do 121 else 122 found=0 123 end if 124 125 ! Now that a symmetric perturbation has been obtained, 126 ! including the expression of the symmetry matrix, see 127 ! if the symmetric perturbations are available 128 if( found==1 ) then 129 130 do idisy1=1,3 131 if(sym1(idir1,idisy1)/=0)then 132 if(pertsy(idisy1,ipesy1)==0)then 133 found=0 134 exit 135 end if 136 end if 137 end do 138 end if 139 140 ! Now, if still found, then it is a symmetric 141 ! of some linear combination of existing perturbations 142 if(found==1)then 143 144 ! DEBUG 145 ! write(std_out,*)' all found ! isym, isign= ',isym,isign 146 ! write(std_out,1010)((sym1(ii,jj),ii=1,3),jj=1,3) 147 ! write(std_out,1010)((sym2(ii,jj),ii=1,3),jj=1,3) 148 ! write(std_out,*)sumr,sumi 149 ! 1010 format(9i4) 150 ! ENDDEBUG 151 152 pertsy(idir1,ipert1)=-1 153 exit ! Exit loop on symmetry operations 154 155 end if 156 157 end if ! End loop on all symmetries + time-reversal 158 end do 159 end do 160 161 ! Now that all symmetries have been examined, 162 ! if still not symmetric of a linear combination 163 ! of basis perturbations, then it is a basis perturbation 164 if(pertsy(idir1,ipert1)/=-1) pertsy(idir1,ipert1)=1 165 ! write(std_out,'(a,3i5)' ) ' irreducible_set_pert :',idir1,ipert1,pertsy(idir1,ipert1) 166 167 end if ! End big loop on all elements 168 end do 169 end do 170 171 end subroutine irreducible_set_pert