Scrambler  1
splines.f90
Go to the documentation of this file.
00001 !#########################################################################
00002 !               
00003 !    Copyright (C) 2003-2012 Department of Physics and Astronomy,
00004 !                            University of Rochester,
00005 !                            Rochester, NY
00006 !
00007 !    splines.f90 is part of AstroBEAR.
00008 !
00009 !    AstroBEAR is free software: you can redistribute it and/or modify    
00010 !    it under the terms of the GNU General Public License as published by 
00011 !    the Free Software Foundation, either version 3 of the License, or    
00012 !    (at your option) any later version.
00013 !
00014 !    AstroBEAR is distributed in the hope that it will be useful, 
00015 !    but WITHOUT ANY WARRANTY; without even the implied warranty of
00016 !    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00017 !    GNU General Public License for more details.
00018 !
00019 !    You should have received a copy of the GNU General Public License
00020 !    along with AstroBEAR.  If not, see <http://www.gnu.org/licenses/>.
00021 !
00022 !#########################################################################
00023 MODULE Splines
00024    USE GlobalDeclarations
00025    USE PhysicsDeclarations
00026 
00027    IMPLICIT NONE
00028 
00029    TYPE SplineDef
00030       REAL(KIND=qPREC), POINTER, DIMENSION(:) :: x=>NULL()
00031       REAL(KIND=qPREC), POINTER, DIMENSION(:) :: y=>NULL()
00032       REAL(KIND=qPREC), POINTER, DIMENSION(:) :: b=>NULL()
00033       REAL(KIND=qPREC), POINTER, DIMENSION(:) :: c=>NULL()
00034       REAL(KIND=qPREC), POINTER, DIMENSION(:) :: d=>NULL()
00035    END type SplineDef
00036 
00037    TYPE pSplineDef
00038       TYPE(SplineDef), POINTER :: p
00039    END type pSplineDef
00040 
00041    ! Example of Usage:
00042    ! Call CreateSpline(Spline, 10)  !10 points
00043    ! Spline%x(1:10)=REAL((/(i,i=1,10)/))   !Initialize x values
00044    ! Spline%y=Spline%x**3+2d0*Spline%x**2   !Initialize y values
00045    ! CALL SolveSpline(Spline)
00046    ! write(*,*) 'value at x=1.5 is ' , GetSplineValue(Spline, 1.5d0)
00047    ! CALL DestroySpline(Spline)
00048 
00049 CONTAINS
00050 
00051 
00052    SUBROUTINE CreateSpline(Spline, size)
00053       TYPE(SplineDef), POINTER :: Spline
00054       INTEGER :: size
00055       ALLOCATE(Spline)
00056       ALLOCATE(Spline%x(size))
00057       ALLOCATE(Spline%y(size))
00058    END SUBROUTINE CreateSpline
00059 
00060    SUBROUTINE DestroySpline(Spline)
00061       TYPE(SplineDef), POINTER :: Spline
00062       IF (ASSOCIATED(Spline%x)) DEALLOCATE(Spline%x)
00063       IF (ASSOCIATED(Spline%y)) DEALLOCATE(Spline%y)
00064       IF (ASSOCIATED(Spline%b)) DEALLOCATE(Spline%b)
00065       IF (ASSOCIATED(Spline%c)) DEALLOCATE(Spline%c)
00066       IF (ASSOCIATED(Spline%d)) DEALLOCATE(Spline%d)
00067       DEALLOCATE(Spline)
00068       NULLIFY(Spline)
00069    END SUBROUTINE DestroySpline
00070 
00071 
00072    !  Calculate the coefficients b(i), c(i), and d(i), i=1,2,...,n for cubic spline interpolation
00073    !  s(x) = y(i) + b(i)*(x-x(i)) + c(i)*(x-x(i))**2 + d(i)*(x-x(i))**3
00074    SUBROUTINE SolveSpline (Spline)
00075       TYPE(SplineDef) :: Spline
00076       INTEGER i, j, n
00077       REAL(KIND=qPREC) :: h
00078       REAL(KIND=qPREC), DIMENSION(:), POINTER :: b,c,d,x,y
00079       n=size(Spline%x)
00080       ALLOCATE(Spline%b(n), Spline%c(n), Spline%d(n))
00081       b=>Spline%b
00082       c=>Spline%c
00083       d=>Spline%d
00084       x=>Spline%x
00085       y=>Spline%y
00086 
00087 
00088       IF ( n < 2 ) RETURN !Need at least two points
00089       IF ( n < 3 ) THEN !Do linear interpolation
00090          b(1) = (y(2)-y(1))/(x(2)-x(1))   
00091          c(1) = 0d0
00092          d(1) = 0d0
00093          b(2) = b(1)
00094          c(2) = 0d0
00095          d(2) = 0d0
00096          RETURN
00097       END IF
00098 
00099       d(1) = x(2) - x(1)
00100       c(2) = (y(2) - y(1))/d(1)
00101       DO i = 2, n-1
00102          d(i) = x(i+1) - x(i)
00103          b(i) = 2d0*(d(i-1) + d(i))
00104          c(i+1) = (y(i+1) - y(i))/d(i)
00105          c(i) = c(i+1) - c(i)
00106       END DO
00107 
00108       b(1) = -d(1)
00109       b(n) = -d(n-1)
00110       c(1) = 0d0
00111       c(n) = 0d0
00112       IF(n /= 3) THEN
00113          c(1) = c(3)/(x(4)-x(2)) - c(2)/(x(3)-x(1))
00114          c(n) = c(n-1)/(x(n)-x(n-2)) - c(n-2)/(x(n-1)-x(n-3))
00115          c(1) = c(1)*d(1)**2/(x(4)-x(1))
00116          c(n) = -c(n)*d(n-1)**2/(x(n)-x(n-3))
00117       END IF
00118 
00119       DO i = 2, n
00120          h = d(i-1)/b(i-1)
00121          b(i) = b(i) - h*d(i-1)
00122          c(i) = c(i) - h*c(i-1)
00123       END DO
00124 
00125       c(n) = c(n)/b(n)
00126       DO j = 1, n-1
00127          i = n-j
00128          c(i) = (c(i) - d(i)*c(i+1))/b(i)
00129       END DO
00130 
00131       b(n) = (y(n) - y(n-1))/d(n-1) + d(n-1)*(c(n-1) + 2d0*c(n))
00132       DO i = 1, n-1
00133          b(i) = (y(i+1) - y(i))/d(i) - d(i)*(c(i+1) + 2d0*c(i))
00134          d(i) = (c(i+1) - c(i))/d(i)
00135          c(i) = 3d0*c(i)
00136       END DO
00137       c(n) = 3d0*c(n)
00138       d(n) = d(n-1)
00139    END SUBROUTINE SolveSpline
00140 
00141 
00143    FUNCTION GetSplineValue(Spline, u)
00144       TYPE(SplineDef) :: Spline
00145       REAL(KIND=qPREC) :: GetSplineValue
00146       INTEGER i, j,k, n
00147       REAL(KIND=qPREC) :: h, u
00148       REAL(KIND=qPREC), DIMENSION(:), POINTER :: b,c,d,x,y
00149       b=>Spline%b
00150       c=>Spline%c
00151       d=>Spline%d
00152       x=>Spline%x
00153       y=>Spline%y
00154       n=size(x)
00155       IF (n == 1) GetSplineValue=y(1)
00156 
00157       IF(u <= x(1)) THEN
00158          GetSplineValue = y(1)
00159          RETURN
00160       END IF
00161       IF(u >= x(n)) THEN
00162          GetSplineValue = y(n)
00163          RETURN
00164       END IF
00165 
00166       i = 1
00167       j = n+1
00168       DO WHILE (j > i+1)
00169          k = (i+j)/2
00170          IF(u < x(k)) THEN
00171             j=k
00172          ELSE
00173             i=k
00174          END IF
00175       END DO
00176 
00177       h = u - x(i)
00178       GetSplineValue = y(i) + h*(b(i) + h*(c(i) + h*d(i)))
00179    END FUNCTION GetSplineValue
00180 
00181 END MODULE Splines
 All Classes Files Functions Variables