ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/SHAPES/calc_shapes.f90
Revision: 1317
Committed: Tue Jun 29 03:05:57 2004 UTC (20 years, 2 months ago) by gezelter
File size: 3886 byte(s)
Log Message:
fixes

File Contents

# User Rev Content
1 gezelter 1314 module shapes
2     implicit none
3     PRIVATE
4    
5     INTEGER, PARAMETER:: CHEBYSHEV_TN = 1
6     INTEGER, PARAMETER:: CHEBYSHEV_UN = 2
7     INTEGER, PARAMETER:: LAGUERRE = 3
8     INTEGER, PARAMETER:: HERMITE = 4
9    
10 gezelter 1317 contains
11 gezelter 1314
12     SUBROUTINE Get_Associated_Legendre(x, l, m, lmax, plm, dlm)
13    
14     ! Purpose: Compute the associated Legendre functions
15     ! Plm(x) and their derivatives Plm'(x)
16     ! Input : x --- Argument of Plm(x)
17     ! l --- Order of Plm(x), l = 0,1,2,...,n
18     ! m --- Degree of Plm(x), m = 0,1,2,...,N
19     ! lmax --- Physical dimension of PLM and DLM
20     ! Output: PLM(l,m) --- Plm(x)
21     ! DLM(l,m) --- Plm'(x)
22    
23     real (kind=8), intent(in) :: x
24     integer, intent(in) :: lmax, l, m
25     real (kind=8), dimension(0:MM,0:N), intent(inout) :: PLM(0:lmax, 0:m)
26     real (kind=8), dimension(0:MM,0:N), intent(inout) :: DLM(0:lmax, 0:m)
27     integer :: i, j
28     real (kind=8) :: xq, xs
29     integer :: ls
30    
31     ! zero out both arrays:
32     DO I = 0, m
33     DO J = 0, l
34     PLM(J,I) = 0.0D0
35     DLM(J,I) = 0.0D0
36     end DO
37     end DO
38    
39     ! start with 0,0:
40     PLM(0,0) = 1.0D0
41    
42     ! x = +/- 1 functions are easy:
43     IF (abs(X).EQ.1.0D0) THEN
44     DO I = 1, m
45     PLM(0, I) = X**I
46     DLM(0, I) = 0.5D0*I*(I+1.0D0)*X**(I+1)
47     end DO
48     DO J = 1, m
49     DO I = 1, l
50     IF (I.EQ.1) THEN
51     DLM(I, J) = 1.0D+300
52     ELSE IF (I.EQ.2) THEN
53     DLM(I, J) = -0.25D0*(J+2)*(J+1)*J*(J-1)*X**(J+1)
54     ENDIF
55     end DO
56     end DO
57     RETURN
58     ENDIF
59    
60     LS = 1
61     IF (abs(X).GT.1.0D0) LS = -1
62     XQ = sqrt(LS*(1.0D0-X*X))
63     XS = LS*(1.0D0-X*X)
64    
65     DO I = 1, l
66     PLM(I, I) = -LS*(2.0D0*I-1.0D0)*XQ*PLM(I-1, I-1)
67     enddo
68    
69     DO I = 0, l
70     PLM(I, I+1)=(2.0D0*I+1.0D0)*X*PLM(I, I)
71     enddo
72    
73     DO I = 0, l
74     DO J = I+2, m
75     PLM(I, J)=((2.0D0*J-1.0D0)*X*PLM(I,J-1) - (I+J-1.0D0)*PLM(I,J-2))/(J-I)
76     end DO
77     end DO
78    
79     DLM(0, 0)=0.0D0
80    
81     DO J = 1, m
82     DLM(0, J)=LS*J*(PLM(0,J-1)-X*PLM(0,J))/XS
83     end DO
84    
85     DO I = 1, l
86     DO J = I, m
87     DLM(I,J) = LS*I*X*PLM(I, J)/XS + (J+I)*(J-I+1.0D0)/XQ*PLM(I-1, J)
88     end DO
89     end DO
90    
91     RETURN
92     END SUBROUTINE Get_Associated_Legendre
93    
94    
95     subroutine Get_Orthogonal_Polynomial(x, m, function_type, pl, dpl)
96    
97     ! Purpose: Compute orthogonal polynomials: Tn(x) or Un(x),
98     ! or Ln(x) or Hn(x), and their derivatives
99     ! Input : function_type --- Function code
100     ! =1 for Chebyshev polynomial Tn(x)
101     ! =2 for Chebyshev polynomial Un(x)
102     ! =3 for Laguerre polynomial Ln(x)
103     ! =4 for Hermite polynomial Hn(x)
104     ! n --- Order of orthogonal polynomials
105     ! x --- Argument of orthogonal polynomials
106     ! Output: PL(n) --- Tn(x) or Un(x) or Ln(x) or Hn(x)
107     ! DPL(n)--- Tn'(x) or Un'(x) or Ln'(x) or Hn'(x)
108    
109     real(kind=8), intent(in) :: x
110     integer, intent(in):: m
111     integer, intent(in):: function_type
112 gezelter 1317 real(kind=8), dimension(0:m), intent(inout) :: pl, dpl
113 gezelter 1314
114 gezelter 1317 real(kind=8) :: a, b, c, y0, y1, dy0, dy1, yn, dyn
115     integer :: k
116 gezelter 1314
117     A = 2.0D0
118     B = 0.0D0
119     C = 1.0D0
120     Y0 = 1.0D0
121     Y1 = 2.0D0*X
122     DY0 = 0.0D0
123     DY1 = 2.0D0
124     PL(0) = 1.0D0
125     PL(1) = 2.0D0*X
126     DPL(0) = 0.0D0
127     DPL(1) = 2.0D0
128     IF (function_type.EQ.CHEBYSHEV_TN) THEN
129     Y1 = X
130     DY1 = 1.0D0
131     PL(1) = X
132     DPL(1) = 1.0D0
133     ELSE IF (function_type.EQ.LAGUERRE) THEN
134     Y1 = 1.0D0-X
135     DY1 = -1.0D0
136     PL(1) = 1.0D0-X
137     DPL(1) = -1.0D0
138     ENDIF
139 gezelter 1317 DO K = 2, m
140 gezelter 1314 IF (function_type.EQ.LAGUERRE) THEN
141     A = -1.0D0/K
142     B = 2.0D0+A
143     C = 1.0D0+A
144     ELSE IF (function_type.EQ.HERMITE) THEN
145     C = 2.0D0*(K-1.0D0)
146     ENDIF
147     YN = (A*X+B)*Y1-C*Y0
148     DYN = A*Y1+(A*X+B)*DY1-C*DY0
149     PL(K) = YN
150     DPL(K) = DYN
151     Y0 = Y1
152     Y1 = YN
153     DY0 = DY1
154     DY1 = DYN
155     end DO
156     RETURN
157    
158     end subroutine Get_Orthogonal_Polynomial
159 gezelter 1317
160     end module shapes