Commit 42591138ec882c4674db5641e0f153c33d0362bc
1 parent
9158e74d6b
Exists in
master
and in
3 other branches
git-svn-id: https://lxsd.femto-st.fr/svn/fvn@26 b657c933-2333-4658-acf2-d3c7c2708721
Showing 1 changed file with 178 additions and 1 deletions Side-by-side Diff
fvnlib.f90
| ... | ... | @@ -6,7 +6,7 @@ |
| 6 | 6 | ! it uses lapack for linear algebra |
| 7 | 7 | ! it uses modified quadpack for integration |
| 8 | 8 | ! |
| 9 | -! William Daniau 2007 | |
| 9 | +! William Daniau 2007->today | |
| 10 | 10 | ! william.daniau@femto-st.fr |
| 11 | 11 | ! |
| 12 | 12 | ! Routines naming scheme : |
| ... | ... | @@ -22,6 +22,9 @@ |
| 22 | 22 | ! if you find it usefull it would be kind to give credits ;-) |
| 23 | 23 | ! |
| 24 | 24 | ! svn version |
| 25 | +! January 2008 : added quadratic interpolation, gamma/factorial function, | |
| 26 | +! a function which return identity matrix, | |
| 27 | +! evaluation of nterm chebyshev series | |
| 25 | 28 | ! September 2007 : added sparse system solving by interfacing umfpack |
| 26 | 29 | ! June 2007 : added some complex trigonometric functions |
| 27 | 30 | ! |
| ... | ... | @@ -51,6 +54,65 @@ |
| 51 | 54 | |
| 52 | 55 | contains |
| 53 | 56 | |
| 57 | +!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! | |
| 58 | +! | |
| 59 | +! Identity Matrix | |
| 60 | +! | |
| 61 | +!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! | |
| 62 | +function fvn_d_ident(n) | |
| 63 | + implicit none | |
| 64 | + integer(kind=4) :: n | |
| 65 | + real(kind=8), dimension(n,n) :: fvn_d_ident | |
| 66 | + | |
| 67 | + real(kind=8),dimension(n*n) :: vect | |
| 68 | + integer(kind=4) :: i | |
| 69 | + | |
| 70 | + vect=0._8 | |
| 71 | + vect(1:n*n:n+1) = 1._8 | |
| 72 | + fvn_d_ident=reshape(vect, shape = (/ n,n /)) | |
| 73 | +end function | |
| 74 | + | |
| 75 | +function fvn_s_ident(n) | |
| 76 | + implicit none | |
| 77 | + integer(kind=4) :: n | |
| 78 | + real(kind=4), dimension(n,n) :: fvn_s_ident | |
| 79 | + | |
| 80 | + real(kind=4),dimension(n*n) :: vect | |
| 81 | + integer(kind=4) :: i | |
| 82 | + | |
| 83 | + vect=0._4 | |
| 84 | + vect(1:n*n:n+1) = 1._4 | |
| 85 | + fvn_s_ident=reshape(vect, shape = (/ n,n /)) | |
| 86 | +end function | |
| 87 | + | |
| 88 | +function fvn_c_ident(n) | |
| 89 | + implicit none | |
| 90 | + integer(kind=4) :: n | |
| 91 | + complex(kind=4), dimension(n,n) :: fvn_c_ident | |
| 92 | + | |
| 93 | + complex(kind=4),dimension(n*n) :: vect | |
| 94 | + integer(kind=4) :: i | |
| 95 | + | |
| 96 | + vect=(0._4,0._4) | |
| 97 | + vect(1:n*n:n+1) = (1._4,0._4) | |
| 98 | + fvn_c_ident=reshape(vect, shape = (/ n,n /)) | |
| 99 | +end function | |
| 100 | + | |
| 101 | +function fvn_z_ident(n) | |
| 102 | + implicit none | |
| 103 | + integer(kind=4) :: n | |
| 104 | + complex(kind=8), dimension(n,n) :: fvn_z_ident | |
| 105 | + | |
| 106 | + complex(kind=8),dimension(n*n) :: vect | |
| 107 | + integer(kind=4) :: i | |
| 108 | + | |
| 109 | + vect=(0._8,0._8) | |
| 110 | + vect(1:n*n:n+1) = (1._8,0._8) | |
| 111 | + fvn_z_ident=reshape(vect, shape = (/ n,n /)) | |
| 112 | +end function | |
| 113 | + | |
| 114 | + | |
| 115 | + | |
| 54 | 116 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
| 55 | 117 | ! |
| 56 | 118 | ! Matrix inversion subroutines |
| 57 | 119 | |
| 58 | 120 | |
| ... | ... | @@ -3124,8 +3186,123 @@ |
| 3124 | 3186 | deallocate(wTi,wTj) |
| 3125 | 3187 | end subroutine |
| 3126 | 3188 | |
| 3189 | +!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! | |
| 3190 | +! | |
| 3191 | +! Special Functions | |
| 3192 | +! | |
| 3193 | +!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! | |
| 3194 | +function fvn_d_lngamma(x) | |
| 3195 | + ! This function returns ln(gamma(x)) | |
| 3196 | + ! adapted from Numerical Recipes | |
| 3197 | + implicit none | |
| 3198 | + real(kind=8) :: x | |
| 3199 | + real(kind=8) :: fvn_d_lngamma | |
| 3127 | 3200 | |
| 3201 | + real(kind=8) :: ser,stp,tmp,y,cof(6) | |
| 3202 | + integer(kind=4) :: i | |
| 3128 | 3203 | |
| 3204 | + cof = (/ 76.18009172947146d0,-86.50532032941677d0, & | |
| 3205 | + 24.01409824083091d0,-1.231739572450155d0, & | |
| 3206 | + .1208650973866179d-2,-.5395239384953d-5 /) | |
| 3207 | + stp = 2.5066282746310005d0 | |
| 3208 | + | |
| 3209 | + tmp=x+5.5d0 | |
| 3210 | + tmp=(x+0.5d0)*log(tmp)-tmp | |
| 3211 | + | |
| 3212 | + ser=1.000000000190015d0 | |
| 3213 | + | |
| 3214 | + y=x | |
| 3215 | + | |
| 3216 | + do i=1,6 | |
| 3217 | + y=y+1.d0 | |
| 3218 | + ser=ser+cof(i)/y | |
| 3219 | + end do | |
| 3220 | + fvn_d_lngamma=tmp+log(stp*ser/x) | |
| 3221 | +end function | |
| 3222 | + | |
| 3223 | +function fvn_d_factorial(n) | |
| 3224 | + ! This function returns factorial(n) as a real(8) | |
| 3225 | + ! adapted from Numerical Recipes | |
| 3226 | + ! real value is calculated for integers lower than 32 | |
| 3227 | + implicit none | |
| 3228 | + integer(kind=4) :: n | |
| 3229 | + real(kind=8) :: fvn_d_factorial | |
| 3230 | + | |
| 3231 | + integer(kind=4) :: j | |
| 3232 | + | |
| 3233 | + fvn_d_factorial=1. | |
| 3234 | + | |
| 3235 | + if (n < 0) then | |
| 3236 | + write(*,*) "Factorial of a negative integer" | |
| 3237 | + stop | |
| 3238 | + end if | |
| 3239 | + | |
| 3240 | + if (n == 0) then | |
| 3241 | + return | |
| 3242 | + end if | |
| 3243 | + | |
| 3244 | + if (n <= 32) then | |
| 3245 | + do j=1,n | |
| 3246 | + fvn_d_factorial=fvn_d_factorial*j | |
| 3247 | + end do | |
| 3248 | + return | |
| 3249 | + else | |
| 3250 | + fvn_d_factorial=exp(fvn_d_lngamma(dble(n)+1.)) | |
| 3251 | + return | |
| 3252 | + end if | |
| 3253 | +end function | |
| 3254 | + | |
| 3255 | +function fvn_d_csevl(x,a,n) | |
| 3256 | + implicit none | |
| 3257 | + ! This function evaluate the n-term chebyshev series a at x | |
| 3258 | + ! directly adapted from http://www.netlib.org/fn | |
| 3259 | + real(kind=8) :: x | |
| 3260 | + real(kind=8), dimension(n) :: a | |
| 3261 | + integer(kind=4) :: n | |
| 3262 | + real(kind=8) :: fvn_d_csevl | |
| 3263 | + | |
| 3264 | + real(kind=8) :: twox, b0, b1, b2 | |
| 3265 | + integer(kind=4) :: i,ni | |
| 3266 | + | |
| 3267 | + twox = 2.0d0*x | |
| 3268 | + b1 = 0.d0 | |
| 3269 | + b0 = 0.d0 | |
| 3270 | + do i=1,n | |
| 3271 | + b2 = b1 | |
| 3272 | + b1 = b0 | |
| 3273 | + ni = n - i + 1 | |
| 3274 | + b0 = twox*b1 - b2 + a(ni) | |
| 3275 | + end do | |
| 3276 | + | |
| 3277 | + fvn_d_csevl = 0.5d0 * (b0-b2) | |
| 3278 | + | |
| 3279 | +end function | |
| 3280 | + | |
| 3281 | +function fvn_s_csevl(x,a,n) | |
| 3282 | + implicit none | |
| 3283 | + ! This function evaluate the n-term chebyshev series a at x | |
| 3284 | + ! directly adapted from http://www.netlib.org/fn | |
| 3285 | + real(kind=4) :: x | |
| 3286 | + real(kind=4), dimension(n) :: a | |
| 3287 | + integer(kind=4) :: n | |
| 3288 | + real(kind=4) :: fvn_s_csevl | |
| 3289 | + | |
| 3290 | + real(kind=4) :: twox, b0, b1, b2 | |
| 3291 | + integer(kind=4) :: i,ni | |
| 3292 | + | |
| 3293 | + twox = 2.0d0*x | |
| 3294 | + b1 = 0.d0 | |
| 3295 | + b0 = 0.d0 | |
| 3296 | + do i=1,n | |
| 3297 | + b2 = b1 | |
| 3298 | + b1 = b0 | |
| 3299 | + ni = n - i + 1 | |
| 3300 | + b0 = twox*b1 - b2 + a(ni) | |
| 3301 | + end do | |
| 3302 | + | |
| 3303 | + fvn_s_csevl = 0.5d0 * (b0-b2) | |
| 3304 | + | |
| 3305 | +end function | |
| 3129 | 3306 | |
| 3130 | 3307 | |
| 3131 | 3308 |