-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmath.f90
1712 lines (1709 loc) · 64.1 KB
/
math.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
module math
!
! Sundry mathematical definitions.
!
use accuracy
use constants
use lapack
!$ use OMP_LIB
implicit none
private
public MathFactorial, MathLogFactorial, MathLegendrePn, MathLegendrePnm
public MathJacobiPn
public MathDoubleFactorial, MathLogDoubleFactorial
public MathBinomial, MathPochhammer
public MathLogGamma, MathSimpleGamma, Math3J, MathLogSum
public MathRotationMatrix, MathYJMRotationMatrix
public MathSetUnitMatrix, MathIsUnitMatrix
public MathGetQuadrature, MathErf
public MathBoysF, MathDawsonF
public MathYLM, MathAllYLM, MathAllYLM2, MathAllYLM2scr
public MathDet3x3, MathInv3x3
public MathInterpolate
public rcsid_math
!
character(len=clen), save :: rcsid_math = "$Id: math.f90,v 1.5 2022/08/02 08:51:57 ps Exp $"
!
! Fortran does not allow choosing a specific function based on the type
! of the result alone. Therefore, we are forced to use some trickery.
!
! For functions which need such resolution, we will add an extra argument,
! with the type matching the desired return type. In order to avoid breaking
! the existing code, we also supply a function without the extra argument;
! this one will make the "default" choice of the return type.
!
! Unfortunately, this makes for a bit of clutter.
!
interface MathBinomial
module procedure MathBinomialIntegerReal
module procedure MathBinomialReal
!*qd module procedure MathBinomialIntegerQuad
end interface MathBinomial
!
interface MathFactorial
module procedure MathFactorialReal
!*qd module procedure MathFactorialQuad
end interface MathFactorial
!
interface MathDoubleFactorial
module procedure MathDoubleFactorialReal
!*qd module procedure MathDoubleFactorialQuad
end interface MathDoubleFactorial
!
interface MathLogFactorial
module procedure MathLogFactorialReal
!*qd module procedure MathLogFactorialQuad
end interface MathLogFactorial
!
interface MathDet3x3
module procedure MathDet3x3Real
module procedure MathDet3x3Complex
end interface MathDet3x3
!
interface MathInv3x3
module procedure MathInv3x3Real
module procedure MathInv3x3Complex
end interface MathInv3x3
!
interface MathIsUnitMatrix
module procedure MathIsUnitMatrixReal
!*qd module procedure MathIsUnitMatrixQuad
end interface MathIsUnitMatrix
!
interface MathSetUnitMatrix
module procedure MathSetUnitMatrixReal
!*qd module procedure MathSetUnitMatrixQuad
end interface MathSetUnitMatrix
!
interface MathBoysF
module procedure MathBoysFReal
!*qd module procedure MathBoysFQuad
end interface MathBoysF
!
integer(ik), parameter :: factorial_slack = 5 ! Extra factorials to produce while filling the cache
integer(ik), save :: factorial_max = -1 ! Largest value of the factorials cached in the table
real(rk), allocatable, save :: factorial_table(:)
integer(ik), save :: factorial_quad_max = -1 ! Largest value of the factorials cached in the table
real(xrk), allocatable, save :: factorial_quad_table(:)
integer(ik), save :: log_factorial_max = -1 ! Largest value of the factorials cached in the table
real(rk), allocatable, save :: log_factorial_table(:)
integer(ik), save :: log_factorial_quad_max = -1 ! Largest value of the factorials cached in the table
real(xrk), allocatable, save :: log_factorial_quad_table(:)
integer(ik), save :: dfactorial_max = -2 ! Largest value of the double factorials cached in the table
real(rk), allocatable, save :: dfactorial_table(:)
integer(ik), save :: dfactorial_quad_max = -2 ! Largest value of the double factorials (quad) cached in the table
real(xrk), allocatable, save :: dfactorial_quad_table(:)
integer(ik), save :: log_dfactorial_max = -2 ! Largest value of the double factorials cached in the table
real(rk), allocatable, save :: log_dfactorial_table(:)
!
contains
!
! External interfaces
!
function MathFactorialReal(n,kind) result(v)
integer(ik), intent(in) :: n
real(rk), intent(in), optional :: kind
real(rk) :: v
!
if (n<0) stop 'math%MathFactorialReal - domain error'
if (n>factorial_max) call fill_factorial_table(n+factorial_slack)
v = factorial_table(n)
end function MathFactorialReal
!
function MathFactorialQuad(n,kind) result(v)
integer(ik), intent(in) :: n
real(xrk), intent(in) :: kind
real(xrk) :: v
!
if (n<0) stop 'math%MathFactorialQuad - domain error'
if (n>factorial_quad_max) call fill_factorial_quad_table(n+factorial_slack)
v = factorial_quad_table(n)
end function MathFactorialQuad
!
function MathLogFactorialReal(n,kind) result(v)
integer(ik), intent(in) :: n
real(rk), intent(in), optional :: kind
real(rk) :: v
!
if (n<0) stop 'math%MathLogFactorial - domain error'
if (n>log_factorial_max) call fill_log_factorial_table(n+factorial_slack)
v = log_factorial_table(n)
end function MathLogFactorialReal
!
function MathLogFactorialQuad(n,kind) result(v)
integer(ik), intent(in) :: n
real(xrk), intent(in) :: kind
real(xrk) :: v
!
if (n<0) stop 'math%MathLogFactorialQuad - domain error'
if (n>log_factorial_quad_max) call fill_log_factorial_quad_table(n+factorial_slack)
v = log_factorial_quad_table(n)
end function MathLogFactorialQuad
!
function MathDoubleFactorialReal(n,kind) result(v)
integer(ik), intent(in) :: n
real(rk), intent(in), optional :: kind
real(rk) :: v
!
if (n<-1) stop 'math%MathDoubleFactorial - domain error'
if (n>dfactorial_max) call fill_dfactorial_table(n+factorial_slack)
v = dfactorial_table(n)
end function MathDoubleFactorialReal
!
function MathDoubleFactorialQuad(n,kind) result(v)
integer(ik), intent(in) :: n
real(xrk), intent(in) :: kind
real(xrk) :: v
!
if (n<-1) stop 'math%MathDoubleFactorialQuad - domain error'
if (n>dfactorial_quad_max) call fill_dfactorial_quad_table(n+factorial_slack)
v = dfactorial_quad_table(n)
end function MathDoubleFactorialQuad
!
function MathLogDoubleFactorial(n) result(v)
integer(ik), intent(in) :: n
real(rk) :: v
!
if (n<-1) stop 'math%MathLogDoubleFactorial - domain error'
if (n>log_dfactorial_max) call fill_log_dfactorial_table(n+factorial_slack)
v = log_dfactorial_table(n)
end function MathLogDoubleFactorial
!
! Pochhammer function: product of (n) integers starting at (a)
!
function MathPochhammer(a,n) result(v)
integer(ik), intent(in) :: a, n
real(rk) :: v
!
integer(ik) :: aa ! Starting integer of the equivalent positive sequence
logical :: minus
!
if (n<0) stop 'math%MathPochhammer - domain error'
!
if (n==0) then
v = 1._rk
return
end if
!
if (a<=0) then
!
! Catch sequences containing zero factors: these are always zero
!
if (a+n-1>=0) then
v = 0._rk
return
end if
aa = -(a+n-1)
minus = mod(n,2)==1
else
aa = a
minus = .false.
end if
!
v = MathLogFactorial(aa+n-1)-MathLogFactorial(aa-1)
v = exp(v)
if (minus) v = -v
end function MathPochhammer
!
! Binomial coefficients
!
function MathBinomialIntegerReal(n,m,kind) result(cnm)
integer(ik), intent(in) :: n, m
real(rk), intent(in), optional :: kind
real(rk) :: cnm
!
if (n<0 .or. m<0 .or. m>n) stop 'MathBinomialIntegerReal - domain error'
cnm = exp(MathLogFactorial(n)-MathLogFactorial(n-m)-MathLogFactorial(m))
end function MathBinomialIntegerReal
!
function MathBinomialIntegerQuad(n,m,kind) result(cnm)
integer(ik), intent(in) :: n, m
real(xrk), intent(in) :: kind
real(xrk) :: cnm
!
if (n<0 .or. m<0 .or. m>n) stop 'MathBinomialIntegerQuad - domain error'
cnm = exp(MathLogFactorial(n,kind)-MathLogFactorial(n-m,kind)-MathLogFactorial(m,kind))
end function MathBinomialIntegerQuad
!
function MathBinomialReal(n,m) result(cnm)
real(rk), intent(in) :: n, m
real(rk) :: cnm
complex(rk), parameter :: one = (1._rk,0._rk)
!
if (n<0._rk .or. m<0._rk .or. m>n) stop 'MathBinomialReal - domain error'
cnm = exp(real(MathLogGamma(n+one)-MathLogGamma(n-m+one)-MathLogGamma(m+one),kind=rk))
end function MathBinomialReal
!
! Evaluate log(exp(a)+exp(b)), where a and b may be too large to fit
! in the exponent range.
!
function MathLogSum(a,b) result(c)
real(rk), intent(in) :: a, b
real(rk) :: c
!
if (a>=b) then
c = a + log(1._rk + exp(b-a))
else
c = b + log(1._rk + exp(a-b))
end if
end function MathLogSum
!
! If all values of Pn or Pnm up to a given order are required, it is much
! more efficient to compute them all at once!
!
! Accuracy of the recursions used to calculate Ledendre functions
! deteriorates in the vicinity of the polynomial roots, especially
! for very high orders (n,m)
!
function MathLegendrePn(n,x) result(pn)
integer(ik), intent(in) :: n ! Order of the Legendre polynomial
real(rk), intent(in) :: x ! Coordinate
real(rk) :: pn ! Value of the Legenre polynomial
!
real(rk) :: tab(0:n)
!
call legendrePn_table(n,x,tab)
pn = tab(n)
end function MathLegendrePn
!
function MathLegendrePnm(n,m,x) result(pnm)
integer(ik), intent(in) :: n, m ! Order of the associated Legendre polynomial
real(rk), intent(in) :: x ! Coordinate, abs(x)<=1
real(rk) :: pnm ! Value of the associated Legenre polynomial
!
real(rk) :: tab(0:n,0:m)
!
call legendrePnm_table(n,m,x,tab)
pnm = tab(n,m)
end function MathLegendrePnm
!
! Use recurrence with respect to degree, see Abramowitz & Stegun 22.7.1
!
function MathJacobiPn(n,alp,bet,x) result(pn)
integer(ik), intent(in) :: n ! Order of the polynomial
real(rk), intent(in) :: alp, bet ! Powers of the weight function, must be > -1
real(rk), intent(in) :: x ! Coordinate, abs(x)<=1
real(rk) :: pn
!
real(rk) :: tab(0:n)
!
call jacobiPn_table(n,alp,bet,x,tab)
pn = tab(n)
end function MathJacobiPn
!
! See Abramowitz & Stegun, 22.18 for the evaluation procedure
! This routine works very well for low orders n, but becomes unstable
! for higher order. It loses 4 decimal digits of significance at n=6,
! deteriorating to complete loss of significance at n=20
!
! function MathJacobiPnUnstable(n,alp,bet,x) result(pn)
! integer(ik), intent(in) :: n ! Order of the polynomial
! real(rk), intent(in) :: alp, bet ! Powers of the weight function, must be > -1
! real(rk), intent(in) :: x ! Coordinate, abs(x)<=1
! real(rk) :: pn
! !
! integer(ik) :: m
! real(rk) :: am ! Auxiliary function; see A&S 22.18
! real(rk) :: fx ! f(x)
! !
! am = 1._rk ! am(n)
! fx = 1._rk - x ! See the second table in A&S 22.18
! recurse_am: do m=n,1,-1
! am = 1._rk - am*bm(m)*fx/cm(m) ! Calculates a(m-1)
! end do recurse_am
! pn = dn()*am
! !
! contains
! real(rk) function dn()
! dn = MathBinomial(n+alp,real(n,kind=rk))
! end function dn
! real(rk) function bm(m)
! integer(ik), intent(in) :: m
! !
! bm = (n-m+1._rk)*(alp+bet+n+m)
! end function bm
! real(rk) function cm(m)
! integer(ik), intent(in) :: m
! !
! cm = 2._rk*m*(alp+m)
! end function cm
! end function MathJacobiPnUnstable
!
! Computes Wigner 3J symbols. The code below is a direct implementation
! of L&L 3J formulae. The accuracy of this routine is reduced relative to
! that is theoretically possible, due to the use of logarithms. The routine
! I had in MNDO99 is more accurate and can handle broader range of J values.
!
function Math3J(j1,j2,j3,m1,m2,m3) result(v)
integer(ik), intent(in) :: j1, j2, j3 ! / J1 J2 J3 \ 3-j
integer(ik), intent(in) :: m1, m2, m3 ! \ M1 M2 M3 /
real(rk) :: v
!
integer(ik) :: ij0, ij1, ij2, ij3, im1a, im1b, im2a, im2b, im3a, im3b
integer(ik) :: t1, t2, t3
integer(ik) :: z, minz, maxz
real(rk) :: logscale, logterm
!
! Before we do anything, check whether this 3J symbol satisfies the
! vector addition constraints
!
ij0 = j1 + j2 + j3 + 1
ij1 = j1 + j2 - j3
ij2 = j1 - j2 + j3
ij3 = - j1 + j2 + j3
im1a = j1 - m1 ; im1b = j1 + m1
im2a = j2 - m2 ; im2b = j2 + m2
im3a = j3 - m3 ; im3b = j3 + m3
if (ij1<0 .or. ij2<0 .or. ij3<0 .or. im1a<0 .or. im1b<0 .or. im2a<0 .or. im2b<0 .or. im3a<0 .or. im3b<0 .or. m1+m2+m3/=0) then
v = 0
return
end if
!
logscale = MathLogFactorial(ij1) + MathLogFactorial(ij2) + MathLogFactorial(ij3) &
+ MathLogFactorial(im1a) + MathLogFactorial(im1b) + MathLogFactorial(im2a) &
+ MathLogFactorial(im2b) + MathLogFactorial(im3a) + MathLogFactorial(im3b) &
- MathLogFactorial(ij0)
logscale = 0.5_rk * logscale
!
t1 = j2 - j3 - m1
t2 = j1 + m2 - j3
t3 = j1 - j2 - m3
minz = max(0_ik,t1,t2)
maxz = min(ij1,im1a,im2b)
v = 0
sum_terms: do z=minz,maxz,1
logterm = logscale - MathLogFactorial(z) - MathLogFactorial(ij1-z) - MathLogFactorial(im1a-z) &
- MathLogFactorial(im2b-z) - MathLogFactorial(z-t1) - MathLogFactorial(z-t2)
if (abs(logterm)>=max_exp) then
write (out,"('Math3J: Intermediate logarithm ',g12.5,' exceeds the real(rk) dynamic range.')") logterm
write (out,"('Math3J: The 3J arguments were: ',6i10)") j1, j2, j3, m1, m2, m3
stop 'math%Math3J - exceeded dynamic range'
end if
if (mod(z+t3,2)==0) then
v = v + exp(logterm)
else
v = v - exp(logterm)
end if
end do sum_terms
!
end function Math3J
!
! Given Euler angles, construct rotation matrix for coordinate axes
! in 3D space. The Euler angles are defined as follows:
! 1. Rotate coordinate axes by alpha around the Z axis
! 2. Rotate axes by beta around the new Y axis
! 3. Rotate axes by gamma around the new Z axis
! This definition of the Euler angles matches the definition from
! section 58 of L&L III.
!
! Note that prior to June 10th, 2010 the definition of all three Euler
! angles in this routine used a wrong sign, corresponding to anti-screw
! rotation sense. Thanks, Mike.
!
! If you would like to rotate an object instead of the axes, take the
! transpose or (equivalently) replace (alpha,beta,gamma) with
! (-gamma,-beta,-alpha).
!
! Some useful special cases are:
!
! a b c
! --- --- ---
! x 0 0 Rotate coordinate system by x around the Z axis
! 0 0 x Rotate coordinate system by x around the Z axis
! 0 x 0 Rotate coordinate system by x around the Y axis
! -pi/2 x pi/2 Rotate coordinate system by x around the X axis
!
subroutine MathRotationMatrix(euler_angles,mat)
real(ark), intent(in) :: euler_angles(3) ! Euler rotation angles: alpha, beta, and gamma
real(ark), intent(out) :: mat(3,3) ! Rotation matrix
!
real(ark) :: a, b, g, rma(3,3), rmb(3,3), rmg(3,3)
!
a = euler_angles(1)
b = euler_angles(2)
g = euler_angles(3)
!
rma(1,:) = (/ cos(a), sin(a), 0._rk /)
rma(2,:) = (/ -sin(a), cos(a), 0._rk /)
rma(3,:) = (/ 0._rk, 0._rk, 1._rk /)
!
rmb(1,:) = (/ cos(b), 0._rk, -sin(b) /)
rmb(2,:) = (/ 0._rk, 1._rk, 0._rk /)
rmb(3,:) = (/ sin(b), 0._rk, cos(b) /)
!
rmg(1,:) = (/ cos(g), sin(g), 0._rk /)
rmg(2,:) = (/ -sin(g), cos(g), 0._rk /)
rmg(3,:) = (/ 0._rk, 0._rk, 1._rk /)
!
mat = matmul(rmg,matmul(rmb,rma))
end subroutine MathRotationMatrix
!
! Rotation matrix for angular momentum eigenfunctions, following L&L III Eq. 58.10
! Both integer and half-integer J values are OK.
!
! The resulting rotation matrix is accurate to 2ulp for multiplicities up to 6,
! with error increasing to 4ulp for multiplicity 20. It loses about 11 decimal places
! of accuracy for multiplicity 81, and overflows IEEE double at higher multiplicities.
!
! Note that the rotation matrix uses somewhat weird conventions: it rotates transposed
! harmonics from the primed coordinate system defined by the Euler angles back into
! the lab system:
!
! Y(L,M) = Sum Y(L,M') D(M',M)
!
! Furthermore, it looks like the expression for the Wigner matrix in the 5th Russian
! edition of L&L is actually incorrect. To get the correct expression, it is necessary
! to change the sign of the Euler beta angle. The code below is a literal implementation
! of L&L 58.10, so don't forget to flip the sign of beta when calling it!
!
subroutine MathYJMRotationMatrix(euler_angles,mult,mat)
real(ark), intent(in) :: euler_angles(3) ! Euler rotation angles: alpha, beta, and gamma
! See comments in MathRotationMatrix
integer(ik), intent(in) :: mult ! Multipliplicity of the angular-momentum state,
! mult = 2*j+1
complex(ark), intent(out) :: mat(:,:) ! Rotation matrix
!
real(ark) :: a, b, g
real(ark) :: cosb2, sinb2
complex(ark) :: expa2, expg2
integer(ik) :: j2, m2, mp2
integer(ik) :: im, imp
!
if (mult<1) then
stop 'math%MathYJMRotationMatrix - multiplicity: domain error'
end if
if (size(mat,dim=1)/=mult .or. size(mat,dim=2)/=mult) then
stop 'math%MathYJMRotationMatrix - rotation matrix dimensions do not match multiplicity'
end if
!
a = euler_angles(1)
b = euler_angles(2)
g = euler_angles(3)
!
! We need to take special care when angle beta approaches n*pi. For these angles,
!
!
sinb2 = sin(0.5_ark*b)
cosb2 = cos(0.5_ark*b)
!
expa2 = exp(cmplx(0._ark,0.5_ark*a,kind=ark))
expg2 = exp(cmplx(0._ark,0.5_ark*g,kind=ark))
!
j2 = mult - 1
mp2 = -j2
row_mp: do imp=1,mult
m2 = -j2
column_m: do im=1,mult
mat(imp,im) = sqrt(hf(j2+mp2)*hf(j2-mp2)/(hf(j2+m2)*hf(j2-m2))) &
* modJacobiPn((j2-mp2)/2,(mp2-m2)/2,(mp2+m2)/2,sinb2,cosb2) &
* expa2**m2 * expg2**mp2
m2 = m2 + 2
end do column_m
mp2 = mp2 + 2
end do row_mp
!
contains
!
! (n/2)! where n must be even
!
function hf(n) result(f)
integer(ik), intent(in) :: n ! Must be even
real(ark) :: f
!
if (mod(n,2)/=0) stop 'math%MathYJMRotationMatrix%hf - domain error'
f = MathFactorial(n/2)
end function hf
!
! Specialized derivative of Jacobi P polynomial:
!
! y^a x^b JacobiP(n,a,b,x^2-y^2)
!
! where x^2+y^2 is equal to 1. Care should be taken in evaluating this function
! when either a or b are negative: standard recursion with respect to degree
! becomes undefined in this case.
!
! As the result, we have to use series expansion around +1/-1 argument to
! evaluate this function.
!
function modJacobiPn(n,a,b,y,x) result(jp)
integer(ik), intent(in) :: n, a, b ! Parameters of the Jacobi P polynomial
real(ark), intent(in) :: y, x ! Coordinate
real(ark) :: jp
!
integer(ik) :: k
!
if (abs(y)<abs(x)) then
!
! Small y, expand JacobiP around z=+1
!
jp = 0
expand_plus1: do k=max(0,-a),n
jp = jp + MathPochhammer(a+k+1,n-k) * MathPochhammer(-n,k) * MathPochhammer(a+b+n+1,k) &
* y**(2*k+a) / MathFactorial(k)
end do expand_plus1
jp = x**b * jp / MathFactorial(n)
! write (out,"(a,3(1x,i3),3(1x,f35.25))") 'A: ',n,a,b,y,x,jp
else
!
! Small x, expand JacobiP around z=-1
!
jp = 0
expand_minus1: do k=max(0,-b),n
jp = jp + MathPochhammer(b+k+1,n-k) * MathPochhammer(-n,k) * MathPochhammer(a+b+n+1,k) &
* x**(2*k+b) / MathFactorial(k)
end do expand_minus1
jp = y**a * jp / MathFactorial(n)
if (mod(n,2)==1) jp = -jp
! write (out,"(a,3(1x,i3),3(1x,f35.25))") 'B: ',n,a,b,y,x,jp
! write (out,"()")
endif
! !
! ! The general case; do not use
! !
! jp = y**a * x**b * MathJacobiPn(n,real(a,kind=rk),real(b,kind=rk),x**2-y**2)
end function modJacobiPn
end subroutine MathYJMRotationMatrix
!
! Initialize unit matrix
!
subroutine MathSetUnitMatrixReal(m)
real(ark), intent(out) :: m(:,:) ! Matrix to initialize
!
integer(ik) :: i
!
if (size(m,dim=1)/=size(m,dim=2)) then
write (out,"('math%MathSetUnitMatrix - argument is not a square matrix?!')")
stop 'math%MathSetUnitMatrix - argument is not a square matrix?!'
end if
!
m = 0
do i=1,size(m,dim=1)
m(i,i) = 1
end do
end subroutine MathSetUnitMatrixReal
subroutine MathSetUnitMatrixQuad(m)
real(xrk), intent(out) :: m(:,:) ! Matrix to initialize
!
integer(ik) :: i
!
if (size(m,dim=1)/=size(m,dim=2)) then
write (out,"('math%MathSetUnitMatrix - argument is not a square matrix?!')")
stop 'math%MathSetUnitMatrix - argument is not a square matrix?!'
end if
!
m = 0
do i=1,size(m,dim=1)
m(i,i) = 1
end do
end subroutine MathSetUnitMatrixQuad
!
! Check matrix for unity
!
function MathIsUnitMatrixReal(m) result(isunit)
real(ark), intent(in) :: m(:,:) ! Matrix to check
logical :: isunit
!
real(kind=kind(m)) :: eps
integer(ik) :: i, j
!
eps = spacing(real(2,kind=kind(m)))
isunit = .false.
if (size(m,dim=1)/=size(m,dim=2)) return
!
scan_rows: do i=1,size(m,dim=2)
if (abs(m(i,i)-1)>eps) return
scan_upper_columns: do j=1,i-1
if (abs(m(j,i))>eps) return
end do scan_upper_columns
scan_lower_columns: do j=i+1,size(m,dim=1)
if (abs(m(j,i))>eps) return
end do scan_lower_columns
end do scan_rows
isunit = .true.
end function MathIsUnitMatrixReal
function MathIsUnitMatrixQuad(m) result(isunit)
real(xrk), intent(in) :: m(:,:) ! Matrix to check
logical :: isunit
!
real(kind=kind(m)) :: eps
integer(ik) :: i, j
!
eps = spacing(real(2,kind=kind(m)))
isunit = .false.
if (size(m,dim=1)/=size(m,dim=2)) return
!
scan_rows: do i=1,size(m,dim=2)
if (abs(m(i,i)-1)>eps) return
scan_upper_columns: do j=1,i-1
if (abs(m(j,i))>eps) return
end do scan_upper_columns
scan_lower_columns: do j=i+1,size(m,dim=1)
if (abs(m(j,i))>eps) return
end do scan_lower_columns
end do scan_rows
isunit = .true.
end function MathIsUnitMatrixQuad
!
! Calculate Gaussian quadrature rules. See: Golub and Welsch, "Calculation of Gaussian
! quadrature rules," mathematics of computation 23 (april, 1969), pp. 221-230
!
! WARNING: Only the Gauss-Hermite rule has been verified against the tables
!
subroutine MathGetQuadrature(type,order,x,w,alpha)
character(len=*), intent(in) :: type ! Desired quadrature type; see below for the allowed types
integer(ik), intent(in) :: order ! Desired quadrature order
real(rk), intent(out) :: x(:) ! Integration points
real(rk), intent(out) :: w(:) ! Integration weights
real(rk), intent(in), optional :: alpha ! Laguerre alpha
!
integer(ik) :: i, alloc
real(rk), allocatable :: a(:), b(:) ! Diagonal and sub-diagonal values for the eigensystem
real(rk), allocatable :: z(:,:) ! Space for the eigenvectors
real(rk) :: mu0 ! Zeroth momentum of the weight function
!
if (size(x)/=order .or. size(w)/=order) then
stop 'math%MathGetQuadrature - inconsistent buffer sizes'
end if
!
allocate (a(order),b(order-1),z(order,order),stat=alloc)
if (alloc/=0) then
write (out,"('Error ',i8,' allocating scratch for order-',i8,' quadrature')") alloc, order
stop 'math%MathGetQuadrature - no memory'
end if
!
! The select case below is the only part of the routine which depends on the
! quadrature required; everything else is the same for all quadratures.
!
select case (type)
case default
write (out,"('MathGetQuadrature called for unknown quadrature type ',a)") trim(type)
stop 'math%MathGetQuadrature - invalid quadrature type'
case ('Legendre')
!
! (-1,+1), 1.0
!
a = 0
legendre: do i=1,order-1
b(i) = i/sqrt(4*i**2 - 1.0_rk)
end do legendre
mu0 = 2._rk
case ('Chebyshev T') ! Chebychev of the 1st kind
!
! (-1,+1), (1-X**2)**-0.5
!
a = 0
b(1) = sqrt(0.5_rk)
b(2:) = 0.5_rk
mu0 = pi
case ('Chebyshev U') ! Chebychev of the 2nd kind
!
! (-1,+1), (1-X**2)**0.5
!
a = 0
b = 0.5_rk
mu0 = 0.5_rk*pi
case ('Hermite')
!
! (-Infinity,+Infinity), EXP(-X**2)
!
a = 0
hermite: do i=1,order-1
b(i) = sqrt(0.5_rk*i)
end do hermite
mu0 = sqrt(pi)
case ('Laguerre')
!
! (0,+Infinity), EXP(-X) * X**ALPHA
!
if (.not.present(alpha)) stop 'MathGetQuadrature(Laguerre) - missing alpha'
laguerre: do i=1,order-1
a(i) = 2*i - 1 + alpha
b(i) = sqrt(i*(i+alpha))
end do laguerre
a(order) = 2*order - 1 + alpha
mu0 = real(MathSimpleGamma(cmplx(1._rk + alpha,0._rk,kind=rk)),kind=rk)
end select
!
call lapack_stev(a,b,z)
!
! Copy integration abscissas and weights, and we are done
!
x = a
w = mu0 * z(1,:)**2
!
deallocate (a,b,z)
end subroutine MathGetQuadrature
!
! Evaluate a given spherical harmonic. If you need to evaluate the same spherical
! harmonic at multiple points, it is preferable to use FLharmonics() in fields.f90
!
function MathYLM(l,m,dir) result(ylm)
integer(ik), intent(in) :: l, m
real(rk), intent(in) :: dir(3) ! (Unnormalized) direction vector
complex(rk) :: ylm
!
complex(rk) :: harm_fact
real(rk) :: r, ct, xymod
complex(rk) :: xy
!
harm_fact = (-1)**((m+abs(m))/2)
harm_fact = harm_fact * (0._rk,1._rk)**l
harm_fact = harm_fact * sqrt((2*l+1)/fourpi)
harm_fact = harm_fact * sqrt(MathFactorial(l-abs(m))/MathFactorial(l+abs(m)))
!
r = sqrt(sum(dir**2))
!
ct = dir(3)/r
if (m>0) then
xy = cmplx(dir(1), dir(2),kind=rk)
else
xy = cmplx(dir(1),-dir(2),kind=rk)
end if
xymod = abs(xy)
if (xymod>0._rk) then
xy = xy / xymod
else
xy = 1._rk
end if
ylm = harm_fact * MathLegendrePnm(l,abs(m),ct) * xy**abs(m)
end function MathYLM
!
! Calculate all spherical harmonics up to angular momentum l_max
! WARNING: MathAllYLM (and MathYLM) starts to lose accuracy beyond
! WARNING: L=30. For L=50, we have only about 5 significant digits
! WARNING: in double precision; at L=80, we have about 2 digits.
! For more accurate results, consider replacing with MathAllYLM2
! from SCID-TDSE.
!
subroutine MathAllYLM(l_max,dir,ylm)
integer(ik), intent(in) :: l_max ! Desired maximum angular momentum
real(rk), intent(in) :: dir(3) ! (Unnormalized) direction vector
complex(rk), intent(out) :: ylm(-l_max:l_max,0:l_max)
!
integer(ik) :: l, m
complex(rk) :: harm_fact
real(rk) :: r, ct, xymod
complex(rk) :: xy_plus, xy_minus
real(rk) :: legendre_tab(0:l_max,0:l_max)
complex(rk) :: xy_powm(-l_max:l_max)
!
if (l_max>=50) stop 'MathAllYLM - accuracy loss'
r = sqrt(sum(dir**2))
ct = dir(3)/r
call legendrePnm_table(l_max,l_max,ct,legendre_tab)
!
xy_plus = cmplx(dir(1), dir(2),kind=rk)
xy_minus = cmplx(dir(1),-dir(2),kind=rk)
xymod = abs(xy_plus)
if (xymod>0._rk) then
xy_plus = xy_plus / xymod
xy_minus = xy_minus / xymod
else
xy_plus = 1._rk
xy_minus = 1._rk
end if
!
xy_pow1: do m=-l_max,-1
xy_powm(m) = xy_minus ** abs(m)
end do xy_pow1
xy_powm(0) = 1._rk
xy_pow2: do m=1,l_max
xy_powm(m) = xy_plus ** abs(m)
end do xy_pow2
!
tab_l: do l=0,l_max
tab_m: do m=-l,l
harm_fact = (-1)**((m+abs(m))/2) * (0._rk,1._rk)**l * sqrt((2*l+1)/fourpi) &
* sqrt(MathFactorial(l-abs(m))/MathFactorial(l+abs(m)))
ylm(m,l) = harm_fact * legendre_tab(l,abs(m)) * xy_powm(m)
end do tab_m
end do tab_l
end subroutine MathAllYLM
!
! Our old spherical harmonics code has a problem with large values of L, where
! we need to multiply through many large and small factors to give an answer on
! the order of 1. The code below integrates evaluation of the prefactor into
! the recurrences for the associated Legendre polynomials, so that not large
! factors arise. By switching to a different recursion, we can also support
! the case of a sub-set of possible M values being required.
!
! Allocation of the two arrays below may come _very_ expensive on Xeon Phi; consider using
! the MathAllYLM2scr entry point directly if you are evaluating a large number of spherical
! harmonics in a loop.
!
subroutine MathAllYLM2(l_max,m_min,m_max,dir,ylm,phase)
integer(ik), intent(in) :: l_max ! Desired maximum angular momentum
integer(ik), intent(in) :: m_min ! Desired minimum Z projection of the angular momentum
integer(ik), intent(in) :: m_max ! Desired maximum Z projection of the angular momentum
real(rk), intent(in) :: dir(3) ! (Unnormalized) direction vector
complex(rk), intent(out) :: ylm(m_min:m_max,0:l_max)
character(len=*), intent(in), optional :: phase ! Phase convention to use. Recognized values are:
! 'L&L' = Landau and Lifshitz v. 3 phase convention (the default)
! 'Arfken' = Arfken and Mathematica phase convention
!
real(rk) :: sqrtn(0:2*l_max+1) ! Square roots of small integers
real(rk) :: rsqrn(1:2*l_max+1) ! Inverse square roots of integers
!
call MathAllYLM2scr(l_max,m_min,m_max,dir,ylm,sqrtn,rsqrn,phase)
end subroutine MathAllYLM2
!
subroutine MathAllYLM2scr(l_max,m_min,m_max,dir,ylm,sqrtn,rsqrn,phase)
integer(ik), intent(in) :: l_max ! Desired maximum angular momentum
integer(ik), intent(in) :: m_min ! Desired minimum Z projection of the angular momentum
integer(ik), intent(in) :: m_max ! Desired maximum Z projection of the angular momentum
real(rk), intent(in) :: dir(3) ! (Unnormalized) direction vector
complex(rk), intent(out) :: ylm(m_min:m_max,0:l_max)
real(rk) :: sqrtn(0:2*l_max+1) ! Square roots of small integers
real(rk) :: rsqrn(1:2*l_max+1) ! Inverse square roots of integers
character(len=*), intent(in), optional :: phase ! Phase convention to use. Recognized values are:
! 'L&L' = Landau and Lifshitz v. 3 phase convention (the default)
! 'Arfken' = Arfken and Mathematica phase convention
!
! I**N
complex(rk), parameter :: ipow(0:3) = (/ (1._rk,0._rk), (0._rk,1._rk), (-1._rk,0._rk), (0._rk,-1._rk) /)
real(rk) :: y00 ! Y00 does not depend on the spherical angles
integer(ik) :: k, l, m, m_low, m_high
! Allocation of the two arrays below may come _very_ expensive on Xeon Phi; consider passing them
! as argumens if forced inlining does not help
real(rk) :: sinth, costh ! Sine and cosine of the spherical theta angle
real(rk) :: r, xymod
complex(rk) :: xpy ! exp(+(0._rk,1._rk)*phi), where phi is the spherical phi angle
complex(rk) :: xmy ! exp(-(0._rk,1._rk)*phi)
complex(rk) :: vlm ! Current YLM in simultaneous L,M recurrence
!
! Sanity testing
!
if (l_max<0 .or. abs(m_min)>l_max .or. abs(m_max)>l_max .or. m_min>m_max) then
write (out,"('MathAllYLM2: l_max= ',i0,' m_min= ',i0,' m_max= ',i0)") l_max, m_min, m_max
! call flush_wrapper(out)
stop 'math%MathAllYLM2 - called with invalid L,M arguments'
end if
!
! Spherical parameters
!
y00 = 1._rk/sqrt(4._rk*pi)
r = sqrt(sum(dir**2))
if (r<=0) then
stop 'math%MathAllYLM2 - direction cannot be zero'
end if
xpy = cmplx(dir(1), dir(2),kind=rk)
xmy = cmplx(dir(1),-dir(2),kind=rk)
xymod = abs(xpy)
if (xymod>0) then
xpy = xpy / xymod
xmy = xmy / xymod
else
xpy = 1
xmy = 1
end if
costh = min(1._rk,max(-1._rk,dir(3)/r)) ! Force sine and cosine into the allowed range,
sinth = min(1._rk,max( 0._rk,xymod /r)) ! to guard against round-off errors
!
! Precompute some of the factors we need
!
sqrtn(0) = 0
fill_square_roots: do k=1,2*l_max+1
sqrtn(k) = sqrt(real(k,kind=rk))
rsqrn(k) = 1._rk / sqrtn(k)
end do fill_square_roots
!
if (m_min<=0 .and. m_max>=0) ylm(0,0) = y00
!
! Simultaneously increase L and M, starting at zero
!
vlm = y00
upward_diagonal_lm: do l=1,l_max
vlm = vlm * (0._rk,1._rk) * sqrtn(2*l+1)*rsqrn(2*l)*xpy*sinth
if (m_min<=l .and. m_max>=l) ylm(l,l) = vlm
end do upward_diagonal_lm
!
! Simultaneously increase L and decrease M, starting at zero
!
vlm = y00
downward_diagonal_lm: do l=1,l_max
vlm = vlm * (0,-1) * sqrtn(2*l+1)*rsqrn(2*l)*xmy*sinth
if (m_min<=-l .and. m_max>=-l) ylm(-l,l) = vlm
end do downward_diagonal_lm
!
! For each M increase L, starting at L=M
!
upward_l: do l=1,l_max
m_low = max(m_min,-(l-1)) ! We also must satisfy |m|<=l-1 for this recursion
m_high = min(m_max, (l-1))
upward_l_m_loop: do m=m_low,m_high
vlm = ylm(m,l-1) * (0._rk,1._rk)*costh*sqrtn(2*l-1)
if (l>=abs(m)+2) then
vlm = vlm + ylm(m,l-2) * rsqrn(2*l-3)*sqrtn(l-1-m)*sqrtn(l-1+m)
end if
ylm(m,l) = vlm * sqrtn(2*l+1)*rsqrn(l-m)*rsqrn(l+m)
end do upward_l_m_loop
end do upward_l
!
if (present(phase)) then
select case (phase)
case default; stop 'math%MathAllYLM2 - Invalid phase convention'
case ('L&L') ! Do nothing
case ('Arfken')
adjust_l: do l=0,l_max
adjust_m: do m=m_min,m_max
! For some reason, Intel Fortran has real trouble generating decent code for the
! commented line below. Let's do a hack here ...
! ylm(m,l) = ylm(m,l) * (0._rk,1._rk)**(-2*m-l)
ylm(m,l) = ylm(m,l) * ipow(modulo(-2*m-l,4))
end do adjust_m
end do adjust_l
end select
end if
end subroutine MathAllYLM2scr
!
! Auxiliary functions
!
subroutine legendrePn_table(nmax,x,pn)
integer(ik), intent(in) :: nmax ! Maximum order of the Legendre polynomials desired
real(rk), intent(in) :: x ! Coordinate at which P values are needed
real(rk), intent(out) :: pn(:) ! Values of LegendreP from n=0 to n=nmax
!
integer(ik) :: n
real(rk) :: invn
!
if (nmax<0) stop 'math%legendreP_table - negative-order polynomial requested'
pn(1) = 1._rk
if (nmax<1) return
pn(2) = x
!
n_recursion: do n=2,nmax
invn = 1._rk / n
pn(1+n) = (2._rk-invn)*x*pn(1+(n-1)) - (1._rk-invn)*pn(1+(n-2))
end do n_recursion
end subroutine legendrePn_table
!
subroutine legendrePnm_table(nmax,mmax,x,pnm)
integer(ik), intent(in) :: nmax ! Maximum order of the Legendre polynomials desired
integer(ik), intent(in) :: mmax ! Maximum order of the Legendre polynomials desired
real(rk), intent(in) :: x ! Coordinate at which P values are needed, abs(x) must be <=1
real(rk), intent(out) :: pnm(:,:) ! Values of LegendreP from n,m=0 to n,m=nmax,mmax
! n is the first subscript; m is the second subscript
!
integer(ik) :: n, m
real(rk) :: sqfac ! sqrt(1-x**2)
!
sqfac = 1._rk - x**2
if (sqfac<0) stop 'math%legendrePnm_table - domain error'
sqfac = sqrt(sqfac)
!
call legendrePn_table(nmax,x,pnm(:,1))
if (mmax<1) return
!
! Special case for m=1: recursion is truncated for n=1, m=1
!