-
Notifications
You must be signed in to change notification settings - Fork 0
/
m_paw_sphharm.F90
3275 lines (2912 loc) · 100 KB
/
m_paw_sphharm.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!!****m* ABINIT/m_paw_sphharm
!! NAME
!! m_paw_sphharm
!!
!! FUNCTION
!! This module contains a set of routines to compute the complex (resp. real)
!! spherical harmonics Ylm (resp. Slm) (and gradients).
!!
!! COPYRIGHT
!! Copyright (C) 2013-2022 ABINIT group (MT, FJ, NH, TRangel)
!! This file is distributed under the terms of the
!! GNU General Public License, see ~abinit/COPYING
!! or http://www.gnu.org/copyleft/gpl.txt .
!!
!! NOTES
!! FOR DEVELOPPERS: in order to preserve the portability of libPAW library,
!! please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt
!!
!! SOURCE
#include "libpaw.h"
MODULE m_paw_sphharm
USE_DEFS
USE_MSG_HANDLING
USE_MEMORY_PROFILING
implicit none
private
!public procedures.
public :: ylmc ! Complex Spherical harmonics for l<=3.
public :: ylmcd ! First derivative of complex Ylm wrt theta and phi up to l<=3
public :: ylm_cmplx ! All (complex) spherical harmonics for lx<=4
public :: initylmr ! Real Spherical Harmonics on a set of vectors
public :: ys ! Matrix element <Yl'm'|Slm>
public :: lxyz ! Matrix element <Yl'm'|L_idir|Ylm>
public :: slxyzs ! Matrix element <Sl'm'|L_idir|Slm>
public :: plm_coeff ! Coefficients depending on Plm used to compute the 2nd der of Ylm
public :: ass_leg_pol ! Associated Legendre Polynomial Plm(x)
public :: plm_d2theta ! d2(Plm (cos(theta)))/d(theta)2 (P_lm= ass. legendre polynomial)
public :: plm_dphi ! m*P_lm(x)/sqrt((1-x^2) (P_lm= ass. legendre polynomial)
public :: plm_dtheta ! -(1-x^2)^1/2*d/dx{P_lm(x)} (P_lm= ass. legendre polynomial)
public :: pl_deriv ! d2(Pl (x)))/d(x)2 where P_l is a legendre polynomial
public :: mkeuler ! For a given symmetry operation, determines the corresponding Euler angles
public :: dble_factorial ! Compute factorial of an integer; returns a double precision real
public :: dbeta ! Calculate the rotation matrix d^l_{m{\prim}m}(beta)
public :: phim ! Computes Phi_m[theta]=Sqrt[2] cos[m theta], if m>0
! Sqrt[2] sin[Abs(m) theta], if m<0
! 1 , if m=0
public :: mat_mlms2jmj ! Change a matrix from the Ylm basis to the J,M_J basis
public :: mat_slm2ylm ! Change a matrix from the Slm to the Ylm basis or from Ylm to Slm
public :: create_slm2ylm ! For a given angular momentum lcor, compute slm2ylm
public :: create_mlms2jmj ! For a given angular momentum lcor, give the rotation matrix msml2jmj
public :: setsym_ylm ! Compute rotation matrices expressed in the basis of real spherical harmonics
public :: gaunt ! Gaunt coeffients for complex Yml
public :: realgaunt ! Compute "real Gaunt coefficients" with "real spherical harmonics"
public :: setnabla_ylm ! Compute rotation matrices expressed in the basis of real spherical harmonics
public :: nablarealgaunt ! Compute the integrals of nablaSlimi.nablaySjmj Slkmk on the PAW spheres
!!***
CONTAINS
!===========================================================
!!***
!!****f* m_paw_sphharm/ylmc
!! NAME
!! ylmc
!!
!! FUNCTION
!! Return a complex spherical harmonic with l <= 3
!!
!! INPUTS
!! il=angular quantum number
!! im=magnetic quantum number
!! kcart=vector in cartesian coordinates defining the value of \theta and \psi
!! where calculate the spherical harmonic
!!
!! OUTPUT
!! ylm= spherical harmonic
!!
!! NOTES
!! Note the use of double precision complex.
!! Case l>3 not implemented.
!!
!! SOURCE
function ylmc(il,im,kcart)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: il,im
complex(dpc) :: ylmc
!arrays
real(dp),intent(in) :: kcart(3)
!Local variables-------------------------------
!scalars
integer,parameter :: LMAX=3
real(dp),parameter :: PPAD=tol8
real(dp) :: cosphi,costh,costhreephi,costwophi,r,rxy,sinphi,sinth,sinthreephi,sintwophi
!complex(dpc) :: new_ylmc
character(len=500) :: msg
! *************************************************************************
if (ABS(im)>ABS(il)) then
write(msg,'(3(a,i0))') 'm is,',im,' however it should be between ',-il,' and ',il
LIBPAW_ERROR(msg)
end if
ylmc = czero
r=SQRT(kcart(1)**2+kcart(2)**2+kcart(3)**2)
if (r<PPAD) r=r+PPAD
!$if (r<tol10) RETURN
rxy=SQRT(kcart(1)**2+kcart(2)**2)
if (rxy<PPAD)rxy=r+PPAD
!
! Determine theta and phi
costh= kcart(3)/r
#if 1
! old buggy coding
sinth= rxy/r
cosphi= kcart(1)/rxy
sinphi= kcart(2)/rxy
#else
sinth=sqrt(abs((one-costh)*(one+costh))) ! abs is needed to prevent very small negative arg
cosphi=one
sinphi=zero
if (sinth>tol10) then
cosphi=kcart(1)/(r*sinth)
sinphi=kcart(2)/(r*sinth)
end if
#endif
costwophi= two*cosphi**2 - one
sintwophi= two*sinphi*cosphi
costhreephi=cosphi*costwophi-sinphi*sintwophi
sinthreephi=cosphi*sintwophi+sinphi*costwophi
select case (il)
case (0)
ylmc= one/SQRT(four_pi)
case (1)
if (ABS(im)==0) then
ylmc = SQRT(three/(four_pi))*costh
else if (ABS(im)==1) then
ylmc = -SQRT(three/(8._dp*pi))*sinth*CMPLX(cosphi,sinphi)
else
msg='wrong im'
LIBPAW_ERROR(msg)
end if
case (2)
if (ABS(im)==0) then
ylmc = SQRT(5.d0/(16.d0*pi))*(three*costh**2-one)
else if (ABS(im)==1) then
ylmc = -SQRT(15.d0/(8.d0*pi))*sinth*costh*cmplx(cosphi,sinphi)
else if (ABS(im)==2) then
ylmc = SQRT(15.d0/(32.d0*pi))*(sinth)**2*CMPLX(costwophi,sintwophi)
else
msg='wrong im'
LIBPAW_ERROR(msg)
end if
case (3)
if (ABS(im)==0) then
ylmc= SQRT(7.d0/(16.d0*pi))*(5.d0*costh**3 -3.d0*costh)
else if (ABS(im)==1) then
ylmc= -SQRT(21.d0/(64.d0*pi))*sinth*(5.d0*costh**2-one)*CMPLX(cosphi,sinphi)
else if (ABS(im)==2) then
ylmc= SQRT(105.d0/(32.d0*pi))*sinth**2*costh*CMPLX(costwophi,sintwophi)
else if (ABS(im)==3) then
ylmc=-SQRT(35.d0/(64.d0*pi))*sinth**3*CMPLX(costhreephi,sinthreephi)
else
msg='wrong im'
LIBPAW_ERROR(msg)
end if
case default
!write(msg,'(a,i6,a,i6)')' The maximum allowed value for l is,',LMAX,' however l=',il
!LIBPAW_ERROR(msg)
end select
!
!=== Treat the case im < 0 ===
if (im < 0) ylmc=(-one)**(im)*CONJG(ylmc)
! FIXME: Use the piece of code below as it works for arbitrary (l,m)
! the implementation above is buggy when the vector is along z!
!
#if 0
! Remember the expression of complex spherical harmonics:
! $Y_{lm}(\theta,\phi)=sqrt{{(2l+1) over (4\pi)} {fact(l-m)/fact(l+m)} } P_l^m(cos(\theta)) e^{i m\phi}$
new_ylmc = SQRT((2*il+1)*dble_factorial(il-ABS(im))/(dble_factorial(il+ABS(im))*four_pi)) * &
& ass_leg_pol(il,ABS(im),costh) * CMPLX(cosphi,sinphi)**ABS(im)
if (im<0) new_ylmc=(-one)**(im)*CONJG(new_ylmc)
if (ABS(new_ylmc-ylmc)>tol6) then
!LIBPAW_WARNING("Check new_ylmc")
!write(std_out,*)"il,im,new_ylmc, ylmc",il,im,new_ylmc,ylmc
!write(std_out,*)"fact",SQRT((2*il+1)*dble_factorial(il-ABS(im))/(dble_factorial(il+ABS(im))*four_pi))
!write(std_out,*)"costh,sinth,ass_leg_pol",costh,sinth,ass_leg_pol(il,ABS(im),costh)
!write(std_out,*)"cosphi,sinphi,e^{imphi}",cosphi,sinphi,CMPLX(cosphi,sinphi)**ABS(im)
end if
ylmc = new_ylmc
#endif
end function ylmc
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/ylmcd
!! NAME
!! ylmcd
!!
!! FUNCTION
!! Computes dth and dphi, the first derivatives of complex Ylm as a function of
!! th and phi (the angles of the spherical coordinates)
!! It works for all spherical harmonics with l <= 3
!!
!! INPUTS
!! il=angular quantum number
!! im=magnetic quantum number
!! kcart=cartesian coordinates of the vector where the first derivatives of Ylm are evaluated
!!
!! OUTPUT
!! dth =derivative of Y_lm with respect to \theta
!! dphi=derivative of Y_lm with respect to \phi
!!
!! NOTES
!! Note the use of double precision complex.
!! Case l>3 not implemented.
!!
!! SOURCE
subroutine ylmcd(il,im,kcart,dth,dphi)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: il,im
complex(dpc),intent(out) :: dphi,dth
!arrays
real(dp),intent(in) :: kcart(3)
!Local variables-------------------------------
!scalars
integer,parameter :: LMAX=3
real(dp),parameter :: PPAD=tol8
real(dp) :: cosphi,costh,costhreephi,costwophi,r,rxy,sinphi,sinth,sinthreephi,sintwophi,c
character(len=500) :: msg
! *************************************************************************
if (ABS(im)>ABS(il))then
write(msg,'(3(a,i0))')' m is,',im,' however it should be between ',-il,' and ',il
LIBPAW_ERROR(msg)
end if
dphi=czero; dth=czero
r=SQRT(kcart(1)**2+kcart(2)**2+kcart(3)**2)
if (r<PPAD) r=r+PPAD
!$if (r<tol10) RETURN
rxy=SQRT(kcart(1)**2+kcart(2)**2)
if (rxy<PPAD) rxy=r+PPAD
! Determine theta and phi
costh= kcart(3)/r
#if 1
! old buggy coding
sinth= rxy/r
cosphi= kcart(1)/rxy
sinphi= kcart(2)/rxy
#else
sinth=sqrt(abs((one-costh)*(one+costh))) ! abs is needed to prevent very small negative arg
cosphi=one
sinphi=zero
if (sinth>tol10) then
cosphi=kcart(1)/(r*sinth)
sinphi=kcart(2)/(r*sinth)
end if
#endif
costwophi= two*cosphi**2 - one
sintwophi= two*sinphi*cosphi
costhreephi=cosphi*costwophi-sinphi*sintwophi
sinthreephi=cosphi*sintwophi+sinphi*costwophi
select case (il)
case (0)
dth = czero
dphi = czero
case (1)
if (ABS(im)==0) then
dth= -SQRT(three/(four_pi))*sinth
dphi= czero
else if (abs(im)==1) then
dth= -SQRT(3.d0/(8.d0*pi))*costh*CMPLX(cosphi,sinphi)
dphi=-SQRT(3.d0/(8.d0*pi))*sinth*CMPLX(-sinphi,cosphi)
end if
case (2)
if (ABS(im)==0) then
dth= -SQRT(5.d0/(16.d0*pi))*6.d0*costh*sinth
dphi= czero
else if (ABS(im)==1) then
dth= -SQRT(15.d0/(8.d0*pi))*(costh**2-sinth**2)*CMPLX(cosphi,sinphi)
dphi= -SQRT(15.d0/(8.d0*pi))*costh*sinth*(0.d0,1.d0)*CMPLX(cosphi,sinphi)
else if (abs(im)==2) then
dth = SQRT(15.d0/(32.d0*pi))*2.d0*costh*sinth*CMPLX(costwophi,sintwophi)
dphi = SQRT(15.d0/(32.d0*pi))*sinth**2*(0.d0,2.d0)*CMPLX(costwophi,sintwophi)
end if
case (3)
if (ABS(im)==0) then
dth = SQRT(7.d0/(16*pi))*(-15.d0*costh**2*sinth + 3.d0**sinth)
dphi= czero
else if (ABS(im)==1) then
c = SQRT(21.d0/(64.d0*pi))
dth= -c* (15.d0*costh**3-11.d0*costh)* CMPLX(cosphi,sinphi)
dphi=-c*sinth*( 5.d0*costh**2-1 )*(0.d0,1.d0)*CMPLX(cosphi,sinphi)
else if (ABS(im)==2) then
c = SQRT(105.d0/(32.d0*pi))
dth =c*(2.d0*sinth*costh**2-sinth**3) *CMPLX(costwophi,sintwophi)
dphi=c*(2.d0*sinth**2*costh)*(0.d0,1.d0)*CMPLX(costwophi,sintwophi)
else if (abs(im)==3) then
dth =-SQRT(35.d0/(64.d0*pi))*3.d0*sinth**2*costh*CMPLX(costhreephi,sinthreephi)
dphi=-SQRT(35.d0/(64.d0*pi))*sinth**3*(0.d0,3.d0)*CMPLX(costhreephi,sinthreephi)
end if
case default
write(msg,'(2(a,i0))')' The maximum allowed value for l is,',LMAX,' however, l=',il
LIBPAW_ERROR(msg)
end select
!
!=== Treat the case im < 0 ===
if (im<0) then
dth = (-one)**(im)*CONJG(dth)
dphi= (-one)**(im)*CONJG(dphi)
end if
end subroutine ylmcd
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/ylm_cmplx
!! NAME
!! ylm_cmplx
!!
!! FUNCTION
!! Calculate all (complex) spherical harmonics for lx<=4
!!
!! INPUTS
!! lx= quantum numbers.
!! xx= cartesian coordinate in the x direction
!! yy= cartesian coordinate in the y direction
!! zz= cartesian coordinate in the z direction
!!
!! cartesian coordinates
!! OUTPUT
!! ylm((lx+1)*(lx+1)) complex spherical harmonics for all l<=lx and all
!! possible values of m.
!!
!! NOTES
!! We are supressing the so-called Condon-Shortley phase
!!
!! SOURCE
subroutine ylm_cmplx(lx,ylm,xx,yy,zz)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: lx
real(dp),intent(in) :: xx,yy,zz
!arrays
complex(dpc),intent(out) :: ylm((lx+1)*(lx+1))
!Local variables-------------------------------
!scalars
integer :: ii,l1,m1,nc,nn
real(dp) :: dc,dl,dm,ds,rr,rrs,rs,sq2,w,x,xs,ya,yi,yr
!arrays
real(dp) :: cosa(lx+1),fact(2*(lx+1)),plm(lx+2,lx+2),qlm(lx+2,lx+2),sgn(lx+1)
real(dp) :: sina(lx+1)
! *************************************************************************
!normalization coefficients
sq2=sqrt(2.0d0)
fact(1)=1.0d0
do ii=2,2*lx+1
fact(ii)=(ii-1)*fact(ii-1)
end do
do l1=1,lx+1
sgn(l1)=(-1.d0)**(l1-1)
do m1=1,l1
qlm(l1,m1)=sqrt((2*l1-1)*fact(l1-m1+1)/&
& (four_pi*fact(l1+m1-1)))
end do
end do
!legendre polynomials
rs=xx**2 + yy**2 + zz**2
if(rs > tol8) then
xs=zz**2/rs
x=zz/sqrt(rs)
w=sqrt(abs(1.0d0 - xs))
else
x=0.0d0
w=1.0d0
end if
plm(1,1)=1.0d0
plm(2,1)=x
plm(2,2)=w
plm(3,2)=3.0d0*x*w
do m1=1,lx
dm=m1-1
if(m1 > 1) then
plm(m1+1,m1)=x*plm(m1,m1) + 2*dm*w*plm(m1,m1-1)
end if
if(m1 < lx) then
do l1=m1+2,lx+1
dl=l1-1
plm(l1,m1)=((2*dl-1)*x*plm(l1-1,m1)&
& - (dl+dm-1)*plm(l1-2,m1))/(dl-dm)
end do
end if
plm(m1+1,m1+1)=(2*dm+1)*w*plm(m1,m1)
end do
!azimuthal angle phase factors
rrs=xx**2 + yy**2
if(rrs > tol8) then
rr=sqrt(rrs)
dc=xx/rr
ds=yy/rr
else
dc=1.0d0
ds=0.0d0
end if
cosa(1)=1.0d0
sina(1)=0.0d0
do m1=2,lx+1
cosa(m1)=dc*cosa(m1-1) - ds*sina(m1-1)
sina(m1)=ds*cosa(m1-1) + dc*sina(m1-1)
end do
!combine factors
do l1=1,lx+1
do m1=2,l1
nn=(l1-1)**2 + (l1-1) + (m1-1) + 1
nc=(l1-1)**2 + (l1-1) - (m1-1) + 1
! note that we are supressing the so-called Condon-Shortley phase
! ya=sgn(m1)*qlm(l1,m1)*plm(l1,m1)
ya=qlm(l1,m1)*plm(l1,m1)
yr=ya*cosa(m1)
yi=ya*sina(m1)
ylm(nc)=sgn(m1)*cmplx(yr,-yi)
ylm(nn)=cmplx(yr,yi)
end do
end do
do l1=1,lx+1
nn=(l1-1)**2 + (l1-1) + 1
ya=qlm(l1,1)*plm(l1,1)
ylm(nn)=cmplx(ya,0.d0)
end do
end subroutine ylm_cmplx
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/initylmr
!! NAME
!! initylmr
!!
!! FUNCTION
!! Calculate the real spherical harmonics Ylm (and gradients)
!! over a set of (r) vectors given in Cartesian coordinates.
!!
!! INPUTS
!! mpsang= 1+maximum angular momentum for nonlocal pseudopotential
!! normchoice=0 the input rr vectors are normalized
!! =1 the norm of the input vector is in nrm() array
!! nrm(npts) = Depending of normchoice, this array contains
!! either the weight of the point or the norm of rr.
!! npts = number of rr vectors
!! option= 1=compute Ylm(R), 2=compute Ylm(R) and dYlm/dRi (cartesian derivatives),
!! 3=compute Ylm(R), dYlm/dRi and d2Ylm/dRidRj (cartesian derivatives)
!! rr(3,npts)= vectors for which ylmr have to be calculated
!! For each point of the spherical mesh, gives the
!! Cartesian coordinates of the corresponding point.
!!
!! OUTPUT
!! if (option=1, 2 or 3)
!! ylm(mpsang*mpsang,npts) = real spherical harmonics for each r point
!! if (option=2 or 3)
!! ylmr_gr(1:3,mpsang*mpsang,npts)= gradients of real spherical harmonics
!! if (option=3)
!! ylmr_gr(4:9,mpsang*mpsang,npts)= first and second gradients of real spherical harmonics
!!
!! NOTES
!! Remember the expression of complex spherical harmonics:
!! $Y_{lm}(%theta ,%phi)=sqrt{{(2l+1) over (4 %pi)} {fact(l-m) over fact(l+m)} } P_l^m(cos(%theta)) func e^{i m %phi}$
!! Remember the expression of real spherical harmonics as linear combination of imaginary spherical harmonics:
!! $Yr_{lm}(%theta ,%phi)=(Re{Y_{l-m}}+(-1)^m Re{Y_{lm}})/sqrt{2}
!! $Yr_{l-m}(%theta ,%phi)=(Im{Y_{l-m}}-(-1)^m Im{Y_{lm}})/sqrt{2}
!!
!! SOURCE
subroutine initylmr(mpsang,normchoice,npts,nrm,option,rr,ylmr,ylmr_gr)
!Arguments ------------------------------------
!scalars
integer,intent(in) :: mpsang,normchoice,npts,option
!arrays
real(dp),intent(in) :: nrm(npts),rr(3,npts)
real(dp),intent(out) :: ylmr(mpsang*mpsang,npts)
real(dp),optional,intent(out) :: ylmr_gr(3*(option/2)+6*(option/3),mpsang*mpsang,npts)
!Local variables ------------------------------
!scalars
integer :: dimgr,ilang,inpt,l0,ll,mm
real(dp) :: cphi,ctheta,fact,onem,rnorm,sphi,stheta,work1,work2,ylmcst,ylmcst2
logical :: compute_ylm,compute_ylm2gr,compute_ylmgr
!arrays
real(dp) :: dphi(3),dtheta(3),iphase(mpsang-1),rphase(mpsang-1)
real(dp),allocatable :: blm(:,:)
!************************************************************************
!What has to be computed ?
compute_ylm = (option==1.or.option==2.or.option==3)
compute_ylmgr =(( option==2.or.option==3).and.present(ylmr_gr))
compute_ylm2gr=(( option==3).and.present(ylmr_gr))
dimgr=3*(option/2)+6*(option/3)
!Initialisation of spherical harmonics
if (compute_ylm ) ylmr (: ,1:npts)=zero
if (compute_ylmgr) ylmr_gr(:,:,1:npts)=zero
!Special case for l=0
if (compute_ylm ) ylmr(1,1:npts)=1._dp/sqrt(four_pi)
if (compute_ylmgr) ylmr_gr(1:dimgr,1,1:npts)=zero
if (mpsang>1) then
! Loop over all rr
do inpt=1,npts
! Load module of rr
rnorm=one
if (normchoice==1) rnorm=nrm(inpt)
! Continue only for r<>0
if (rnorm>tol10) then
! Determine theta and phi
cphi=one
sphi=zero
ctheta=rr(3,inpt)/rnorm
! MM030519 : abs is needed to prevent very small negative arg
stheta=sqrt(abs((one-ctheta)*(one+ctheta)))
if (stheta>tol10) then
cphi=rr(1,inpt)/(rnorm*stheta)
sphi=rr(2,inpt)/(rnorm*stheta)
end if
do mm=1,mpsang-1
rphase(mm)=dreal(dcmplx(cphi,sphi)**mm)
iphase(mm)=aimag(dcmplx(cphi,sphi)**mm)
end do
! Determine gradients of theta and phi
if (compute_ylmgr) then
dtheta(1)=ctheta*cphi
dtheta(2)=ctheta*sphi
dtheta(3)=-stheta
dphi(1)=-sphi
dphi(2)=cphi
dphi(3)=zero
end if
! COMPUTE Ylm(R)
if (compute_ylm) then
! Loop over angular momentum l
do ilang=2,mpsang
ll=ilang-1
l0=ll**2+ll+1
fact=1._dp/real(ll*(ll+1),dp)
ylmcst=sqrt(real(2*ll+1,dp)/four_pi)
! Special case m=0
ylmr(l0,inpt)=ylmcst*ass_leg_pol(ll,0,ctheta)
! Compute for m>0
onem=one
do mm=1,ll
onem=-onem
work1=ylmcst*sqrt(fact)*onem*ass_leg_pol(ll,mm,ctheta)*sqrt(2._dp)
ylmr(l0+mm,inpt)=work1*rphase(mm)
ylmr(l0-mm,inpt)=work1*iphase(mm)
if (mm/=ll) fact=fact/real((ll+mm+1)*(ll-mm),dp)
end do ! End loop over m
end do ! End loop over l
end if
! COMPUTE dYlm/dRi
if (compute_ylmgr) then
! Loop over angular momentum l
do ilang=2,mpsang
ll=ilang-1
l0=ll**2+ll+1
fact=1._dp/real(ll*(ll+1),dp)
ylmcst=sqrt(real(2*ll+1,dp)/four_pi)/rnorm
! Special case m=0
work1=ylmcst*plm_dtheta(ll,0,ctheta)
ylmr_gr(1:3,l0,inpt)=work1*dtheta(1:3)
! Compute for m>0
onem=one
do mm=1,ll
onem=-onem
work1=ylmcst*sqrt(fact)*onem*plm_dtheta(ll,mm,ctheta)*sqrt(2._dp)
work2=ylmcst*sqrt(fact)*onem*plm_dphi (ll,mm,ctheta)*sqrt(2._dp)
ylmr_gr(1:3,l0+mm,inpt)=rphase(mm)*work1*dtheta(1:3)-iphase(mm)*work2*dphi(1:3)
ylmr_gr(1:3,l0-mm,inpt)=iphase(mm)*work1*dtheta(1:3)+rphase(mm)*work2*dphi(1:3)
if (mm/=ll) fact=fact/real((ll+mm+1)*(ll-mm),dp)
end do ! End loop over m
end do ! End loop over l
end if
! COMPUTE d2Ylm/dRidRj
if (compute_ylm2gr) then
LIBPAW_ALLOCATE(blm,(5,mpsang*mpsang))
call plm_coeff(blm,mpsang,ctheta)
! Loop over angular momentum l
do ilang=2,mpsang
ll=ilang-1
l0=ll**2+ll+1
fact=1._dp/real(ll*(ll+1),dp)
ylmcst=sqrt(real(2*ll+1,dp)/four_pi)/(rnorm**2)
! Special case m=0
ylmr_gr(4,l0,inpt)=ylmcst*(-blm(3,l0)*sphi*sphi+blm(4,l0)*cphi*cphi)
ylmr_gr(5,l0,inpt)=ylmcst*(-blm(3,l0)*cphi*cphi+blm(4,l0)*sphi*sphi)
ylmr_gr(6,l0,inpt)=ylmcst*blm(1,l0)
ylmr_gr(7,l0,inpt)=ylmcst*blm(2,l0)*sphi
ylmr_gr(8,l0,inpt)=ylmcst*blm(2,l0)*cphi
ylmr_gr(9,l0,inpt)=ylmcst*(blm(3,l0)+blm(4,l0))*sphi*cphi
! Compute for m>0
onem=one
do mm=1,ll
onem=-onem;ylmcst2=ylmcst*sqrt(fact)*sqrt(two)
ylmr_gr(4,l0+mm,inpt)=ylmcst2*((-blm(3,l0+mm)*sphi*sphi+blm(4,l0+mm)*cphi*cphi)*rphase(mm)-&
& blm(5,l0+mm)*2.d0*cphi*sphi*mm*iphase(mm))
ylmr_gr(4,l0-mm,inpt)=ylmcst2*((-blm(3,l0+mm)*sphi*sphi+blm(4,l0+mm)*cphi*cphi)*iphase(mm)+&
& blm(5,l0+mm)*2.d0*cphi*sphi*mm*rphase(mm))
ylmr_gr(5,l0+mm,inpt)=ylmcst2*((-blm(3,l0+mm)*cphi*cphi+blm(4,l0+mm)*sphi*sphi)*rphase(mm)+&
& blm(5,l0+mm)*2.d0*cphi*sphi*mm*iphase(mm))
ylmr_gr(5,l0-mm,inpt)=ylmcst2*((-blm(3,l0+mm)*cphi*cphi+blm(4,l0+mm)*sphi*sphi)*iphase(mm)-&
& blm(5,l0+mm)*2.d0*cphi*sphi*mm*rphase(mm))
ylmr_gr(6,l0+mm,inpt)=ylmcst2*blm(1,l0+mm)*rphase(mm)
ylmr_gr(6,l0-mm,inpt)=ylmcst2*blm(1,l0+mm)*iphase(mm)
ylmr_gr(7,l0+mm,inpt)=ylmcst2*(blm(2,l0+mm)*sphi*rphase(mm)+&
& mm*iphase(mm)*cphi*onem*plm_dtheta(ll,mm,ctheta))
ylmr_gr(7,l0-mm,inpt)=ylmcst2*(blm(2,l0+mm)*sphi*iphase(mm)-&
& mm*rphase(mm)*cphi*onem*plm_dtheta(ll,mm,ctheta))
ylmr_gr(8,l0+mm,inpt)=ylmcst2*(blm(2,l0+mm)*cphi*rphase(mm)-&
& mm*iphase(mm)*sphi*onem*plm_dtheta(ll,mm,ctheta))
ylmr_gr(8,l0-mm,inpt)=ylmcst2*(blm(2,l0+mm)*cphi*iphase(mm)+&
& mm*rphase(mm)*sphi*onem*plm_dtheta(ll,mm,ctheta))
ylmr_gr(9,l0+mm,inpt)=ylmcst2*((blm(3,l0+mm)+blm(4,l0+mm))*sphi*cphi*rphase(mm)-&
& blm(5,l0+mm)*(sphi*sphi-cphi*cphi)*mm*iphase(mm))
ylmr_gr(9,l0-mm,inpt)=ylmcst2*((blm(3,l0+mm)+blm(4,l0+mm))*sphi*cphi*iphase(mm)+&
& blm(5,l0+mm)*(sphi*sphi-cphi*cphi)*mm*rphase(mm))
if (mm/=ll) fact=fact/real((ll+mm+1)*(ll-mm),dp)
end do ! End loop over m
end do ! End loop over l
LIBPAW_DEALLOCATE(blm)
end if
! End condition r<>0
end if
! End loop over rr
end do
! End condition l<>0
end if
end subroutine initylmr
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/ys
!! NAME
!! ys
!!
!! FUNCTION
!! Computes the matrix element <Y_(l2,m2)|S_(l1,m1)>
!!
!! INPUTS
!! integer :: l2,m2,l1,m1
!!
!! OUTPUT
!! complex(dpc) :: ys_val
!!
!! NOTES
!! Ylm is the standard complex-valued spherical harmonic, Slm is the real spherical harmonic
!! used througout abinit.
!!
!! SOURCE
subroutine ys(l2,m2,l1,m1,ys_val)
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: l1,l2,m1,m2
complex(dpc),intent(out) :: ys_val
!Local variables ---------------------------------------
!scalars
integer :: mp1
! *********************************************************************
! See Blanco et al., J. Mol Struct. 419, 19-27 (1997) Eq. 19
! <Y_l2,m2|S_l1,m1> is given by C^l_{m1,m2} where
! l1 == l2 and |m1| == |m2|, 0 otherwise
ys_val = czero
if ( l2 /= l1 ) return
if ( abs(m2) /= abs(m1) ) return
mp1=(-1)**abs(m1)
if(m1.EQ.0) then
ys_val=cone
else if((m1.GT.0).AND.(m2.GT.0)) then
ys_val=mp1*sqrthalf
else if((m1.GT.0).AND.(m2.LT.0)) then
ys_val=sqrthalf
else if((m1.LT.0).AND.(m2.GT.0)) then
ys_val=-j_dpc*mp1*sqrthalf
else if((m1.LT.0).AND.(m2.LT.0)) then
ys_val=j_dpc*sqrthalf
else
ys_val=czero
end if
end subroutine ys
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/lxyz.F90
!! NAME
!! lxyz
!!
!! FUNCTION
!! Computes the matrix element <Yl'm'|L_idir|Ylm>
!!
!! INPUTS
!! integer :: lp,mp,idir,ll,mm
!!
!! OUTPUT
!! complex(dpc) :: lidir
!!
!! NOTES
!! Ylm is the standard complex-valued spherical harmonic,
!! idir is the direction in space of L
!!
!! SOURCE
subroutine lxyz(lp,mp,idir,ll,mm,lidir)
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: idir,ll,lp,mm,mp
complex(dpc),intent(out) :: lidir
!Local variables ---------------------------------------
!scalars
complex(dpc) :: jme, jmme, jpme
! *********************************************************************
lidir = czero
if ( lp /= ll ) return
jpme=czero; jmme=czero; jme=czero
if (mp==mm) then
jme=cone*mm
else if (mp==mm+1) then
jpme=-cone*sqrt(half*((ll*(ll+1))-mm*(mm+1)))
else if (mp==mm-1) then
jmme= cone*sqrt(half*((ll*(ll+1))-mm*(mm-1)))
end if
select case (idir)
case (1) ! Lx
lidir = -sqrthalf*(jpme - jmme)
case (2) ! Ly
lidir = j_dpc*sqrthalf*(jpme + jmme)
case (3) ! Lz
lidir = jme
end select
end subroutine lxyz
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/slxyzs
!! NAME
!! slxyzs
!!
!! FUNCTION
!! computes the matrix element <Sl'm'|L_idir|Slm>
!!
!! INPUTS
!! integer :: lp,mp,idir,ll,mm
!!
!! OUTPUT
!! complex(dpc) :: sls_val
!!
!! NOTES
!! Slm is the real spherical harmonic used througout abinit,
!! L_idir is a component of the angular momentum operator.
!! The subroutine computes <S_l'm'|L_idir|S_lm>
!!
!! SOURCE
subroutine slxyzs(lp,mp,idir,ll,mm,sls_val)
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: idir,ll,lp,mm,mp
complex(dpc),intent(out) :: sls_val
!Local variables ---------------------------------------
!scalars
integer :: mpp,mppp
complex(dpc) :: lidir,sy_val,ys_val
! *********************************************************************
sls_val = czero
if ( lp /= ll ) return
do mpp = -ll, ll
call ys(ll,mpp,ll,mp,sy_val)
do mppp = -ll, ll
call lxyz(ll,mpp,idir,ll,mppp,lidir)
call ys(ll,mppp,ll,mm,ys_val)
sls_val = sls_val + conjg(sy_val)*lidir*ys_val
end do
end do
end subroutine slxyzs
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/plm_coeff
!! NAME
!! plm_coeff
!!
!! FUNCTION
!! Compute coefficients depending on Plm and its derivatives where P_lm is a legendre polynomial.
!! They are used to compute the second derivatives of spherical harmonics
!!
!! INPUTS
!! mpsang=1+ maximum l quantum number
!! xx= input value
!!
!! OUTPUT
!! blm(5,mpsang*mpsang)=coefficients depending on Plm and its derivatives where P_lm is a legendre polynome
!!
!! SOURCE
subroutine plm_coeff(blm,mpsang,xx)
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: mpsang
real(dp),intent(in) :: xx
!arrays
real(dp),intent(out) :: blm(5,mpsang*mpsang)
!Local variables ---------------------------------------
!scalars
integer :: il,ilm,ilm0,ilm1,im
real(dp) :: dplm_dt,d2plm_dt2,llp1,onemx2,plm,sqrx,xsqrx,xx2,yy
logical :: is_one
character(len=500) :: msg
!arrays
real(dp) :: pl_d2(mpsang),plm_d2t(mpsang*mpsang)
!************************************************************************
if (abs(xx).gt.1.d0) then
msg = ' plm_coeff : xx > 1 !'
LIBPAW_ERROR(msg)
end if
blm=zero
is_one=(abs(abs(xx)-one)<=tol12)
xx2=xx**2
onemx2=abs(one-xx2)
sqrx=sqrt(onemx2)
xsqrx=xx*sqrt(onemx2)
call plm_d2theta(mpsang,plm_d2t,xx)
if (is_one) then
yy=sign(one,xx)
call pl_deriv(mpsang,pl_d2,yy)
end if
do il=0,mpsang-1
llp1=dble(il*(il+1))
ilm0=il*il+il+1
do im=0,il
ilm=ilm0+im;ilm1=ilm0-im
plm =(-1)**im*ass_leg_pol(il,im,xx)
dplm_dt =(-1)**im*plm_dtheta(il,im,xx)
d2plm_dt2= plm_d2t(ilm)
blm(1,ilm)= two*xsqrx *dplm_dt+onemx2*d2plm_dt2
blm(2,ilm)= (one-two*xx2)*dplm_dt-xsqrx *d2plm_dt2
blm(3,ilm)=llp1*plm+ d2plm_dt2
blm(4,ilm)= -two*xsqrx *dplm_dt+xx2 *d2plm_dt2
if (is_one) then
if (im==1) then
blm(5,ilm)=llp1*plm+d2plm_dt2
end if
if (im==2) then
blm(5,ilm)=d2plm_dt2-three*pl_d2(il+1)
end if
else
if(im>0) then
blm(5,ilm)=plm/onemx2-dplm_dt*xx/sqrx
end if
end if
if (im>0) then
blm(1,ilm1)=blm(1,ilm)
blm(2,ilm1)=blm(2,ilm)
blm(3,ilm1)=blm(3,ilm)
blm(4,ilm1)=blm(4,ilm)
blm(5,ilm1)=blm(5,ilm)
end if
end do
end do
end subroutine plm_coeff
!!***
!----------------------------------------------------------------------
!!****f* m_paw_sphharm/ass_leg_pol
!! NAME
!! ass_leg_pol
!!
!! FUNCTION
!! Compute the associated Legendre Polynomial Plm(x),
!! using a stable recursion formula.
!! Here m and l are integers satisfying 0<=m<=l,
!! while x lies in the range -1<=x<=1
!!
!! INPUTS
!! l,m= l,m numbers
!! xarg=argument of the polynom
!!
!! OUTPUT
!!
!! SOURCE
function ass_leg_pol(l,m,xarg)
!Arguments ------------------------------------
!scalars