src/res_opt.f

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

Source module last modified on Wed, 10 May 2006, 16:36;
HTML image of Fortran source automatically generated by for2html on Mon, 29 May 2006, 15:06.


#//////////////////////////////////////////////////////////////////////
#////  $Id: res_opt.f,v 1.4 2006/05/10 14:36:26 saroun Exp $
#////
#////  R E S T R A X   (c) J. Saroun, J. Kulda, 1995-2005
#////
#////  Subroutines for TAS optimization
#////
#//////////////////////////////////////////////////////////////////////


#**************************************************************
      SUBROUTINE GETDERIV1(a,delta,functn,der1,der2)
#     1st and 2nd derivative for 1 parameter
#**************************************************************
      implicit none

      real*8 aj,y1,y2,y0,der1,der2
      real*4 a(1),delta
      real*4 functn
      external functn
      der1=0
      der2=0
      if(delta.ne.0) then
         aj = a(1)
10       a(1) = aj+delta
         y1 = functn(a(1))
         a(1) = aj-delta
         y2 = functn(a(1))
         if(abs(y1+y2).gt.1d-10.and.abs((y1-y2)/(y1+y2)).gt.0.4) then
            delta=delta/2.
            goto 10
         endif
         a(1) = aj            
         y0 = functn(a(1))
         der1=(y1-y2)/2/delta
         der2=(y1+y2-2*y0)/delta/delta
      endif           
      end

#**************************************************************
      real*8 FUNCTION GETDERIV(ider,j,k,a,deltaa,n,functn)
# get derivatives of FUNCTN
# IDER=1,2 for 1st,2nd derivative
# J,K indexes of A(N) array pointing to the dependent variables
# DELTAA(N) ... increments for numerical calculation of derivatives
#**************************************************************
      implicit none

      integer*4 ider,j,k,n
      real*8 aj,ak,y1,y2,y0,w1,w2,delta,deltb
      real*4 a(n),deltaa(n)
      real*4 functn
      external functn

      if (ider.eq.0) then
        GETDERIV=functn(a)
      else if (ider.eq.1) then
        if(deltaa(j).ne.0) then
          aj = a(j)
10        delta = deltaa(j)
          a(j) = aj+delta
          y1 = functn(a)
          a(j) = aj-delta
          y2 = functn(a)
          if(abs(y1+y2).gt.1d-10.and.abs((y1-y2)/(y1+y2)).gt.0.4) then
            deltaa(j)=deltaa(j)/2.
            goto 10
          endif
          GETDERIV=(y1-y2)/2/delta
          a(j)=aj          
        else
          GETDERIV=0      
        endif
      else if (ider.eq.2) then
        if(deltaa(j).ne.0.and.deltaa(k).ne.0) then
         aj = a(j)
         ak = a(k)
         delta = deltaa(j)
         deltb = deltaa(k)
         if (j.eq.k) then
           a(j) = aj+delta
           y1 = functn(a)
           a(j) = aj-delta
           y2 = functn(a)
           a(j) = aj            
           y0 = functn(a)
           GETDERIV=(y1+y2-2*y0)/delta/delta
         else
           a(k) = ak+deltb
           a(j) = aj+delta
           y1 = functn(a)
           a(j) = aj-delta
           y2 = functn(a)
           w1=(y1-y2)/2/delta
           a(k) = ak-deltb
           a(j) = aj+delta
           y1 = functn(a)
           a(j) = aj-delta
           y2 = functn(a)
           w2=(y1-y2)/2/delta
           GETDERIV=(w1-w2)/2/deltb
         endif
         a(k) = ak
         a(j) = aj
        else
         GETDERIV=0
        endif 
      endif     
      
#
      return
      end
#
#------------------------------------------------------------
      SUBROUTINE LMOPT(functn,a,ami,ama,np,tol,isil)
#  find minimum of FUNCTN function by Levenberg-Marquardt algorithm
#------------------------------------------------------------
      implicit none
      
      INCLUDE 'const.inc'
      
      integer*4 maxit,j,k,np,ins(mpar),it,ierr,isil
      real*4 tol,flamda
      real*8 array(mpar,mpar),alpha(mpar,mpar),beta(mpar),an(mpar,mpar)           
      real*4 a(np),ami(np),ama(np),b(mpar),deltaa(mpar)
      real*4 chisqr,chisq1,db,ksi
      real*8 GETDERIV
      real*4 functn
      common /error/ierr
      external functn
1     format( 'CHISQR: ',g12.5,2x, 'PAR: ',6(g12.5,2x))
2     format( 'x',$)

#/// initializations
      ierr=0
      it=0       ! iter. counter
      maxit=100  ! max. 100 iterations
      do j=1,np
        deltaa(j)=abs(a(j)*tol)  ! derive increments from TOL
      enddo  
      flamda=0.1

#/// Get initial Chi^2
      chisqr = functn(a)

#/// Start main cycle      
41    chisq1 = chisqr  

      if (isil.le.0) write(*,1) chisqr,(a(j),j=1,min(np,6))
 
#/// Calculate derivatives: 
      if(np.eq.1) then       
        call GETDERIV1(a(1),deltaa(1),functn,beta(1),alpha(1,1))
      else
        do j=1,np
            beta(j)=GETDERIV(1,j,j,a,deltaa,np,functn)
            do k=1,j
              alpha(j,k)=GETDERIV(2,j,k,a,deltaa,np,functn)
            end do
        end do
      endif  

#/// symetrize ALPHA
      do j=1,np
        do k=1,j
          alpha(k,j) = alpha(j,k)
        end do
      end do

#/// check zeros on digonal
      do  j=1,np
          ins(j) = 0
          if(abs(alpha(j,j)).lt.1.e-19) ins(j)=1
      end do

#/// create Hessian matrix
71    do j=1,np
        do k=1,np
          if((ins(j)+ins(k)).eq.0) then
            an(j,k)=sqrt(abs(alpha(j,j))*abs(alpha(k,k)))
            array(j,k) = alpha(j,k)/an(j,k)
          else
            an(j,k)=1.d0
            array(j,k) = 0.d0
          endif
        end do 
        array(j,j) = 1.d0+flamda
      end do
      
#/// invert Hessian matrix
      call INVERT(np,array,mpar,array,mpar)

#/// increment with limits checking
      ksi=1.0 
10    do j=1,np
        db = 0.0
        do k=1,np
          if((ins(j)+ins(k)).eq.0) then
            db=db-ksi*beta(k)*array(j,k)/an(j,k)
          endif
        end do
        if (a(j)+db.gt.ama(j).or.a(j)+db.lt.ami(j)) then
          ksi=ksi/2.0
#      write(*,1) KSI,A(J),A(J)+DB
          goto 10
        endif  
        b(j)=a(j)+db       
      end do
      
#// get new CHISQR            
      chisqr = functn(b)  

#// check stop conditions
      if(chisqr.lt.1e-20) goto 110  ! exact match ... speed-up linear fit
      if ((chisq1-chisqr).le.0.d0) then 
        flamda = 5.*flamda
        if(flamda.le.1000.) then  ! try again with larger FLAMDA
           goto 71
        else
           ierr=-1
           goto 111   ! exit if FLAMBDA is too large
        endif      
      endif

#//  accept new values and continue 
110   do j=1,np
         a(j) = b(j)
      end do
      flamda = flamda/5.
      it=it+1

111   if (chisqr.lt.1e-20) ierr=1  ! exit on exact match
      if (abs(abs(chisq1/chisqr)-1.).lt.0.01) ierr=2  ! exit on < 1/100 change
      if (it.ge.maxit) ierr=-2  ! exit on maximum iteration number

      if (ierr.eq.0) goto 41

      end

#--------------------------------------------------------------------------------
# $Log: res_opt.f,v $
# Revision 1.4  2006/05/10 14:36:26  saroun
# Added ISIL argument to control printed output (ISIL=1 ... silent)
#
# Revision 1.3  2005/07/19 20:38:08  saroun
# *** empty log message ***
#
# Revision 1.2  2005/07/16 16:46:06  saroun
# Improved TAS optimization (LMOPT, OPINSTR, ...)
# Added limits of free parameters in LMOPT, etc
# Added demo with optimization in ./demo/opt
#