Statistiques
| Révision :

root / src / Hupdate_Bofill.f90 @ 12

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

1 1 pfleura2
! This program reads an old hessian, old forces, new forces
2 1 pfleura2
! and print the updated hessian and old coordinates.
3 1 pfleura2
! Instead of the Murtagh-Sargent update, we use the Bofill update
4 1 pfleura2
! which is an average of Murtagh-Sargent and Powell-symmetryc-Broyden
5 1 pfleura2
! cf JCP 1999, vol 111 (19), p 8773
6 1 pfleura2
7 1 pfleura2
! This subroutine slightly differs from previous ones as
8 1 pfleura2
! it is written with the Gradient as input and not the forces
9 1 pfleura2
! Using the relation F=-Grad, one can easily go from older
10 1 pfleura2
! versions to this one (which more closely mimics the JCP article)
11 1 pfleura2
12 1 pfleura2
SUBROUTINE Hupdate_Bofill(Nfree,ds,dgrad,Hess)
13 1 pfleura2
14 12 pfleura2
!----------------------------------------------------------------------
15 12 pfleura2
!  Copyright 2003-2014 Ecole Normale Supérieure de Lyon,
16 12 pfleura2
!  Centre National de la Recherche Scientifique,
17 12 pfleura2
!  Université Claude Bernard Lyon 1. All rights reserved.
18 12 pfleura2
!
19 12 pfleura2
!  This work is registered with the Agency for the Protection of Programs
20 12 pfleura2
!  as IDDN.FR.001.100009.000.S.P.2014.000.30625
21 12 pfleura2
!
22 12 pfleura2
!  Authors: P. Fleurat-Lessard, P. Dayal
23 12 pfleura2
!  Contact: optnpath@gmail.com
24 12 pfleura2
!
25 12 pfleura2
! This file is part of "Opt'n Path".
26 12 pfleura2
!
27 12 pfleura2
!  "Opt'n Path" is free software: you can redistribute it and/or modify
28 12 pfleura2
!  it under the terms of the GNU Affero General Public License as
29 12 pfleura2
!  published by the Free Software Foundation, either version 3 of the License,
30 12 pfleura2
!  or (at your option) any later version.
31 12 pfleura2
!
32 12 pfleura2
!  "Opt'n Path" is distributed in the hope that it will be useful,
33 12 pfleura2
!  but WITHOUT ANY WARRANTY; without even the implied warranty of
34 12 pfleura2
!
35 12 pfleura2
!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
36 12 pfleura2
!  GNU Affero General Public License for more details.
37 12 pfleura2
!
38 12 pfleura2
!  You should have received a copy of the GNU Affero General Public License
39 12 pfleura2
!  along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>.
40 12 pfleura2
!
41 12 pfleura2
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr,
42 12 pfleura2
! for commercial licensing opportunities.
43 12 pfleura2
!----------------------------------------------------------------------
44 12 pfleura2
45 1 pfleura2
  IMPLICIT NONE
46 1 pfleura2
47 1 pfleura2
  INTEGER, PARAMETER :: KINT=KIND(1)
48 1 pfleura2
  INTEGER, PARAMETER :: KREAL=KIND(1.D0)
49 1 pfleura2
50 1 pfleura2
  INTEGER(KINT), INTENT(IN) :: NFREE
51 1 pfleura2
  REAL(KREAL), INTENT(INOUT) :: Hess(Nfree,Nfree)
52 1 pfleura2
  REAL(KREAL), INTENT(IN) :: ds(Nfree)
53 1 pfleura2
  REAL(KREAL), INTENT(IN) :: dGrad(Nfree)
54 1 pfleura2
55 1 pfleura2
56 1 pfleura2
  REAL(KREAL), ALLOCATABLE  :: HessOld(:,:) , HupMS(:,:)
57 1 pfleura2
  REAL(KREAL), ALLOCATABLE  :: HupPSB(:,:), Hup(:,:)
58 1 pfleura2
  REAL(KREAL), ALLOCATABLE  :: HupInt(:)
59 1 pfleura2
  REAL(KREAL), ALLOCATABLE  :: HS(:), Num(:,:), Num2(:,:)
60 1 pfleura2
  REAL(KREAL) :: Phi, Num3, Den, Den2,NormDS
61 1 pfleura2
  INTEGER(KINT) :: i,j, idx
62 1 pfleura2
63 1 pfleura2
  LOGICAL debug,valid
64 1 pfleura2
  EXTERNAL valid
65 1 pfleura2
66 1 pfleura2
  debug=valid("Hupdate").OR.Valid("Hupdate_bofill")
67 1 pfleura2
68 1 pfleura2
  if (debug) WRITE(*,*) "=========================== Entering Hupdate_Bofill ===================="
69 1 pfleura2
70 1 pfleura2
  ALLOCATE(HessOld(nfree,nfree), HupMS(nfree,nfree),  &
71 1 pfleura2
       HupPSB(nfree,nfree), Hup(nfree,nfree), &
72 1 pfleura2
       HupInt(nfree), HS(nfree),  &
73 1 pfleura2
       Num(nfree,nfree), Num2(nfree,nfree) )
74 1 pfleura2
75 1 pfleura2
  HessOld=Hess
76 1 pfleura2
77 1 pfleura2
  if (debug) THEN
78 1 pfleura2
     WRITE(*,*) 'DBG Hupdate HessOld'
79 1 pfleura2
     Idx=min(12,Nfree)
80 1 pfleura2
     DO I=1,NFree
81 1 pfleura2
        WRITE(*,'(12(1X,F8.4))') HessOld(I,1:Idx)
82 1 pfleura2
     END DO
83 1 pfleura2
     WRITE(*,*) 'DBG Hupdate ds'
84 1 pfleura2
     WRITE(*,'(12(1X,F8.4))') ds(1:Idx)
85 1 pfleura2
     WRITE(*,*) 'DBG Hupdate dGrad'
86 1 pfleura2
     WRITE(*,'(12(1X,F8.4))') dGrad(1:Idx)
87 1 pfleura2
  END IF
88 1 pfleura2
89 1 pfleura2
  NormDS=sqrt(DOT_PRODUCT(DS,DS))
90 1 pfleura2
  if (NormDS>=1e-6) THEN
91 1 pfleura2
92 1 pfleura2
     HS=0
93 1 pfleura2
     DO I=1,nfree
94 1 pfleura2
        DO J=1,nfree
95 1 pfleura2
           HS(I)=HS(I)+HessOld(i,j)*Ds(j)
96 1 pfleura2
        END DO
97 1 pfleura2
     END DO
98 1 pfleura2
     HupInt=DGrad-HS
99 1 pfleura2
100 1 pfleura2
     ! WE calculate the MS update
101 1 pfleura2
     Num=0.
102 1 pfleura2
     Den=0.
103 1 pfleura2
     DO I=1,nfree
104 1 pfleura2
        DO J=1,nfree
105 1 pfleura2
           Num(I,J)=HUpInt(i)*HupInt(j)
106 1 pfleura2
        END DO
107 1 pfleura2
        Den=Den+Hupint(i)*Ds(i)
108 1 pfleura2
     END DO
109 1 pfleura2
110 1 pfleura2
     HupMS=Num/Den
111 1 pfleura2
112 1 pfleura2
     ! Now we calculate the PSB update
113 1 pfleura2
     Num=0
114 1 pfleura2
     Num2=0
115 1 pfleura2
     Num3=0
116 1 pfleura2
     Den=0
117 1 pfleura2
     DO I=1,nfree
118 1 pfleura2
        DO J=1,nfree
119 1 pfleura2
           Num(I,J)=HUpInt(i)*Ds(j)+Ds(i)*HupInt(j)
120 1 pfleura2
           Num2(I,J)=DS(i)*DS(j)
121 1 pfleura2
        END DO
122 1 pfleura2
        Den=Den+Ds(i)*Ds(i)
123 1 pfleura2
        Num3=Num3+DS(i)*HupInt(i)
124 1 pfleura2
     END DO
125 1 pfleura2
     HupPSB=Num/Den-Num3*Num2/(Den)**2
126 1 pfleura2
127 1 pfleura2
     ! We calculate Phi
128 1 pfleura2
     Den2=0.
129 1 pfleura2
     Do I=1,nfree
130 1 pfleura2
        Den2=Den2+HupInt(i)*HupInt(i)
131 1 pfleura2
     END DO
132 1 pfleura2
     Phi=Num3**2/(Den*Den2)
133 1 pfleura2
134 1 pfleura2
     Hup=Phi*HupMS+(1-Phi)*HupPSB
135 1 pfleura2
136 1 pfleura2
     Hess=HessOld+Hup
137 1 pfleura2
  ELSE
138 1 pfleura2
     Hess=HessOld
139 1 pfleura2
  END IF
140 1 pfleura2
  !     WRITE(*,*) 'New hess'
141 1 pfleura2
  !     DO I=1,nfree
142 1 pfleura2
  !        WRITE(*,fmt) Hess(i,:)
143 1 pfleura2
  !     END DO
144 1 pfleura2
145 1 pfleura2
  DEALLOCATE(HessOld, HupMS, HupPSB, Hup, &
146 1 pfleura2
       HupInt, HS, Num, Num2 )
147 1 pfleura2
148 1 pfleura2
  if (debug) WRITE(*,*) "=========================== Hupdate Over ===================="
149 1 pfleura2
150 1 pfleura2
END SUBROUTINE Hupdate_Bofill