TABLE OF CONTENTS
ABINIT/m_angles [ Modules ]
NAME
m_angles
FUNCTION
COPYRIGHT
Copyright (C) 2008-2018 ABINIT group (BA, 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 .
PARENTS
CHILDREN
SOURCE
20 #if defined HAVE_CONFIG_H 21 #include "config.h" 22 #endif 23 24 #include "abi_common.h" 25 26 module m_angles 27 28 use defs_basis 29 use m_errors 30 use m_profiling_abi 31 32 use m_special_funcs, only : factorial 33 34 implicit none 35 36 private
m_angles/create_mlms2jmj [ Functions ]
[ Top ] [ m_angles ] [ Functions ]
NAME
create_mlms2jmj
FUNCTION
For a given angular momentum lcor, give the rotation matrix msml2jmj
INPUTS
lcor= angular momentum
SIDE EFFECTS
mlms2jmj= rotation matrix
PARENTS
CHILDREN
SOURCE
336 subroutine create_mlms2jmj(lcor,mlmstwojmj) 337 338 339 !This section has been created automatically by the script Abilint (TD). 340 !Do not modify the following lines by hand. 341 #undef ABI_FUNC 342 #define ABI_FUNC 'create_mlms2jmj' 343 !End of the abilint section 344 345 implicit none 346 347 !Arguments --------------------------------------------- 348 !scalars 349 integer,intent(in) :: lcor 350 !arrays 351 complex(dpc),intent(out) :: mlmstwojmj(2*(2*lcor+1),2*(2*lcor+1)) 352 353 !Local variables --------------------------------------- 354 !scalars 355 integer :: jc1,jj,jm,ll,ml1,ms1 356 real(dp) :: invsqrt2lp1,xj,xmj 357 character(len=500) :: message 358 !arrays 359 integer, allocatable :: ind_msml(:,:) 360 complex(dpc),allocatable :: mat_mlms2(:,:) 361 362 !********************************************************************* 363 364 !--------------- Built indices + allocations 365 ll=lcor 366 mlmstwojmj=czero 367 ABI_ALLOCATE(ind_msml,(2,-ll:ll)) 368 ABI_ALLOCATE(mat_mlms2,(2*(2*lcor+1),2*(2*lcor+1))) 369 mlmstwojmj=czero 370 jc1=0 371 do ms1=1,2 372 do ml1=-ll,ll 373 jc1=jc1+1 374 ind_msml(ms1,ml1)=jc1 375 end do 376 end do 377 !--------------- built mlmstwojmj 378 !do jj=ll,ll+1 ! the physical value of j are ll-0.5,ll+0.5 379 !xj(jj)=jj-0.5 380 if(ll==0)then 381 message=' ll should not be equal to zero !' 382 MSG_BUG(message) 383 end if 384 jc1=0 385 invsqrt2lp1=one/sqrt(float(2*lcor+1)) 386 do jj=ll,ll+1 387 xj=float(jj)-half 388 do jm=-jj,jj-1 389 xmj=float(jm)+half 390 jc1=jc1+1 391 if(nint(xj+0.5)==ll+1) then 392 if(nint(xmj+0.5)==ll+1) then 393 mlmstwojmj(ind_msml(2,ll),jc1)=1.0 ! J=L+0.5 and m_J=L+0.5 394 else if(nint(xmj-0.5)==-ll-1) then 395 mlmstwojmj(ind_msml(1,-ll),jc1)=1.0 ! J=L+0.5 and m_J=-L-0.5 396 else 397 mlmstwojmj(ind_msml(2,nint(xmj-0.5)),jc1)=invsqrt2lp1*(sqrt(float(ll)+xmj+0.5)) 398 mlmstwojmj(ind_msml(1,nint(xmj+0.5)),jc1)=invsqrt2lp1*(sqrt(float(ll)-xmj+0.5)) 399 end if 400 end if 401 if(nint(xj+0.5)==ll) then 402 mlmstwojmj(ind_msml(1,nint(xmj+0.5)),jc1)=invsqrt2lp1*(sqrt(float(ll)+xmj+0.5)) 403 mlmstwojmj(ind_msml(2,nint(xmj-0.5)),jc1)=-invsqrt2lp1*(sqrt(float(ll)-xmj+0.5)) 404 end if 405 end do 406 end do 407 408 ABI_DEALLOCATE(ind_msml) 409 ABI_DEALLOCATE(mat_mlms2) 410 411 end subroutine create_mlms2jmj
m_angles/create_slm2ylm [ Functions ]
[ Top ] [ m_angles ] [ Functions ]
NAME
create_slm2ylm
FUNCTION
For a given angular momentum lcor, compute slm2ylm
INPUTS
lcor= angular momentum, size of the matrix is 2(2*lcor+1)
OUTPUT
slm2ylm(2lcor+1,2lcor+1) = rotation matrix.
NOTES
useful only in ndij==4
PARENTS
CHILDREN
SOURCE
270 subroutine create_slm2ylm(lcor,slmtwoylm) 271 272 273 !This section has been created automatically by the script Abilint (TD). 274 !Do not modify the following lines by hand. 275 #undef ABI_FUNC 276 #define ABI_FUNC 'create_slm2ylm' 277 !End of the abilint section 278 279 implicit none 280 281 !Arguments --------------------------------------------- 282 !scalars 283 integer,intent(in) :: lcor 284 !arrays 285 complex(dpc),intent(out) :: slmtwoylm(2*lcor+1,2*lcor+1) 286 287 !Local variables --------------------------------------- 288 !scalars 289 integer :: jm,ll,mm,im 290 real(dp),parameter :: invsqrt2=one/sqrt2 291 real(dp) :: onem 292 !arrays 293 294 ! ********************************************************************* 295 296 ll=lcor 297 slmtwoylm=czero 298 do im=1,2*ll+1 299 mm=im-ll-1;jm=-mm+ll+1 300 onem=dble((-1)**mm) 301 if (mm> 0) then 302 slmtwoylm(im,im)= cmplx(onem*invsqrt2,zero,kind=dp) 303 slmtwoylm(jm,im)= cmplx(invsqrt2, zero,kind=dp) 304 end if 305 if (mm==0) then 306 slmtwoylm(im,im)=cone 307 end if 308 if (mm< 0) then 309 slmtwoylm(im,im)= cmplx(zero, invsqrt2,kind=dp) 310 slmtwoylm(jm,im)=-cmplx(zero,onem*invsqrt2,kind=dp) 311 end if 312 end do 313 314 end subroutine create_slm2ylm
m_angles/dbeta [ Functions ]
[ Top ] [ m_angles ] [ Functions ]
NAME
dbeta
FUNCTION
Calculate the rotation matrix d^l_{m{\prim}m}(beta) using Eq. 4.14 of M.E. Rose, Elementary Theory of Angular Momentum, John Wiley & Sons, New-York, 1957
INPUTS
cosbeta= cosinus of beta (=Euler angle) ll= index l mm= index m mp= index m_prime
OUTPUT
dbeta= rotation matrix
NOTES
- This file comes from the file crystal_symmetry.f by N.A.W. Holzwarth and A. Tackett for the code pwpaw - Assume l relatively small so that factorials do not cause roundoff error
PARENTS
setsymrhoij
CHILDREN
SOURCE
185 function dbeta(cosbeta,ll,mp,mm) 186 187 188 !This section has been created automatically by the script Abilint (TD). 189 !Do not modify the following lines by hand. 190 #undef ABI_FUNC 191 #define ABI_FUNC 'dbeta' 192 !End of the abilint section 193 194 implicit none 195 196 !Arguments --------------------------------------------- 197 !scalars 198 integer,intent(in) :: ll,mm,mp 199 real(dp) :: dbeta 200 real(dp),intent(in) :: cosbeta 201 202 !Local variables ------------------------------ 203 !scalars 204 integer,parameter :: mxterms=200 205 integer :: ii,ina,inb,inc,ml,ms 206 real(dp) :: arg,cosbetab2,pref,sinbetab2,sum,tt 207 208 !************************************************************************ 209 dbeta=zero 210 211 !Special cases 212 if (abs(cosbeta-1._dp).lt.tol10) then 213 if (mp.eq.mm) dbeta=1 214 else if (abs(cosbeta+1._dp).lt.tol10) then 215 if (mp.eq.-mm) dbeta=(-1)**(ll+mm) 216 else 217 218 ! General case 219 cosbetab2=sqrt((1+cosbeta)*0.5_dp) 220 sinbetab2=sqrt((1-cosbeta)*0.5_dp) 221 ml=max(mp,mm) 222 ms=min(mp,mm) 223 if (ml.ne.mp) sinbetab2=-sinbetab2 224 tt=-(sinbetab2/cosbetab2)**2 225 pref=sqrt((factorial(ll-ms)*factorial(ll+ml))& 226 & /(factorial(ll+ms)*factorial(ll-ml)))& 227 & /factorial(ml-ms)*(cosbetab2**(2*ll+ms-ml))& 228 & *((-sinbetab2)**(ml-ms)) 229 sum=1._dp 230 arg=1._dp 231 ina=ml-ll 232 inb=-ms-ll 233 inc=ml-ms+1 234 do ii=1,mxterms 235 if (ina.eq.0.or.inb.eq.0) exit 236 arg=(arg*ina*inb*tt)/(ii*inc) 237 sum=sum+arg 238 ina=ina+1 239 inb=inb+1 240 inc=inc+1 241 end do 242 dbeta=pref*sum 243 end if 244 245 end function dbeta
m_angles/mkeuler [ Functions ]
[ Top ] [ m_angles ] [ Functions ]
NAME
mkeuler
FUNCTION
For a given symmetry operation, determines the corresponding Euler angles
INPUTS
rot(3,3)= symmetry matrix
OUTPUT
cosalp= cos(alpha) with alpha=Euler angle 1 cosbeta= cos(beta) with beta =Euler angle 2 cosgam= cos(gamma) with gamma=Euler angle 3 isn= error code (0 if the routine exit normally) sinalp= sin(alpha) with alpha=Euler angle 1 singam= sin(gamma) with gamma=Euler angle 3
NOTES
This file comes from the file crystal_symmetry.f by N.A.W. Holzwarth and A. Tackett for the code pwpaw
PARENTS
setsymrhoij
CHILDREN
SOURCE
77 subroutine mkeuler(rot,cosbeta,cosalp,sinalp,cosgam,singam,isn) 78 79 80 !This section has been created automatically by the script Abilint (TD). 81 !Do not modify the following lines by hand. 82 #undef ABI_FUNC 83 #define ABI_FUNC 'mkeuler' 84 !End of the abilint section 85 86 implicit none 87 88 !Arguments --------------------------------------------- 89 !scalars 90 integer,intent(out) :: isn 91 real(dp),intent(out) :: cosalp,cosbeta,cosgam,sinalp,singam 92 !arrays 93 real(dp),intent(in) :: rot(3,3) 94 95 !Local variables --------------------------------------- 96 !scalars 97 integer :: ier 98 real(dp) :: check,sinbeta 99 character(len=500) :: message 100 101 ! ********************************************************************* 102 103 do isn= -1,1,2 104 cosbeta=real(isn)*rot(3,3) 105 if(abs(1._dp-cosbeta*cosbeta)<tol10) then 106 sinbeta=zero 107 else 108 sinbeta=sqrt(1._dp-cosbeta*cosbeta) 109 end if 110 if (abs(sinbeta).gt.tol10) then 111 cosalp=isn*rot(3,1)/sinbeta 112 sinalp=isn*rot(3,2)/sinbeta 113 cosgam=-isn*rot(1,3)/sinbeta 114 singam=isn*rot(2,3)/sinbeta 115 else 116 cosalp=isn*rot(1,1)/cosbeta 117 sinalp=isn*rot(1,2)/cosbeta 118 cosgam=one 119 singam=zero 120 end if 121 122 ! Check matrix: 123 ier=0 124 check=cosalp*cosbeta*cosgam-sinalp*singam 125 if (abs(check-isn*rot(1,1))>tol8) ier=ier+1 126 check=sinalp*cosbeta*cosgam+cosalp*singam 127 if (abs(check-isn*rot(1,2))>tol8) ier=ier+1 128 check=-sinbeta*cosgam 129 if (abs(check-isn*rot(1,3))>tol8) ier=ier+1 130 check=-cosalp*cosbeta*singam-sinalp*cosgam 131 if (abs(check-isn*rot(2,1))>tol8) ier=ier+1 132 check=-sinalp*cosbeta*singam+cosalp*cosgam 133 if (abs(check-isn*rot(2,2))>tol8) ier=ier+1 134 check=sinbeta*singam 135 if (abs(check-isn*rot(2,3))>tol8) ier=ier+1 136 check=cosalp*sinbeta 137 if (abs(check-isn*rot(3,1))>tol8) ier=ier+1 138 check=sinalp*sinbeta 139 if (abs(check-isn*rot(3,2))>tol8) ier=ier+1 140 if (ier.eq.0) return 141 end do 142 143 isn=0 144 write(message, '(7a)' )& 145 & 'Error during determination of symetries !',ch10,& 146 & 'Action: check your input file:',ch10,& 147 & ' unit cell vectors and/or atoms positions',ch10,& 148 & ' have to be given with a better precision...' 149 MSG_ERROR(message) 150 151 end subroutine mkeuler