root / src / Hinvup_DFP.f90
Historique | Voir | Annoter | Télécharger (5,31 ko)
1 | 1 | pfleura2 | SUBROUTINE Hinvup_DFP(Nfree,ds,dgrad,Hess) |
---|---|---|---|
2 | 1 | pfleura2 | |
3 | 12 | pfleura2 | !---------------------------------------------------------------------- |
4 | 12 | pfleura2 | ! Copyright 2003-2014 Ecole Normale Supérieure de Lyon, |
5 | 12 | pfleura2 | ! Centre National de la Recherche Scientifique, |
6 | 12 | pfleura2 | ! Université Claude Bernard Lyon 1. All rights reserved. |
7 | 12 | pfleura2 | ! |
8 | 12 | pfleura2 | ! This work is registered with the Agency for the Protection of Programs |
9 | 12 | pfleura2 | ! as IDDN.FR.001.100009.000.S.P.2014.000.30625 |
10 | 12 | pfleura2 | ! |
11 | 12 | pfleura2 | ! Authors: P. Fleurat-Lessard, P. Dayal |
12 | 12 | pfleura2 | ! Contact: optnpath@gmail.com |
13 | 12 | pfleura2 | ! |
14 | 12 | pfleura2 | ! This file is part of "Opt'n Path". |
15 | 12 | pfleura2 | ! |
16 | 12 | pfleura2 | ! "Opt'n Path" is free software: you can redistribute it and/or modify |
17 | 12 | pfleura2 | ! it under the terms of the GNU Affero General Public License as |
18 | 12 | pfleura2 | ! published by the Free Software Foundation, either version 3 of the License, |
19 | 12 | pfleura2 | ! or (at your option) any later version. |
20 | 12 | pfleura2 | ! |
21 | 12 | pfleura2 | ! "Opt'n Path" is distributed in the hope that it will be useful, |
22 | 12 | pfleura2 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of |
23 | 12 | pfleura2 | ! |
24 | 12 | pfleura2 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
25 | 12 | pfleura2 | ! GNU Affero General Public License for more details. |
26 | 12 | pfleura2 | ! |
27 | 12 | pfleura2 | ! You should have received a copy of the GNU Affero General Public License |
28 | 12 | pfleura2 | ! along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>. |
29 | 12 | pfleura2 | ! |
30 | 12 | pfleura2 | ! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr, |
31 | 12 | pfleura2 | ! for commercial licensing opportunities. |
32 | 12 | pfleura2 | !---------------------------------------------------------------------- |
33 | 12 | pfleura2 | |
34 | 1 | pfleura2 | IMPLICIT NONE |
35 | 1 | pfleura2 | |
36 | 1 | pfleura2 | INTEGER, PARAMETER :: KINT=KIND(1) |
37 | 1 | pfleura2 | INTEGER, PARAMETER :: KREAL=KIND(1.D0) |
38 | 1 | pfleura2 | |
39 | 1 | pfleura2 | INTEGER(KINT), INTENT(IN) :: NFREE |
40 | 1 | pfleura2 | REAL(KREAL), INTENT(INOUT) :: Hess(Nfree,Nfree) |
41 | 1 | pfleura2 | REAL(KREAL), INTENT(IN) :: ds(Nfree) |
42 | 1 | pfleura2 | REAL(KREAL), INTENT(IN) :: dGrad(Nfree) |
43 | 1 | pfleura2 | |
44 | 1 | pfleura2 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
45 | 1 | pfleura2 | ! |
46 | 1 | pfleura2 | ! This subroutine updates the inverse of the Hessian |
47 | 1 | pfleura2 | ! using a DFP formula. |
48 | 1 | pfleura2 | ! |
49 | 1 | pfleura2 | ! input: |
50 | 1 | pfleura2 | ! - nfree: number of coordinates |
51 | 1 | pfleura2 | ! - ds: GeomOld-Geom (expressed in the optimization coord.) |
52 | 1 | pfleura2 | ! - dgrad: GradOld-Grad (expressed in the optimization coord.) |
53 | 1 | pfleura2 | ! |
54 | 1 | pfleura2 | ! input/output: |
55 | 1 | pfleura2 | ! - Hess: old Hessian in input, replaced by updated Hessian on output |
56 | 1 | pfleura2 | ! |
57 | 1 | pfleura2 | ! Notations for the formula: |
58 | 1 | pfleura2 | ! s=Geom-GeomOld |
59 | 1 | pfleura2 | ! y=Grad-GradOld |
60 | 1 | pfleura2 | ! W= (H)-1 ie the inverse matrix that we are updating |
61 | 1 | pfleura2 | ! W+=W + s(sT)/(s,y)-[Wy(yT)W]/(y,Wy) |
62 | 1 | pfleura2 | ! |
63 | 1 | pfleura2 | ! We decompose the computation in two steps: |
64 | 1 | pfleura2 | ! first we compute the update, and then we add it to Hess to get the |
65 | 1 | pfleura2 | ! new Hessian. |
66 | 1 | pfleura2 | ! This allows us to introduce a scaling factor if needed later. |
67 | 1 | pfleura2 | ! |
68 | 1 | pfleura2 | !!!!!!!!!!!!!!!! |
69 | 1 | pfleura2 | ! |
70 | 1 | pfleura2 | ! Due to the dual properties, this routine is also called to compute |
71 | 1 | pfleura2 | ! the BFGS update for the Hessian |
72 | 1 | pfleura2 | ! |
73 | 1 | pfleura2 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
74 | 1 | pfleura2 | |
75 | 1 | pfleura2 | |
76 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: HessOld(:,:) , Wy(:) |
77 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: Hup(:,:) |
78 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: Num(:,:), Num2(:,:) |
79 | 1 | pfleura2 | REAL(KREAL) :: Den, Den2,NormDS |
80 | 2 | pfleura2 | INTEGER(KINT) :: i, j, Idx |
81 | 1 | pfleura2 | |
82 | 1 | pfleura2 | LOGICAL :: debug |
83 | 1 | pfleura2 | |
84 | 1 | pfleura2 | interface |
85 | 1 | pfleura2 | function valid(string) result (isValid) |
86 | 1 | pfleura2 | logical :: isValid |
87 | 1 | pfleura2 | character(*), intent(in) :: string |
88 | 1 | pfleura2 | end function valid |
89 | 1 | pfleura2 | |
90 | 1 | pfleura2 | end interface |
91 | 1 | pfleura2 | |
92 | 1 | pfleura2 | debug=valid("Hinvup_DFP").OR. Valid('HUPDATE') |
93 | 1 | pfleura2 | |
94 | 1 | pfleura2 | if (debug) WRITE(*,*) "=========================== Entering Hinvup_DFP ====================" |
95 | 1 | pfleura2 | |
96 | 1 | pfleura2 | ALLOCATE(HessOld(nfree,nfree), Hup(nfree,nfree), & |
97 | 1 | pfleura2 | Num(nfree,nfree), Num2(nfree,nfree),Wy(Nfree) ) |
98 | 1 | pfleura2 | |
99 | 1 | pfleura2 | |
100 | 1 | pfleura2 | |
101 | 1 | pfleura2 | HessOld=Hess |
102 | 1 | pfleura2 | |
103 | 1 | pfleura2 | if (debug) THEN |
104 | 1 | pfleura2 | WRITE(*,*) 'DBG Hinvup_DFP HessOld' |
105 | 1 | pfleura2 | Idx=min(12,Nfree) |
106 | 1 | pfleura2 | DO I=1,NFree |
107 | 1 | pfleura2 | WRITE(*,'(12(1X,F8.4))') HessOld(I,1:Idx) |
108 | 1 | pfleura2 | END DO |
109 | 1 | pfleura2 | WRITE(*,*) 'DBG Hinvup_DFP ds' |
110 | 1 | pfleura2 | WRITE(*,'(12(1X,F8.4))') ds(1:Idx) |
111 | 1 | pfleura2 | WRITE(*,*) 'DBG Hinvup_DFP dGrad' |
112 | 1 | pfleura2 | WRITE(*,'(12(1X,F8.4))') dGrad(1:Idx) |
113 | 1 | pfleura2 | END IF |
114 | 1 | pfleura2 | |
115 | 1 | pfleura2 | |
116 | 1 | pfleura2 | |
117 | 1 | pfleura2 | NormDS=sqrt(DOT_PRODUCT(DS,DS)) |
118 | 1 | pfleura2 | if (NormDS>=1e-6) THEN |
119 | 1 | pfleura2 | IF (debug) WRITE(*,*) 'NormDS=',NormDS |
120 | 1 | pfleura2 | IF (debug) WRITE(*,*) 'DS=',DS |
121 | 1 | pfleura2 | IF (debug) WRITE(*,*) 'DGrad=',DGrad |
122 | 1 | pfleura2 | ! We calculate the denominator (DGrad,DS) |
123 | 1 | pfleura2 | den=DOT_PRODUCT(DGrad,DS) |
124 | 1 | pfleura2 | |
125 | 1 | pfleura2 | if (debug) WRITE(*,*) "(y,s)=(yT)s",den |
126 | 1 | pfleura2 | |
127 | 1 | pfleura2 | ! We calculate Wy |
128 | 1 | pfleura2 | Wy=0 |
129 | 1 | pfleura2 | DO I=1,nfree |
130 | 1 | pfleura2 | DO J=1,NFree |
131 | 1 | pfleura2 | Wy(I)=Wy(I)+HessOld(i,j)*DGrad(j) |
132 | 1 | pfleura2 | END DO |
133 | 1 | pfleura2 | END DO |
134 | 1 | pfleura2 | |
135 | 1 | pfleura2 | if (debug) WRITE(*,*) "Wy=",Wy |
136 | 1 | pfleura2 | |
137 | 1 | pfleura2 | ! We calculate Num=s(sT) |
138 | 1 | pfleura2 | Num=0. |
139 | 1 | pfleura2 | DO I=1,Nfree |
140 | 1 | pfleura2 | DO J=1,Nfree |
141 | 1 | pfleura2 | Num(I,J)=ds(i)*Ds(j) |
142 | 1 | pfleura2 | END DO |
143 | 1 | pfleura2 | END DO |
144 | 1 | pfleura2 | |
145 | 1 | pfleura2 | if (debug) THEN |
146 | 1 | pfleura2 | WRITE(*,*) "s(sT)" |
147 | 1 | pfleura2 | DO I=1,NFree |
148 | 1 | pfleura2 | WRITE(*,*) Num(I,:) |
149 | 1 | pfleura2 | END DO |
150 | 1 | pfleura2 | END if |
151 | 1 | pfleura2 | |
152 | 1 | pfleura2 | ! We calculate Num2=Wy(WyT) |
153 | 1 | pfleura2 | Num2=0 |
154 | 1 | pfleura2 | Do I=1,NFree |
155 | 1 | pfleura2 | DO J=1,NFree |
156 | 1 | pfleura2 | Num2(I,J)=Wy(i)*Wy(j) |
157 | 1 | pfleura2 | END DO |
158 | 1 | pfleura2 | END DO |
159 | 1 | pfleura2 | |
160 | 1 | pfleura2 | if (debug) THEN |
161 | 1 | pfleura2 | WRITE(*,*) "Wy(WyT)" |
162 | 1 | pfleura2 | DO I=1,NFree |
163 | 1 | pfleura2 | WRITE(*,*) Num2(I,:) |
164 | 1 | pfleura2 | END DO |
165 | 1 | pfleura2 | END if |
166 | 1 | pfleura2 | |
167 | 1 | pfleura2 | |
168 | 1 | pfleura2 | ! we calculate den=(sT)y |
169 | 1 | pfleura2 | den=DOT_PRODUCT(ds,dgrad) |
170 | 1 | pfleura2 | |
171 | 1 | pfleura2 | ! we calculate den2=(yT)Wy |
172 | 1 | pfleura2 | den2=DOT_PRODUCT(dgrad,Wy) |
173 | 1 | pfleura2 | |
174 | 1 | pfleura2 | if (debug) THEN |
175 | 1 | pfleura2 | WRITE(*,*) "(sT)y,(yT)Wy ",den,den2 |
176 | 1 | pfleura2 | END if |
177 | 1 | pfleura2 | |
178 | 1 | pfleura2 | |
179 | 1 | pfleura2 | !We add everything |
180 | 1 | pfleura2 | Hup=Num/den-Num2/den2 |
181 | 1 | pfleura2 | Hess=HessOld+Hup |
182 | 1 | pfleura2 | |
183 | 1 | pfleura2 | if (debug) THEN |
184 | 1 | pfleura2 | WRITE(*,*) 'DBG Hinvup_DFP Hup' |
185 | 1 | pfleura2 | Idx=min(12,Nfree) |
186 | 1 | pfleura2 | DO I=1,NFree |
187 | 1 | pfleura2 | WRITE(*,'(12(1X,F8.4))') Hup(I,1:Idx) |
188 | 1 | pfleura2 | END DO |
189 | 1 | pfleura2 | WRITE(*,*) 'DBG Hinvup_DFP Hess new' |
190 | 1 | pfleura2 | Idx=min(12,Nfree) |
191 | 1 | pfleura2 | DO I=1,NFree |
192 | 1 | pfleura2 | WRITE(*,'(12(1X,F8.4))') Hess(I,1:Idx) |
193 | 1 | pfleura2 | END DO |
194 | 1 | pfleura2 | END IF |
195 | 1 | pfleura2 | |
196 | 1 | pfleura2 | ELSE |
197 | 1 | pfleura2 | Hess=HessOld |
198 | 1 | pfleura2 | END IF |
199 | 1 | pfleura2 | |
200 | 1 | pfleura2 | DEALLOCATE(HessOld, Wy, Hup, Num, Num2 ) |
201 | 1 | pfleura2 | |
202 | 1 | pfleura2 | if (debug) WRITE(*,*) "=========================== Hinvup_DFP Over ====================" |
203 | 1 | pfleura2 | |
204 | 1 | pfleura2 | END SUBROUTINE Hinvup_DFP |