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