TABLE OF CONTENTS
ABINIT/setsymrhoij [ Functions ]
NAME
setsymrhoij
FUNCTION
PAW only Compute rotation matrices expressed in the basis of real spherical harmonics This coefficients are used later to symmetrize rhoij quantities (augmentation occupancies) and other similar quantities.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (NH, FJ, MT) 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 translations for reciprocal space ($\textrm{bohr}^{-1}$) lmax=value of lmax mentioned at the second line of the psp file nsym=number of symmetry elements in space group pawprtvol=control print volume and debugging output for PAW rprimd(3,3)=dimensional primitive translations in real space (bohr) sym(3,3,nsym)=symmetries of group in terms of operations on primitive translations
OUTPUT
zarot(2*lmax+1,2*lmax+1,lmax+1,nsym)=coefficients of the transformation of real spherical harmonics under the symmetry operations
NOTES
Typical use: sym(:,:,:) is symrec(:,:,:) (rotations in reciprocal space) because we need symrel^-1 (=transpose[symrec]) to symmetrize quantities. - This file comes from the file crystal_symmetry.f by N.A.W. Holzwarth and A. Tackett for the code pwpaw - Uses sign & phase convension of M. E. Rose, Elementary Theory of Angular Momentum, John Wiley & Sons,. inc. 1957) zalpha = exp(-i*alpha) zgamma = exp (-i*gamma) - Assumes each transformation can be expressed in terms of 3 Euler angles with or without inversion Reference for evaluation of rotation matrices in the basis of real SH: Blanco M.A., Florez M. and Bermejo M. Journal of Molecular Structure: THEOCHEM, Volume 419, Number 1, 8 December 1997 , pp. 19-27(9) http://www.unioviedo.es/qcg/art/Theochem419-19-ov-BF97-rotation-matrices.pdf
PARENTS
bethe_salpeter,dfpt_looppert,gstate,initberry,initorbmag,respfn screening,sigma,wfk_analyze
CHILDREN
mkeuler,wrtout
SOURCE
59 #if defined HAVE_CONFIG_H 60 #include "config.h" 61 #endif 62 63 #include "abi_common.h" 64 65 subroutine setsymrhoij(gprimd,lmax,nsym,pawprtvol,rprimd,sym,zarot) 66 67 use defs_basis 68 use m_errors 69 use m_profiling_abi 70 71 use m_special_funcs, only : phim 72 use m_angles, only : mkeuler, dbeta 73 74 !This section has been created automatically by the script Abilint (TD). 75 !Do not modify the following lines by hand. 76 #undef ABI_FUNC 77 #define ABI_FUNC 'setsymrhoij' 78 use interfaces_14_hidewrite 79 !End of the abilint section 80 81 implicit none 82 83 !Arguments --------------------------------------------- 84 !scalars 85 integer,intent(in) :: lmax,nsym,pawprtvol 86 !arrays 87 integer,intent(in) :: sym(3,3,nsym) 88 real(dp),intent(in) :: gprimd(3,3),rprimd(3,3) 89 real(dp),intent(out) :: zarot(2*lmax+1,2*lmax+1,lmax+1,nsym) 90 91 !Local variables ------------------------------ 92 !scalars 93 integer :: i1,ii,il,irot,isn,j1,jj,k1,ll,mm,mp 94 real(dp) :: cosalp,cosbeta,cosgam,sinalp,singam 95 character(len=1000) :: message 96 !arrays 97 real(dp) :: prod(3,3),rot(3,3) 98 !************************************************************************ 99 100 DBG_ENTER("COLL") 101 102 if (abs(pawprtvol)>=3) then 103 write(message,'(8a,i4)') ch10,& 104 & ' PAW TEST:',ch10,& 105 & ' ==== setsymrhoij: rotation matrices in the basis ============',ch10,& 106 & ' ==== of real spherical harmonics ============',ch10,& 107 & ' > Number of symmetries (nsym)=',nsym 108 call wrtout(std_out,message,'COLL') 109 end if 110 111 zarot=zero 112 113 do irot=1,nsym 114 115 if (abs(pawprtvol)>=3) then 116 write(message,'(a,i2,a,9i2,a)') ' >For symmetry ',irot,' (',sym(:,:,irot),')' 117 call wrtout(std_out,message,'COLL') 118 end if 119 120 ! === l=0 case === 121 zarot(1,1,1,irot)=one 122 123 ! === l>0 case === 124 if (lmax>0) then 125 ! Calculate the rotations in the cartesian basis 126 rot=zero;prod=zero 127 do k1=1,3 128 do j1=1,3 129 do i1=1,3 130 prod(i1,j1)=prod(i1,j1)+sym(i1,k1,irot)*rprimd(j1,k1) 131 end do 132 end do 133 end do 134 do j1=1,3 135 do i1=1,3 136 do k1=1,3 137 rot(i1,j1)=rot(i1,j1)+gprimd(i1,k1)*prod(k1,j1) 138 end do 139 if(abs(rot(i1,j1))<tol10) rot(i1,j1)=zero 140 end do 141 end do 142 call mkeuler(rot,cosbeta,cosalp,sinalp,cosgam,singam,isn) 143 do ll=1,lmax 144 il=(isn)**ll 145 do mp=-ll,ll 146 jj=mp+ll+1 147 do mm=-ll,ll 148 ii=mm+ll+1 149 150 ! Formula (47) from the paper of Blanco et al 151 zarot(ii,jj,ll+1,irot)=il& 152 & *(phim(cosalp,sinalp,mm)*phim(cosgam,singam,mp)*sign(1,mp)& 153 *(dbeta(cosbeta,ll,abs(mp),abs(mm))& 154 & +(-1._dp)**mm*dbeta(cosbeta,ll,abs(mm),-abs(mp)))*half& 155 & -phim(cosalp,sinalp,-mm)*phim(cosgam,singam,-mp)*sign(1,mm)& 156 *(dbeta(cosbeta,ll,abs(mp),abs(mm))& 157 & -(-1._dp)**mm*dbeta(cosbeta,ll,abs(mm),-abs(mp)))*half) 158 end do 159 end do 160 end do 161 end if ! lmax case 162 163 if (abs(pawprtvol)>=3) then 164 if(lmax>0) then 165 write(message,'(2a,3(3(2x,f7.3),a))') & 166 & ' Rotation matrice for l=1:',ch10,& 167 & (zarot(1,jj,2,irot),jj=1,3),ch10,& 168 & (zarot(2,jj,2,irot),jj=1,3),ch10,& 169 & (zarot(3,jj,2,irot),jj=1,3) 170 call wrtout(std_out,message,'COLL') 171 end if 172 if(lmax>1) then 173 write(message,'(2a,5(5(2x,f7.3),a))') & 174 & ' Rotation matrice for l=2:',ch10,& 175 & (zarot(1,jj,3,irot),jj=1,5),ch10,& 176 & (zarot(2,jj,3,irot),jj=1,5),ch10,& 177 & (zarot(3,jj,3,irot),jj=1,5),ch10,& 178 & (zarot(4,jj,3,irot),jj=1,5),ch10,& 179 & (zarot(5,jj,3,irot),jj=1,5) 180 call wrtout(std_out,message,'COLL') 181 end if 182 if(lmax>2) then 183 write(message,'(2a,7(7(2x,f7.3),a))') & 184 & ' Rotation matrice for l=3:',ch10,& 185 & (zarot(1,jj,4,irot),jj=1,7),ch10,& 186 & (zarot(2,jj,4,irot),jj=1,7),ch10,& 187 & (zarot(3,jj,4,irot),jj=1,7),ch10,& 188 & (zarot(4,jj,4,irot),jj=1,7),ch10,& 189 & (zarot(5,jj,4,irot),jj=1,7),ch10,& 190 & (zarot(6,jj,4,irot),jj=1,7),ch10,& 191 & (zarot(7,jj,4,irot),jj=1,7) 192 call wrtout(std_out,message,'COLL') 193 end if 194 end if 195 196 end do ! isym loop 197 198 DBG_EXIT("COLL") 199 200 end subroutine setsymrhoij