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