src/res_fit.f

Fortran project RESTRAX, source module src/res_fit.f.

Source module last modified on Fri, 12 May 2006, 17:22;
HTML image of Fortran source automatically generated by for2html on Mon, 29 May 2006, 15:06.


#
# $Id: res_fit.f,v 1.4 2006/05/12 15:22:12 saroun Exp $  
#**************************************************************
      SUBROUTINE MATIN(array,norder,dt)
#**************************************************************
#  
#     ** INVERSION OF A SYMMETRIC MATRIX FROM BEVINGTON **
#   

      INCLUDE 'const.inc'


      double precision array,amax,save
      dimension array(mpar,mpar),ik(mpar),jk(mpar)
10    dt = 1.
11    do 100 k=1,norder
      amax = .0
21    do 30 i=k,norder
      do 30 j=k,norder
23    if(dabs(amax)-dabs(array(i,j))) 24,24,30
24    amax = array(i,j)
      ik(k) = i
      jk(k) = j
30    continue
31    if(amax) 41,32,41
32    dt = .0
      goto 140
41    i = ik(k)
      if(i-k) 21,51,43
43    do 50 j=1,norder
      save = array(k,j)
      array(k,j) = array(i,j)
50    array(i,j) = -save
51    j = jk(k)
      if(j-k) 21,61,53
53    do 60 i=1,norder
      save = array(i,k)
      array(i,k) = array(i,j)
60    array(i,j) = -save
61    do 70 i=1,norder
      if(i-k) 63,70,63
63    array(i,k) = -array(i,k)/amax
70    continue
71    do 80 i=1,norder
      if(i-k) 74,80,74
74    do 77 j = 1,norder
      if(j-k) 75,77,75
75    array(i,j) = array(i,j)+array(i,k)*array(k,j)
77    continue
80    continue
81    do 90 j=1,norder
      if(j-k) 83,90,83
83    array(k,j) = array(k,j)/amax
90    continue
      array(k,k) = 1./amax
100   dt = dt*amax
101   do 130 l=1,norder
      k = norder - l +1
      j = ik(k)
      if(j-k) 111,111,105
105   do 110 i=1,norder
      save = array(i,k)
      array(i,k) = -array(i,j)
110   array(i,j) = save
111   i = jk(k)
      if(i-k) 130,130,113
113   do 120 j=1,norder
      save = array(k,j)
      array(k,j) = -array(i,j)
120   array(i,j) = save
130   continue
140   return
      end


#**************************************************************
#      SUBROUTINE FDERIV(FUNCTN,X,A,DELTAA,NPT,NPTS,NTERMS,DERIV)
      SUBROUTINE FDERIV(deltaa,deriv)
# 14/4/2006, J.S.    
# rewritten to deal with common arrays defined in restrax.inc  
#/// Modified on 20.4.99 by J. Saroun:
#/// NPT contains partition of X,Y array to individual spectra
#**************************************************************
      implicit none
      INCLUDE 'const.inc'
      INCLUDE 'restrax.inc'

      real*4 deltaa(mpar),deriv(maxd,mpar)
      real*4 yfit1(maxd),yfit2(maxd),aj,delta
      integer*4 i,j    
      real*4 eval,HIST


      do j=1,nfpar
      if(deltaa(j).ne.0) then
        aj = fpar(j)
        delta = deltaa(j)
        fpar(j) = aj+delta
        eval = hist(xhist,yfit1,npt,npt(mdat),fpar,nfpar)
        fpar(j) = aj-delta
        eval = hist(xhist,yfit2,npt,npt(mdat),fpar,nfpar)

        do i = 1,npt(mdat)
          deriv(i,j) = (yfit1(i)-yfit2(i))/(2.*delta)
        enddo

        fpar(j) = aj
      else
        do i = 1,npt(mdat)
          deriv(i,j) = .0
        enddo
      endif
      enddo

      return
      end



#------------------------------------------------------------
      SUBROUTINE CURF3T(deltaa,flamda,nfree)
# 14/4/2006, J.S.    
# rewritten to deal with common arrays defined in restrax.inc  
#
#      SUBROUTINE CURF3T(FUNCTN,X,Y,SIGMAY,IPT,NPT,NPTS,NTERMS,MODE,
#     *         A,DELTAA,SIGMAA,COVAR,FLAMDA,YFIT,CHISQR,DCHISQ)
#
#  ****** CURF1T IS CURFIT MODIFIED TO ALLOW FOR COVARIANCE
#  ****** MATRIX OUTPUT
#  ****** CURF3T IS CURFIT evaluating FUNCTN in a single call
#
#/// Modified on 20.4.99 by J. Saroun:
#
#/// IPT contains indexes for data points, which are used as mask.
#/// NPT contains partition of X,Y array to individual spectra
#/// If IPT(i)<=0 , the point is excluded from fitting
#/// DCHISQ(i) now contains partial Chi^2 for individual spectra 
#------------------------------------------------------------

      implicit none
      
      INCLUDE 'const.inc'
      INCLUDE 'restrax.inc'
      
      real*4 deltaa(mpar),flamda,dt
      integer*4 nfree
      real*4 weight(maxd),alpha(mpar,mpar),beta(mpar),deriv(maxd,mpar),b(mpar)
      real*8 array(mpar,mpar)
      integer*4 ins(mpar),np
      real*4 eval,hist,fchisq,dy,chisq1,chi
      integer*4 i,j,k
   
      np=npt(mdat)
  
#///  set weigths:
  
      do  i=1,np
        if (ipt(i).le.0) then      ! weight=0 to ignore the point if IPT(i)=0
           weight(i)=0
        else 
          select case(mode)  
          case(-1)
            if(abs(spy(i)).le.1.) then
              weight(i) = 1.
            else
              weight(i) = 1./abs(spy(i))
            endif
          case(0)
            weight(i) = 1.
          case default
            weight(i) = 1./(spz(i)**2)
          end select
        endif
      enddo

#/// get derivatives:  
      eval = hist(xhist,rhist,nhist,nhist(mdat),fpar,nfpar)
      call FDERIV(deltaa,deriv)

#// get initial Chisq

#      CHISQ1 = FCHISQ(SPY,SPZ,IPT,NPT(MDAT),NFREE,MODE,RHIST,DCHISQ)  
# modified J.S., 5/5/2006:
# this call is not needed, it is sure that CHISQR has been calculated before
      chisq1=chisqr

#// construct the matrix
      do j=1,nfpar
         beta(j) = .0
         do  k=1,j
           alpha(j,k) = .0
         end do
      end do   
      do i=1,np
        dy = spy(i)-rhist(i)
        do j=1,nfpar
          beta(j) = beta(j)+weight(i)*dy*deriv(i,j)
          do k=1,j
            alpha(j,k) = alpha(j,k)+weight(i)*deriv(i,j)*deriv(i,k)
          end do
        end do  
      end do
#// make it symmetric
      do j=2,nfpar
        do k=1,j-1
          alpha(k,j) = alpha(j,k)
        end do
      end do

#// check if some diagonal element = 0, then keep it fixed 
      do  j=1,nfpar
          ins(j) = 0
          if(abs(alpha(j,j)).lt.1.e-19) ins(j)=1
      end do

#// start here the cycle with FLAMDA adjustment
      do while (flamda.le.1000.and.chisq1.gt.0)

#// get the normalized matrix, ARRAY
        do j=1,nfpar
          do k=1,nfpar
            if((ins(j)+ins(k)).eq.0) then
              array(j,k) = alpha(j,k)/sqrt(alpha(j,j)*alpha(k,k))
              covar(j,k) = array(j,k)
            else  ! set COVAR = delta(i,j) for i or j fixed 
              array(j,k) = .0
              covar(j,k) = .0
              if(j.eq.k) covar(j,k)=1.
            endif
          end do 
          array(j,j) = 1.+flamda
        end do

#// invert ARRAY      
        call MATIN(array,nfpar,dt)

#// calculate increments, B(i) 
        do j=1,nfpar
          b(j) = .0
          do k=1,nfpar
            if((ins(j)+ins(k)).eq.0) then
              b(j) = b(j)+beta(k)*array(j,k)/sqrt(alpha(j,j)*alpha(k,k))
            endif
          end do 
        end do

#// set B(i)= new estimated values
        do j=1,nfpar
          b(j) = b(j)+fpar(j)
        end do
      
#// check new ChiSq     
        eval = hist(xhist,rhist,nhist,nhist(mdat),b,nfpar)
        chisqr = FCHISQ(spy,spz,ipt,npt(mdat),nfree,mode,rhist,dchisq)

#// if ChiSq increases, try again with larger FLAMDA
        if (chisq1-chisqr.lt.0) then
           flamda = 10.*flamda 
        else
           chisq1=-1.  ! signal to stop cycle
        endif       
      enddo
      
#// finalization of an iteration step:
#// update FPAR(i) values and error estimates
      do j=1,nfpar
         fpar(j) = b(j)
         if(ins(j).eq.0) then
           chi = chisqr
           if(chisqr.lt.1) chi=1.
           sigmaa(j) = sqrt(abs(array(j,j))*chi*(1.+flamda)/alpha(j,j))
         else
           sigmaa(j) = .0
         endif
      enddo

#// decrease FLAMDA for next iteration
      flamda = flamda/10.

      end



#------------------------------------------------------------
#
      real*4 FUNCTION FCHISQ(y,sigmay,ipt,npts,nfree,mode,yfit,dchisq)
#/// Modified on 20.4.99 by J. Saroun:
#
#/// IPT contains indexes for data points, which are used as mask.
#/// If IPT(i)<=0 , the point is excluded from fitting
#/// DCHISQ(i) now contains partial Chi^2 for individual spectra 
#------------------------------------------------------------
#  
#       EVALUATES REDUCED CHISQUARE FOR FIT TO DATA (BEVINGTON)
#  
      implicit none
      INCLUDE 'const.inc'

      integer*4 npts,nfree,mode
      real*8 chisq, weight
      real*4 y(npts),sigmay(npts),yfit(npts)
      integer*4 ipt(npts),nf(0:mdat),nsum,j,i
      real*4 dchisq(mdat),z,res
      
      chisq=0.0
      res=0.0
      nsum=0
      do j=1,mdat
         dchisq(j)=0.0
         nf(j)=0.0
      enddo   
      if (nfree.le.0) goto 50      
      do i=1,npts
        j=ipt(i)
        if (j.le.0) then      ! weight=0 to ignore the point if IPT(i)=0
           weight=0.0
        else
          select case(mode)
          case(-1)
            if(abs(y(i)).gt.0.e0) then
              weight = 1./y(i)
            else
              weight = 1.
            endif
          case(0)
            weight = 1.      
          case default
            weight = 1./(sigmay(i)**2)
          end select
# 1     format('FCHISQ: i,Y,FIT: ',I3,2x,2(1x,G12.6))
#        write(*,*) I,Y(I),YFIT(I)
          z=weight*(y(i)-yfit(i))**2
          dchisq(j)=dchisq(j)+z
          nf(j)=nf(j)+1
          nsum=nsum+1
          chisq = chisq+z
        endif   
      enddo
      
      res = chisq/nfree
      do j=1,mdat
         if (nf(j).gt.0) dchisq(j)=dchisq(j)/nf(j)/nfree*nsum
      enddo   
50    FCHISQ=res
      return
      end
#
#
#------------------------------------------------------------
#
      SUBROUTINE GETPEAKPARAM(x,y,n,suma,center,fwhm,wspread)
# get peak parameters
#------------------------------------------------------------
      implicit none
      INCLUDE 'const.inc'
      
      integer*4 n,i
      real*4 x(*),y(*)
      
      real*8 s0,s1,s2,suma,center,fwhm,wspread,z
      suma=0.
      s0=0
      s1=0
      s2=0
      do i=1,n
        z=x(i)
        s0=s0+y(i)
        s1=s1+y(i)*z
        s2=s2+y(i)*z**2
      enddo
      suma=s0
      center=s1/s0
      wspread=sqrt8ln2*sqrt(s2/s0-center**2)
      call CALFWHM(x,y,n,fwhm)
       
      end

     
#-------------------------------------------------------------------
      SUBROUTINE CALFWHM(x,y,n,fwhm)
# Calculate fwhm      
#-------------------------------------------------------------------
      implicit none

      integer*4 n,i,imax,i1,i2    
      real*4 x(*),y(*)
      real*8 ymax,x1,x2,fwhm
      ymax=y(1)
      imax=1
      do i=1,n         
         if(ymax.le.y(i)) then
           imax=i
           ymax=y(i)
         endif  
      enddo
      i1=1
      do while ((y(i1).lt.ymax/2.d0).and.(i1.lt.n))
        i1=i1+1
      enddo
      i2=n
      do while ((y(i2).lt.ymax/2.d0).and.(i2.gt.0))
        i2=i2-1
      enddo
      if((i1.gt.1).and.(i1.lt.n).and.(i2.gt.1).and.(i2.lt.n).and.
     *   (i1.le.i2)) then
         x1=x(i1-1)+(ymax/2.d0-y(i1-1))/(y(i1)-y(i1-1))*(x(i1)-x(i1-1))
         x2=x(i2)+(ymax/2.d0-y(i2))/(y(i2+1)-y(i2))*(x(i2+1)-x(i2))
         fwhm=x2-x1
      else
         fwhm=0
      endif
      end

# $Log: res_fit.f,v $
# Revision 1.4  2006/05/12 15:22:12  saroun
# bug fix: DET structure passed to MATIN as the last argument instead of a dummy real
#