Statistiques
| Révision :

root / src / Hupdate_MS.f90 @ 12

Historique | Voir | Annoter | Télécharger (4,23 ko)

1 12 pfleura2
SUBROUTINE Hupdate_MS(Nfree,ds,dgrad,Hess)
2 12 pfleura2
3 1 pfleura2
! Hessian - or its inverse - update
4 1 pfleura2
! We use the Murtagh-Sargent update, also called Symmetric Rank 1
5 1 pfleura2
!
6 1 pfleura2
! input:
7 1 pfleura2
!  - nfree: number of coordinates
8 1 pfleura2
!  - ds: GeomOld-Geom (expressed in the optimization coord.)
9 1 pfleura2
!  - dgrad: GradOld-Grad  (expressed in the optimization coord.)
10 1 pfleura2
!
11 1 pfleura2
!  input/output:
12 1 pfleura2
!  - Hess: old Hessian in input, replaced by updated Hessian on output
13 1 pfleura2
!             (or inverse Hessian)
14 1 pfleura2
!
15 1 pfleura2
!
16 1 pfleura2
!!!!!!!!!!!!!!!!!!!!!!
17 1 pfleura2
! The update formula for the Hessian is:
18 1 pfleura2
! (following Fletcher, approximate Hessian is denoted by B)
19 1 pfleura2
! B(k+1)=B+[(y-Bs){(y-Bs)T}]/[{(y-Bs)T}s]
20 1 pfleura2
!
21 1 pfleura2
! it is self dual, ie for the inverse Hessian, denote by H:
22 1 pfleura2
! H(k+1)=H+[(s-Hy){(s-Hy)T}]/[{(s-Hy)T}y]
23 1 pfleura2
!
24 1 pfleura2
!
25 1 pfleura2
! We decompose the computation in two steps:
26 1 pfleura2
! first we compute the update, and then we add it to Hess to get the
27 1 pfleura2
! new Hessian.
28 1 pfleura2
! This allows us to introduce a scaling factor if needed later.
29 1 pfleura2
!
30 1 pfleura2
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
31 1 pfleura2
32 1 pfleura2
33 12 pfleura2
!----------------------------------------------------------------------
34 12 pfleura2
!  Copyright 2003-2014 Ecole Normale Supérieure de Lyon,
35 12 pfleura2
!  Centre National de la Recherche Scientifique,
36 12 pfleura2
!  Université Claude Bernard Lyon 1. All rights reserved.
37 12 pfleura2
!
38 12 pfleura2
!  This work is registered with the Agency for the Protection of Programs
39 12 pfleura2
!  as IDDN.FR.001.100009.000.S.P.2014.000.30625
40 12 pfleura2
!
41 12 pfleura2
!  Authors: P. Fleurat-Lessard, P. Dayal
42 12 pfleura2
!  Contact: optnpath@gmail.com
43 12 pfleura2
!
44 12 pfleura2
! This file is part of "Opt'n Path".
45 12 pfleura2
!
46 12 pfleura2
!  "Opt'n Path" is free software: you can redistribute it and/or modify
47 12 pfleura2
!  it under the terms of the GNU Affero General Public License as
48 12 pfleura2
!  published by the Free Software Foundation, either version 3 of the License,
49 12 pfleura2
!  or (at your option) any later version.
50 12 pfleura2
!
51 12 pfleura2
!  "Opt'n Path" is distributed in the hope that it will be useful,
52 12 pfleura2
!  but WITHOUT ANY WARRANTY; without even the implied warranty of
53 12 pfleura2
!
54 12 pfleura2
!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
55 12 pfleura2
!  GNU Affero General Public License for more details.
56 12 pfleura2
!
57 12 pfleura2
!  You should have received a copy of the GNU Affero General Public License
58 12 pfleura2
!  along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>.
59 12 pfleura2
!
60 12 pfleura2
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr,
61 12 pfleura2
! for commercial licensing opportunities.
62 12 pfleura2
!----------------------------------------------------------------------
63 12 pfleura2
64 1 pfleura2
  IMPLICIT NONE
65 1 pfleura2
66 1 pfleura2
  INTEGER, PARAMETER :: KINT=KIND(1)
67 1 pfleura2
  INTEGER, PARAMETER :: KREAL=KIND(1.D0)
68 1 pfleura2
69 1 pfleura2
  INTEGER(KINT), INTENT(IN) :: NFREE
70 1 pfleura2
  REAL(KREAL), INTENT(INOUT) :: Hess(Nfree,Nfree)
71 1 pfleura2
  REAL(KREAL), INTENT(IN) :: ds(Nfree)
72 1 pfleura2
  REAL(KREAL), INTENT(IN) :: dGrad(Nfree)
73 1 pfleura2
74 1 pfleura2
75 2 pfleura2
  REAL(KREAL), ALLOCATABLE  :: HessOld(:, :)
76 1 pfleura2
  REAL(KREAL), ALLOCATABLE  :: Hup(:,:)
77 1 pfleura2
  REAL(KREAL), ALLOCATABLE  :: HupInt(:)
78 1 pfleura2
  REAL(KREAL), ALLOCATABLE  :: HS(:), Num(:,:)
79 2 pfleura2
  REAL(KREAL) :: Den, NormDS
80 1 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
93 1 pfleura2
94 1 pfleura2
95 1 pfleura2
  debug=valid("Hupdate").OR.Valid("Hupdate_MS")
96 1 pfleura2
97 1 pfleura2
  if (debug) Call Header("Entering Hupdate_MS")
98 1 pfleura2
99 1 pfleura2
  ALLOCATE(HessOld(nfree,nfree),  &
100 1 pfleura2
       Hup(nfree,nfree), &
101 1 pfleura2
       HupInt(nfree), HS(nfree),  &
102 1 pfleura2
       Num(nfree,nfree))
103 1 pfleura2
104 1 pfleura2
  HessOld=Hess
105 1 pfleura2
106 1 pfleura2
  if (debug) THEN
107 1 pfleura2
     WRITE(*,*) 'DBG Hupdate HessOld'
108 1 pfleura2
     Idx=min(12,Nfree)
109 1 pfleura2
     DO I=1,NFree
110 1 pfleura2
        WRITE(*,'(12(1X,F8.4))') HessOld(I,1:Idx)
111 1 pfleura2
     END DO
112 1 pfleura2
     WRITE(*,*) 'DBG Hupdate ds'
113 1 pfleura2
     WRITE(*,'(12(1X,F8.4))') ds(1:Idx)
114 1 pfleura2
     WRITE(*,*) 'DBG Hupdate dGrad'
115 1 pfleura2
     WRITE(*,'(12(1X,F8.4))') dGrad(1:Idx)
116 1 pfleura2
  END IF
117 1 pfleura2
118 1 pfleura2
  NormDS=sqrt(DOT_PRODUCT(DS,DS))
119 1 pfleura2
  if (NormDS>=1e-6) THEN
120 1 pfleura2
121 1 pfleura2
     HS=0
122 1 pfleura2
     DO I=1,nfree
123 1 pfleura2
        DO J=1,nfree
124 1 pfleura2
           HS(I)=HS(I)+HessOld(i,j)*Ds(j)
125 1 pfleura2
        END DO
126 1 pfleura2
     END DO
127 1 pfleura2
128 1 pfleura2
     if (debug) WRITE(*,*) "HS:",HS
129 1 pfleura2
     HupInt=DGrad-HS
130 1 pfleura2
131 1 pfleura2
     ! WE calculate the MS update
132 1 pfleura2
     Num=0.
133 1 pfleura2
     Den=0.
134 1 pfleura2
     DO I=1,nfree
135 1 pfleura2
        DO J=1,nfree
136 1 pfleura2
           Num(I,J)=HUpInt(i)*HupInt(j)
137 1 pfleura2
        END DO
138 1 pfleura2
        Den=Den+Hupint(i)*Ds(i)
139 1 pfleura2
     END DO
140 1 pfleura2
     if (debug) WRITE(*,*) "Den:",Den,1./Den
141 1 pfleura2
142 1 pfleura2
     Hup=Num/Den
143 1 pfleura2
144 1 pfleura2
     Hess=HessOld+Hup
145 1 pfleura2
  ELSE
146 1 pfleura2
     Hess=HessOld
147 1 pfleura2
  END IF
148 1 pfleura2
  !     WRITE(*,*) 'New hess'
149 1 pfleura2
  !     DO I=1,nfree
150 1 pfleura2
  !        WRITE(*,fmt) Hess(i,:)
151 1 pfleura2
  !     END DO
152 1 pfleura2
153 1 pfleura2
  DEALLOCATE(HessOld,Hup, HupInt, HS, Num )
154 1 pfleura2
155 1 pfleura2
  if (debug) Call Header("Hupdate_MS Over")
156 1 pfleura2
157 1 pfleura2
END SUBROUTINE Hupdate_MS