Statistiques
| Révision :

root / src / Hupdate_all.f90 @ 4

Historique | Voir | Annoter | Télécharger (3,47 ko)

1 1 equemene
 subroutine hupdate_all(n, ds, dgrad, Hess)
2 1 equemene
3 4 pfleura2
  use Path_module, only :  HUpdate,Hinv
4 1 equemene
  use Io_module, only : IoOut
5 1 equemene
6 1 equemene
  IMPLICIT NONE
7 1 equemene
8 1 equemene
  INTEGER, PARAMETER :: KINT=KIND(1)
9 1 equemene
  INTEGER, PARAMETER :: KREAL=KIND(1.0D0)
10 1 equemene
11 1 equemene
  INTEGER(KINT), INTENT(IN) :: n
12 1 equemene
  real(KREAL), INTENT(IN) :: ds(n), dgrad(n)
13 1 equemene
  REAL(KREAL),INTENT(INOUT) :: Hess(n:n)
14 1 equemene
15 1 equemene
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16 1 equemene
!
17 4 pfleura2
! This subroutine is a driver for the update of the Hessian
18 1 equemene
! in the Quasi-Newton optimization step.
19 1 equemene
! We propose mostly update of the inverse Hessian.
20 1 equemene
!
21 1 equemene
! input:
22 1 equemene
!  - n: number of coordinates
23 1 equemene
!  - ds: GeomOld-Geom (expressed in the optimization coord.)
24 1 equemene
!  - dgrad: GradOld-Grad  (expressed in the optimization coord.)
25 1 equemene
!
26 1 equemene
!  input/output:
27 1 equemene
!  - Hess: old Hessian in input, replaced by updated Hessian on output
28 1 equemene
!
29 1 equemene
!!!!!!!!!!!!!!!!
30 1 equemene
!
31 1 equemene
! For now, we propose :
32 4 pfleura2
! - For update of inverse Hessian:
33 4 pfleura2
!      - BFGS (Broyden, Fletcher, Goldfarb, and Shanno)
34 4 pfleura2
!      - DFP (Davidson - Fletcher - Powell )
35 4 pfleura2
!      - MS (Murtagh-Sargent)
36 1 equemene
!
37 1 equemene
! - For update of hessian:
38 1 equemene
!      - Bofill (combination of Murtagh-Sargent and Powell-symmetric-Broyden
39 4 pfleura2
!      - MS (Murtagh-Sargent)
40 4 pfleura2
!      - BFGS
41 4 pfleura2
!      - DFP
42 4 pfleura2
!----------------------------------------------------
43 4 pfleura2
! We use the fact that some formula are the dual of other
44 4 pfleura2
! that is, one can convert the formula for update of the Hessian
45 4 pfleura2
! into the update of the inverse Hessian by replacing the Hessian
46 4 pfleura2
! (denoted by B in the Fletcher book and in the Nocedal&Wrigth book)
47 4 pfleura2
! by the inverse hessian (denoted by H), and replacing ds by dgrad in the
48 4 pfleura2
! formula.
49 4 pfleura2
! See for example:
50 4 pfleura2
! 1) Numerical Optimization, Springer 2nd Ed., 2006
51 4 pfleura2
! by Jorge Nocedal & Stephen J. Wright
52 4 pfleura2
! 2) ractical Methods of Optimization, John Wiley & Sons, second ed., 1987.
53 4 pfleura2
! by R. Fletcher
54 1 equemene
!
55 4 pfleura2
!
56 1 equemene
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
57 1 equemene
58 1 equemene
59 1 equemene
60 1 equemene
  interface
61 1 equemene
     function valid(string) result (isValid)
62 1 equemene
       logical                  :: isValid
63 1 equemene
       character(*), intent(in) :: string
64 1 equemene
     end function valid
65 1 equemene
66 1 equemene
  end interface
67 1 equemene
68 1 equemene
  ! ======================================================================
69 1 equemene
  LOGICAL :: Debug
70 1 equemene
71 1 equemene
  ! ======================================================================
72 1 equemene
73 4 pfleura2
  debug = valid ('HESS') .OR. Valid('HUPDATE')
74 1 equemene
75 1 equemene
  if (debug) WRITE(*,*) "================================= Entering Hupdate_all ==============="
76 1 equemene
77 1 equemene
  SELECT CASE (Hinv)
78 1 equemene
     CASE (.TRUE.)
79 4 pfleura2
        SELECT CASE (HUpdate)
80 1 equemene
           CASE ("BFGS")
81 1 equemene
              Call Hinvup_BFGS(n,ds,dgrad,hess)
82 4 pfleura2
           CASE ("MS")
83 4 pfleura2
! we use the fact that MS is self-dual
84 4 pfleura2
              Call Hupdate_MS(n,dgrad,ds,hess)
85 4 pfleura2
           CASE ("DFP")
86 4 pfleura2
              Call Hinvup_DFP(n,ds,dgrad,hess)
87 1 equemene
           CASE DEFAULT
88 4 pfleura2
              WRITE(*,*) Trim(HUpdate) // " not known: using BFGS"
89 1 equemene
              Call Hinvup_BFGS(n,ds,dgrad,hess)
90 1 equemene
         END SELECT
91 1 equemene
     CASE (.FALSE.)
92 4 pfleura2
        SELECT CASE (HUpdate)
93 1 equemene
           CASE ("BOFILL")
94 1 equemene
              Call Hupdate_Bofill(n,ds,dgrad,hess)
95 1 equemene
           CASE ("MS")
96 1 equemene
              Call Hupdate_MS(n,ds,dgrad,hess)
97 4 pfleura2
           CASE ("DFP")
98 4 pfleura2
! we use the fact that DFP is the dual of BFGS
99 4 pfleura2
              Call Hinvup_BFGS(n,dgrad,ds,hess)
100 4 pfleura2
           CASE ("BFGS")
101 4 pfleura2
! we use the fact that DFP is the dual of BFGS
102 4 pfleura2
              Call Hinvup_DFP(n,dgrad,ds,hess)
103 1 equemene
           CASE DEFAULT
104 4 pfleura2
              WRITE(*,*) Trim(HUpdate) // " not known: using Bofill"
105 1 equemene
              Call Hupdate_Bofill(n,ds,dgrad,hess)
106 1 equemene
           END SELECT
107 1 equemene
    END SELECT
108 1 equemene
109 4 pfleura2
  if (debug) WRITE(*,*) "================================= Exiting Hupdate_all ==============="
110 1 equemene
end subroutine hupdate_all