src/sim_opt.f

Fortran project SIMRES, source module src/sim_opt.f.

Source module last modified on Sun, 27 Mar 2005, 23:48;
HTML image of Fortran source automatically generated by for2html on Mon, 23 May 2005, 21:29.


#//////////////////////////////////////////////////////////////////////
#////
#////  R E S T R A X   4.4
#////
#////  Subroutines for 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
#      write(*,*) 'GETDERIV1 ',DELTA
      if(delta.ne.0) then
         aj = a(1)
10       a(1) = aj+delta
#         write(*,*) 'GETDERIV1  before' 
         y1 = functn(a(1))
#         write(*,*) 'GETDERIV1 after'
         a(1) = aj-delta
         y2 = functn(a(1))
         if(abs((y1-y2)/(y1+y2)).gt.0.4) then
            delta=delta/2.
            write(*,*) delta,y1,y2,abs((y1-y2)/(y1+y2))
            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)
#
#**************************************************************
      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)/(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,np,tol,inca,isil)
#  find minimum of FUNCTN function by Levenberg-Marquardt algorithm
#  A(i)    ...  parameters
#  INAC(i) ...  minimum increments for numerical derivatives of A
#  TOL     ... tolerance indicator
#  ISIL    ... silence level (no output for ISIL>1)
#------------------------------------------------------------
      implicit none
      
      integer*4 m,np,isil
      parameter (m=16)
      real*4 a(np),tol,inca(np)
      
      integer*4 maxit,j,k,ins(m),it,ierr
      real*8 array(m,m),alpha(m,m),an(m,m),beta(m)            
      real*4 b(m),deltaa(m),chisqr,chisq1,flamda
      real*8 GETDERIV
      real*4 functn
      common /error/ierr
      external functn

1     format(/,a10, ': ',i3,2x,8(g10.4,2x),2x)
#       write(*,*) ISIL
#       pause
      if (m.lt.np) goto 999
      ierr=0
      maxit=20   ! do no more than 20 iterations
      it=0

#/// set increments
      do j=1,np
        deltaa(j)=abs(a(j)*tol)
        if(deltaa(j).lt.inca(j)) deltaa(j)=inca(j)
      enddo  
      flamda=0.1
#      write(*,1) 'Delta',(DELTAA(J),J=1,NP),NP
#      pause       
      
#/// get initial CHISQR
      chisqr = functn(a(1))  
      
#/// start iterations      
41    chisq1=chisqr

      
      if (isil.lt.1)  write(*,1)  'IT,CHI,PAR',it,chisqr,(a(j),j=1,np)      

#/// get 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)
        end do
        do j=1,np
            do k=1,j
              alpha(j,k)=GETDERIV(2,j,k,a,deltaa,np,functn)
            end do
        end do
      endif  
#      write(*,1) 'dF/dRo',(BETA(J),J=1,NP),(DELTAA(J),J=1,NP)       
#      pause       

#/// 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.
            array(j,k) = 0.d0
          endif
        end do 
        array(j,j) = 1.d0+flamda
      end do
      
      
#      write(*,*)       
#        WRITE(*,1) 'BETA: ',CHISQR,FLAMDA,(BETA(K),K=1,NP)
#      DO J=1,NP
#        WRITE(*,1) 'ALPHA: ',(ALPHA(J,K),K=1,NP)
#      ENDDO  
#       DO J=1,NP
#         WRITE(*,1) 'ARRAY: ',(ARRAY(J,K),K=1,NP)
#       ENDDO  
#      write(*,*)
      
#// invert ARRAY
      call INVERT(np,array,m,array,m)
#      write(*,*) 'INVERT OK'

#// calculate new values
      do j=1,np
        b(j) = 0.d0
        do k=1,np
          if((ins(j)+ins(k)).eq.0) then
            b(j) = b(j)-beta(k)*array(j,k)/an(j,k)
          endif
        end do 
      end do
      
      do j=1,np
         b(j) = b(j)+a(j)
      end do

#// get new CHISQR            
      chisqr = functn(b(1))  

#// 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 FLAMBDA
           goto 71
        else
           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

      if (chisqr.lt.1e-20) goto 111   ! exit on exact match
      if (abs(abs(chisq1/chisqr)-1.).lt.0.001) goto 111  ! exit on < 1/1000 change
      if (it.ge.maxit) goto 112   ! exit on maximum iteration number

      goto 41

111   continue
      if (isil.lt.1)  write(*,1)  'IT,CHI,PAR',it,chisqr,(a(j),j=1,np)      
      return
112   if (isil.lt.2) then
        write(*,*)  'LMOPT: iteration > ',maxit, ', fit not finished'
      endif 
      return
999   if (isil.lt.2) then      
        write(*,*)  'ERROR in LMOPT: max. dimension ',m, '<',np 
      endif
      end