root / src / Hupdate_all.f90
Historique | Voir | Annoter | Télécharger (5,59 ko)
1 |
subroutine hupdate_all(n, ds, dgrad, Hess) |
---|---|
2 |
|
3 |
!---------------------------------------------------------------------- |
4 |
! Copyright 2003-2014 Ecole Normale Supérieure de Lyon, |
5 |
! Centre National de la Recherche Scientifique, |
6 |
! Université Claude Bernard Lyon 1. All rights reserved. |
7 |
! |
8 |
! This work is registered with the Agency for the Protection of Programs |
9 |
! as IDDN.FR.001.100009.000.S.P.2014.000.30625 |
10 |
! |
11 |
! Authors: P. Fleurat-Lessard, P. Dayal |
12 |
! Contact: optnpath@gmail.com |
13 |
! |
14 |
! This file is part of "Opt'n Path". |
15 |
! |
16 |
! "Opt'n Path" is free software: you can redistribute it and/or modify |
17 |
! it under the terms of the GNU Affero General Public License as |
18 |
! published by the Free Software Foundation, either version 3 of the License, |
19 |
! or (at your option) any later version. |
20 |
! |
21 |
! "Opt'n Path" is distributed in the hope that it will be useful, |
22 |
! but WITHOUT ANY WARRANTY; without even the implied warranty of |
23 |
! |
24 |
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
25 |
! GNU Affero General Public License for more details. |
26 |
! |
27 |
! You should have received a copy of the GNU Affero General Public License |
28 |
! along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>. |
29 |
! |
30 |
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr, |
31 |
! for commercial licensing opportunities. |
32 |
!---------------------------------------------------------------------- |
33 |
|
34 |
use Path_module, only : HUpdate,Hinv |
35 |
use Io_module, only : IoOut |
36 |
|
37 |
IMPLICIT NONE |
38 |
|
39 |
INTEGER, PARAMETER :: KINT=KIND(1) |
40 |
INTEGER, PARAMETER :: KREAL=KIND(1.0D0) |
41 |
|
42 |
INTEGER(KINT), INTENT(IN) :: n |
43 |
real(KREAL), INTENT(IN) :: ds(n), dgrad(n) |
44 |
REAL(KREAL),INTENT(INOUT) :: Hess(n:n) |
45 |
|
46 |
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
47 |
! |
48 |
! This subroutine is a driver for the update of the Hessian |
49 |
! in the Quasi-Newton optimization step. |
50 |
! We propose mostly update of the inverse Hessian. |
51 |
! |
52 |
! input: |
53 |
! - n: number of coordinates |
54 |
! - ds: GeomOld-Geom (expressed in the optimization coord.) |
55 |
! - dgrad: GradOld-Grad (expressed in the optimization coord.) |
56 |
! |
57 |
! input/output: |
58 |
! - Hess: old Hessian in input, replaced by updated Hessian on output |
59 |
! |
60 |
!!!!!!!!!!!!!!!! |
61 |
! |
62 |
! For now, we propose : |
63 |
! - For update of inverse Hessian: |
64 |
! - BFGS (Broyden, Fletcher, Goldfarb, and Shanno) |
65 |
! - DFP (Davidson - Fletcher - Powell ) |
66 |
! - MS (Murtagh-Sargent) |
67 |
! |
68 |
! - For update of hessian: |
69 |
! - Bofill (combination of Murtagh-Sargent and Powell-symmetric-Broyden |
70 |
! - MS (Murtagh-Sargent) |
71 |
! - BFGS |
72 |
! - DFP |
73 |
!---------------------------------------------------- |
74 |
! We use the fact that some formula are the dual of other |
75 |
! that is, one can convert the formula for update of the Hessian |
76 |
! into the update of the inverse Hessian by replacing the Hessian |
77 |
! (denoted by B in the Fletcher book and in the Nocedal&Wrigth book) |
78 |
! by the inverse hessian (denoted by H), and replacing ds by dgrad in the |
79 |
! formula. |
80 |
! See for example: |
81 |
! 1) Numerical Optimization, Springer 2nd Ed., 2006 |
82 |
! by Jorge Nocedal & Stephen J. Wright |
83 |
! 2) ractical Methods of Optimization, John Wiley & Sons, second ed., 1987. |
84 |
! by R. Fletcher |
85 |
! |
86 |
! |
87 |
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
88 |
|
89 |
|
90 |
|
91 |
interface |
92 |
function valid(string) result (isValid) |
93 |
logical :: isValid |
94 |
character(*), intent(in) :: string |
95 |
end function valid |
96 |
|
97 |
|
98 |
SUBROUTINE die(routine, msg, file, line, unit) |
99 |
|
100 |
Use VarTypes |
101 |
Use io_module |
102 |
|
103 |
implicit none |
104 |
!--------------------------------------------------------------- Input Variables |
105 |
character(len=*), intent(in) :: routine, msg |
106 |
character(len=*), intent(in), optional :: file |
107 |
integer(KINT), intent(in), optional :: line, unit |
108 |
|
109 |
END SUBROUTINE die |
110 |
|
111 |
SUBROUTINE Warning(routine, msg, file, line, unit) |
112 |
|
113 |
Use VarTypes |
114 |
Use io_module |
115 |
|
116 |
implicit none |
117 |
|
118 |
character(len=*), intent(in) :: routine, msg |
119 |
character(len=*), intent(in), optional :: file |
120 |
integer(KINT), intent(in), optional :: line, unit |
121 |
|
122 |
END SUBROUTINE Warning |
123 |
|
124 |
|
125 |
end interface |
126 |
|
127 |
! ====================================================================== |
128 |
LOGICAL :: Debug |
129 |
|
130 |
! ====================================================================== |
131 |
|
132 |
debug = valid ('HESS') .OR. Valid('HUPDATE') |
133 |
|
134 |
if (debug) WRITE(*,*) "================================= Entering Hupdate_all ===============" |
135 |
|
136 |
SELECT CASE (Hinv) |
137 |
CASE (.TRUE.) |
138 |
SELECT CASE (HUpdate) |
139 |
CASE ("BFGS") |
140 |
Call Hinvup_BFGS(n,ds,dgrad,hess) |
141 |
CASE ("MS") |
142 |
! we use the fact that MS is self-dual |
143 |
Call Hupdate_MS(n,dgrad,ds,hess) |
144 |
CASE ("DFP") |
145 |
Call Hinvup_DFP(n,ds,dgrad,hess) |
146 |
CASE DEFAULT |
147 |
Call Warning("Hupdate_all Hinv=T",Trim(HUpdate) // " not known: using BFGS" & |
148 |
,Unit=IoOut) |
149 |
Call Hinvup_BFGS(n,ds,dgrad,hess) |
150 |
END SELECT |
151 |
CASE (.FALSE.) |
152 |
SELECT CASE (HUpdate) |
153 |
CASE ("BOFILL") |
154 |
Call Hupdate_Bofill(n,ds,dgrad,hess) |
155 |
CASE ("MS") |
156 |
Call Hupdate_MS(n,ds,dgrad,hess) |
157 |
CASE ("DFP") |
158 |
! we use the fact that DFP is the dual of BFGS |
159 |
Call Hinvup_BFGS(n,dgrad,ds,hess) |
160 |
CASE ("BFGS") |
161 |
! we use the fact that DFP is the dual of BFGS |
162 |
Call Hinvup_DFP(n,dgrad,ds,hess) |
163 |
CASE DEFAULT |
164 |
Call Warning("Hupdate_all Hinv=F",Trim(HUpdate) // " not known: using BoFill" & |
165 |
,Unit=IoOut) |
166 |
Call Hupdate_Bofill(n,ds,dgrad,hess) |
167 |
END SELECT |
168 |
END SELECT |
169 |
|
170 |
if (debug) WRITE(*,*) "================================= Exiting Hupdate_all ===============" |
171 |
end subroutine hupdate_all |