Statistiques
| Révision :

root / src / Step_GEDIIS_All.f90 @ 12

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

1 1 pfleura2
     ! Geom = input parameter vector (Geometry), Grad = input gradient vector, HEAT is Energy(Geom)
2 1 pfleura2
      SUBROUTINE Step_GEDIIS_All(NGeomF,IGeom,Step,Geom,Grad,HEAT,Hess,NCoord,allocation_flag,Tangent)
3 12 pfleura2
4 12 pfleura2
!----------------------------------------------------------------------
5 12 pfleura2
! This routine was adapted from the public domain mopac6 diis.f
6 12 pfleura2
!  source file (c) 2009, Stewart Computational Chemistry.
7 12 pfleura2
!  <http://www.openmopac.net/Downloads/Downloads.html>
8 12 pfleura2
!
9 12 pfleura2
!----------------------------------------------------------------------
10 12 pfleura2
!  Copyright 2003-2014 Ecole Normale Supérieure de Lyon,
11 12 pfleura2
!  Centre National de la Recherche Scientifique,
12 12 pfleura2
!  Université Claude Bernard Lyon 1. All rights reserved.
13 12 pfleura2
!
14 12 pfleura2
!  This work is registered with the Agency for the Protection of Programs
15 12 pfleura2
!  as IDDN.FR.001.100009.000.S.P.2014.000.30625
16 12 pfleura2
!
17 12 pfleura2
!  Authors: P. Fleurat-Lessard, P. Dayal
18 12 pfleura2
!  Contact: optnpath@gmail.com
19 12 pfleura2
!
20 12 pfleura2
! This file is part of "Opt'n Path".
21 12 pfleura2
!
22 12 pfleura2
!  "Opt'n Path" is free software: you can redistribute it and/or modify
23 12 pfleura2
!  it under the terms of the GNU Affero General Public License as
24 12 pfleura2
!  published by the Free Software Foundation, either version 3 of the License,
25 12 pfleura2
!  or (at your option) any later version.
26 12 pfleura2
!
27 12 pfleura2
!  "Opt'n Path" is distributed in the hope that it will be useful,
28 12 pfleura2
!  but WITHOUT ANY WARRANTY; without even the implied warranty of
29 12 pfleura2
!
30 12 pfleura2
!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
31 12 pfleura2
!  GNU Affero General Public License for more details.
32 12 pfleura2
!
33 12 pfleura2
!  You should have received a copy of the GNU Affero General Public License
34 12 pfleura2
!  along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>.
35 12 pfleura2
!
36 12 pfleura2
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr,
37 12 pfleura2
! for commercial licensing opportunities.
38 12 pfleura2
!----------------------------------------------------------------------
39 12 pfleura2
40 1 pfleura2
    use Io_module
41 12 pfleura2
    use Path_module, only :  Vfree
42 1 pfleura2
      IMPLICIT NONE
43 1 pfleura2
44 1 pfleura2
      INTEGER(KINT) :: NGeomF,IGeom
45 1 pfleura2
      INTEGER(KINT), INTENT(IN) :: NCoord
46 2 pfleura2
      REAL(KREAL) :: Geom(NCoord), Grad(NCoord), Hess(NCoord*NCoord), Step(NCoord)
47 1 pfleura2
    REAL(KREAL) :: HEAT ! HEAT= Energy
48 1 pfleura2
    LOGICAL :: allocation_flag
49 1 pfleura2
    REAL(KREAL), INTENT(INOUT) :: Tangent(Ncoord)
50 1 pfleura2
51 1 pfleura2
      ! MRESET = maximum number of iterations.
52 1 pfleura2
      INTEGER(KINT), PARAMETER :: MRESET=15, M2=(MRESET+1)*(MRESET+1) !M2 = 256
53 1 pfleura2
      REAL(KREAL), ALLOCATABLE, SAVE :: GeomSet(:,:), GradSet(:,:) ! NGeomF,MRESET*NCoord
54 1 pfleura2
      REAL(KREAL), ALLOCATABLE, SAVE :: GSAVE(:,:) !NGeomF,NCoord
55 1 pfleura2
    REAL(KREAL), ALLOCATABLE, SAVE :: ESET(:,:)
56 1 pfleura2
    REAL(KREAL) :: ESET_tmp(MRESET), B(M2), BS(M2), BST(M2), B_tmp(M2) ! M2=256
57 1 pfleura2
      LOGICAL :: DEBUG, PRINT, ci_lt_zero
58 1 pfleura2
      INTEGER(KINT), ALLOCATABLE, SAVE :: MSET(:) ! mth Iteration
59 1 pfleura2
    LOGICAL, ALLOCATABLE, SAVE :: FRST(:)
60 1 pfleura2
    REAL(KREAL) :: ci(MRESET), ci_tmp(MRESET) ! MRESET = maximum number of iterations.
61 2 pfleura2
      INTEGER(KINT) :: NGEDIIS, MPLUS, INV, ITERA, MM, cis_zero
62 2 pfleura2
      INTEGER(KINT) :: I, J, K, JJ, JNV, II, IONE, IJ, IX, JX, KX
63 2 pfleura2
    INTEGER(KINT) :: current_size_B_mat, MyPointer, Isch, NFree, Idx
64 2 pfleura2
      REAL(KREAL) :: XMax, XNorm, DET, THRES, tmp, ER_star, ER_star_tmp, Norm
65 1 pfleura2
    REAL(KREAL), PARAMETER :: eps=1e-12
66 1 pfleura2
      REAL(KREAL), PARAMETER :: crit=1e-8
67 1 pfleura2
    REAL(KREAL), ALLOCATABLE :: Tanf(:) ! NCoord
68 1 pfleura2
    REAL(KREAL), ALLOCATABLE :: HFree(:) ! NFree*NFree
69 1 pfleura2
    REAL(KREAL), ALLOCATABLE :: Htmp(:,:) ! NCoord,NFree
70 1 pfleura2
    REAL(KREAL), ALLOCATABLE :: Grad_free(:), Step_free(:) ! NFree
71 1 pfleura2
    REAL(KREAL), ALLOCATABLE :: Geom_free(:), Geom_new_free(:) ! NFree
72 1 pfleura2
    REAL(KREAL), ALLOCATABLE, SAVE :: GeomSet_free(:,:), GradSet_free(:,:)
73 1 pfleura2
74 1 pfleura2
      DEBUG=.TRUE.
75 1 pfleura2
      PRINT=.FALSE.
76 1 pfleura2
77 1 pfleura2
      IF (PRINT)  WRITE(*,'(/,''      BEGIN Step_GEDIIS_ALL   '')')
78 1 pfleura2
79 1 pfleura2
      ! Initialization
80 1 pfleura2
      IF (allocation_flag) THEN
81 1 pfleura2
      ! allocation_flag will be set to False in SPACE_GEDIIS, so no need to modify it here
82 1 pfleura2
         IF (ALLOCATED(GeomSet)) THEN
83 1 pfleura2
            IF (PRINT)  WRITE(*,'(/,''    In allocation_flag, GEDIIS_ALL  Dealloc  '')')
84 1 pfleura2
            DEALLOCATE(GeomSet,GradSet,GSave,GeomSet_free,GradSet_free)
85 1 pfleura2
            RETURN
86 1 pfleura2
         ELSE
87 1 pfleura2
            IF (PRINT)  WRITE(*,'(/,''     In allocation_flag,  GEDIIS_ALL  Alloc  '')')
88 1 pfleura2
            ALLOCATE(GeomSet(NGeomF,MRESET*NCoord),GradSet(NGeomF,MRESET*NCoord),GSAVE(NGeomF,NCoord))
89 1 pfleura2
      ALLOCATE(GeomSet_free(NGeomF,MRESET*NCoord),GradSet_free(NGeomF,MRESET*NCoord))
90 1 pfleura2
      ALLOCATE(MSET(NGeomF),FRST(NGeomF),ESET(NGeomF,MRESET))
91 1 pfleura2
      DO I=1,NGeomF
92 1 pfleura2
         FRST(I) = .TRUE.
93 1 pfleura2
         GeomSet(I,:) = 0.d0
94 1 pfleura2
         GradSet(I,:) = 0.d0
95 1 pfleura2
         GSAVE(I,:)=0.d0
96 1 pfleura2
         GeomSet_free(I,:) = 0.d0
97 1 pfleura2
         GradSet_free(I,:) = 0.d0
98 1 pfleura2
      END DO
99 1 pfleura2
      MSET(:)=0
100 1 pfleura2
         END IF
101 1 pfleura2
     allocation_flag = .FALSE.
102 1 pfleura2
      END IF ! IF (allocation_flag) THEN
103 1 pfleura2
104 1 pfleura2
    ! ADDED FROM HERE:
105 1 pfleura2
    Call FreeMv(NCoord,Vfree) ! VFree(Ncoord,Ncoord), as of now, an Identity matrix.
106 1 pfleura2
      ! we orthogonalize Vfree to the tangent vector of this geom only if Tangent/=0.d0
107 1 pfleura2
      Norm=sqrt(dot_product(Tangent,Tangent))
108 1 pfleura2
      IF (Norm.GT.eps) THEN
109 1 pfleura2
         ALLOCATE(Tanf(NCoord))
110 1 pfleura2
111 1 pfleura2
         ! We normalize Tangent
112 1 pfleura2
         Tangent=Tangent/Norm
113 1 pfleura2
114 1 pfleura2
         ! We convert Tangent into Vfree only displacements. This is useless for now (2007.Apr.23)
115 1 pfleura2
     ! as Vfree=Id matrix but it will be usefull as soon as we introduce constraints.
116 1 pfleura2
         DO I=1,NCoord
117 1 pfleura2
            Tanf(I)=dot_product(reshape(Vfree(:,I),(/NCoord/)),Tangent)
118 1 pfleura2
         END DO
119 1 pfleura2
         Tangent=0.d0
120 1 pfleura2
         DO I=1,NCoord
121 1 pfleura2
            Tangent=Tangent+Tanf(I)*Vfree(:,I)
122 1 pfleura2
         END DO
123 1 pfleura2
         ! first we subtract Tangent from vfree
124 1 pfleura2
         DO I=1,NCoord
125 1 pfleura2
            Norm=dot_product(reshape(vfree(:,I),(/NCoord/)),Tangent)
126 1 pfleura2
            Vfree(:,I)=Vfree(:,I)-Norm*Tangent
127 1 pfleura2
     END DO
128 1 pfleura2
129 8 pfleura2
         Idx=0
130 1 pfleura2
         ! Schmidt orthogonalization of the Vfree vectors
131 1 pfleura2
         DO I=1,NCoord
132 1 pfleura2
            ! We subtract the first vectors, we do it twice as the Schmidt procedure is not numerically stable.
133 1 pfleura2
            DO Isch=1,2
134 1 pfleura2
               DO J=1,Idx
135 1 pfleura2
                  Norm=dot_product(reshape(Vfree(:,I),(/NCoord/)),reshape(Vfree(:,J),(/NCoord/)))
136 1 pfleura2
                  Vfree(:,I)=Vfree(:,I)-Norm*Vfree(:,J)
137 1 pfleura2
               END DO
138 1 pfleura2
            END DO
139 1 pfleura2
            Norm=dot_product(reshape(Vfree(:,I),(/NCoord/)),reshape(Vfree(:,I),(/NCoord/)))
140 1 pfleura2
            IF (Norm.GE.crit) THEN
141 1 pfleura2
               Idx=Idx+1
142 1 pfleura2
               Vfree(:,Idx)=Vfree(:,I)/sqrt(Norm)
143 1 pfleura2
            END IF
144 1 pfleura2
         END DO
145 1 pfleura2
146 1 pfleura2
         IF (Idx/= NCoord-1) THEN
147 1 pfleura2
            WRITE(*,*) "Pb in orthogonalizing Vfree to tangent for geom",IGeom
148 1 pfleura2
            WRITE(IOOut,*) "Pb in orthogonalizing Vfree to tangent for geom",IGeom
149 1 pfleura2
            STOP
150 1 pfleura2
         END IF
151 1 pfleura2
152 1 pfleura2
         DEALLOCATE(Tanf)
153 1 pfleura2
         NFree=Idx
154 1 pfleura2
      ELSE ! Tangent =0, matches IF (Norm.GT.eps) THEN
155 1 pfleura2
         if (debug) WRITE(*,*) "Tangent=0, using full displacement"
156 1 pfleura2
         NFree=NCoord
157 1 pfleura2
      END IF !IF (Norm.GT.eps) THEN
158 1 pfleura2
159 1 pfleura2
      if (debug) WRITE(*,*) 'DBG Step_GEDIIS_All, IGeom, NFree=', IGeom, NFree
160 1 pfleura2
161 1 pfleura2
      ! We now calculate the new step
162 1 pfleura2
      ! we project the hessian onto the free vectors
163 1 pfleura2
      ALLOCATE(HFree(NFree*NFree),Htmp(NCoord,NFree),Grad_free(NFree))
164 1 pfleura2
    ALLOCATE(Geom_free(NFree),Step_free(NFree),Geom_new_free(NFree))
165 1 pfleura2
      DO J=1,NFree
166 1 pfleura2
        DO I=1,NCoord
167 1 pfleura2
           Htmp(I,J)=0.d0
168 1 pfleura2
           DO K=1,NCoord
169 1 pfleura2
              Htmp(I,J)=Htmp(I,J)+Hess(((I-1)*NCoord)+K)*Vfree(K,J)
170 1 pfleura2
           END DO
171 1 pfleura2
        END DO
172 1 pfleura2
      END DO
173 1 pfleura2
      DO J=1,NFree
174 1 pfleura2
        DO I=1,NFree
175 1 pfleura2
       HFree(I+((J-1)*NFree))=0.d0
176 1 pfleura2
           DO K=1,NCoord
177 1 pfleura2
        HFree(I+((J-1)*NFree))=HFree(I+((J-1)*NFree))+Vfree(K,I)*Htmp(K,J)
178 1 pfleura2
           END DO
179 1 pfleura2
        END DO
180 1 pfleura2
      END DO
181 1 pfleura2
182 1 pfleura2
      DO I=1,NFree
183 1 pfleura2
         Grad_free(I)=dot_product(reshape(Vfree(:,I),(/NCoord/)),Grad)
184 1 pfleura2
     Geom_free(I)=dot_product(reshape(Vfree(:,I),(/NCoord/)),Geom)
185 1 pfleura2
      END DO
186 1 pfleura2
      !ADDED ENDS HERE.***********************************************
187 1 pfleura2
188 1 pfleura2
      ! SPACE_GEDIIS SIMPLY LOADS THE CURRENT VALUES OF Geom AND Grad INTO THE ARRAYS GeomSet
189 1 pfleura2
      ! AND GradSet, MSET is set to zero and then 1 in SPACE_GEDIIS_All at first iteration.
190 1 pfleura2
      CALL SPACE_GEDIIS_All(NGeomF,IGeom,MRESET,MSET,Geom,Grad,HEAT,NCoord,GeomSet,GradSet,ESET,FRST)
191 1 pfleura2
192 1 pfleura2
      IF (PRINT)  WRITE(*,'(/,''       GEDIIS after SPACE_GEDIIS_ALL   '')')
193 1 pfleura2
194 1 pfleura2
    DO J=1,MSet(IGeom)
195 1 pfleura2
         DO K=1,NFree
196 1 pfleura2
            GradSet_free(IGeom,((J-1)*NFree)+K)=dot_product(reshape(Vfree(:,K),(/NCoord/)),&
197 1 pfleura2
               GradSet(IGeom,((J-1)*NCoord)+1:((J-1)*NCoord)+NCoord))
198 1 pfleura2
        GeomSet_free(IGeom,((J-1)*NFree)+K)=dot_product(reshape(Vfree(:,K),(/NCoord/)),&
199 1 pfleura2
                GeomSet(IGeom,((J-1)*NCoord)+1:((J-1)*NCoord)+NCoord))
200 1 pfleura2
         END DO
201 1 pfleura2
    END DO
202 1 pfleura2
203 1 pfleura2
      ! INITIALIZE SOME VARIABLES AND CONSTANTS:
204 1 pfleura2
      NGEDIIS = MSET(IGeom) !MSET=mth iteration
205 1 pfleura2
      MPLUS = MSET(IGeom) + 1
206 1 pfleura2
      MM = MPLUS * MPLUS
207 1 pfleura2
208 1 pfleura2
      ! CONSTRUCT THE GEDIIS MATRIX:
209 1 pfleura2
      ! B_ij calculations from <B_ij=(g_i-g_j)(R_i-R_j)>
210 1 pfleura2
      JJ=0
211 1 pfleura2
      INV=-NFree
212 1 pfleura2
      DO I=1,MSET(IGeom)
213 1 pfleura2
         INV=INV+NFree
214 1 pfleura2
         JNV=-NFree
215 1 pfleura2
         DO J=1,MSET(IGeom)
216 1 pfleura2
            JNV=JNV+NFree
217 1 pfleura2
            JJ = JJ + 1
218 1 pfleura2
            B(JJ)=0.D0
219 1 pfleura2
      DO K=1, NFree
220 1 pfleura2
         B(JJ) = B(JJ) + (((GradSet_free(IGeom,INV+K)-GradSet_free(IGeom,JNV+K))* &
221 1 pfleura2
                 (GeomSet_free(IGeom,INV+K)-GeomSet_free(IGeom,JNV+K)))/2.D0)
222 1 pfleura2
      END DO
223 1 pfleura2
         END DO
224 1 pfleura2
      END DO
225 1 pfleura2
226 1 pfleura2
     ! The following shifting is required to correct indices of B_ij elements in the GEDIIS matrix.
227 1 pfleura2
   ! The correction is needed because the last coloumn of the matrix contains all 1 and one zero.
228 1 pfleura2
      DO I=MSET(IGeom)-1,1,-1
229 1 pfleura2
         DO J=MSET(IGeom),1,-1
230 1 pfleura2
            B(I*MSET(IGeom)+J+I) = B(I*MSET(IGeom)+J)
231 1 pfleura2
         END DO
232 1 pfleura2
    END DO
233 1 pfleura2
234 1 pfleura2
      ! For the last row and last column of GEDIIS matrix:
235 1 pfleura2
      DO I=1,MPLUS
236 1 pfleura2
         B(MPLUS*I) = 1.D0
237 1 pfleura2
         B(MPLUS*MSET(IGeom)+I) = 1.D0
238 1 pfleura2
      END DO
239 1 pfleura2
      B(MM) = 0.D0
240 1 pfleura2
241 1 pfleura2
    DO I=1, MPLUS
242 1 pfleura2
       !WRITE(*,'(10(1X,F20.4))') B((I-1)*MPLUS+1:I*(MPLUS))
243 1 pfleura2
    END DO
244 1 pfleura2
245 1 pfleura2
      ! ELIMINATE ERROR VECTORS WITH THE LARGEST NORM:
246 1 pfleura2
   80 CONTINUE
247 1 pfleura2
      DO I=1,MM !MM = (MSET(IGeom)+1) * (MSET(IGeom)+1)
248 1 pfleura2
         BS(I) = B(I) !just a copy of the original B (GEDIIS) matrix
249 1 pfleura2
      END DO
250 1 pfleura2
251 1 pfleura2
      IF (NGEDIIS .NE. MSET(IGeom)) THEN
252 1 pfleura2
        DO II=1,MSET(IGeom)-NGEDIIS
253 1 pfleura2
           XMAX = -1.D10
254 1 pfleura2
           ITERA = 0
255 1 pfleura2
           DO I=1,MSET(IGeom)
256 1 pfleura2
              XNORM = 0.D0
257 1 pfleura2
              INV = (I-1) * MPLUS
258 1 pfleura2
              DO J=1,MSET(IGeom)
259 1 pfleura2
                 XNORM = XNORM + ABS(B(INV + J))
260 1 pfleura2
              END DO
261 1 pfleura2
              IF (XMAX.LT.XNORM .AND. XNORM.NE.1.0D0) THEN
262 1 pfleura2
                 XMAX = XNORM
263 1 pfleura2
                 ITERA = I
264 1 pfleura2
                 IONE = INV + I
265 1 pfleura2
              ENDIF
266 1 pfleura2
           END DO
267 1 pfleura2
268 1 pfleura2
           DO I=1,MPLUS
269 1 pfleura2
              INV = (I-1) * MPLUS
270 1 pfleura2
              DO J=1,MPLUS
271 1 pfleura2
                 JNV = (J-1) * MPLUS
272 1 pfleura2
                 IF (J.EQ.ITERA) B(INV + J) = 0.D0
273 1 pfleura2
                 B(JNV + I) = B(INV + J)
274 1 pfleura2
              END DO
275 1 pfleura2
       END DO
276 1 pfleura2
           B(IONE) = 1.0D0
277 1 pfleura2
        END DO
278 1 pfleura2
    END IF ! matches IF (NGEDIIS .NE. MSET(IGeom)) THEN
279 1 pfleura2
280 1 pfleura2
      ! SCALE GEDIIS MATRIX BEFORE INVERSION:
281 1 pfleura2
      DO I=1,MPLUS
282 1 pfleura2
         II = MPLUS * (I-1) + I ! B(II)=diagonal elements of B matrix
283 1 pfleura2
         GSAVE(IGeom,I) = 1.D0 / DSQRT(1.D-20+DABS(B(II)))
284 1 pfleura2
     !Print *, 'GSAVE(',IGeom,',',I,')=', GSAVE(IGeom,I)
285 1 pfleura2
      END DO
286 1 pfleura2
      GSAVE(IGeom,MPLUS) = 1.D0
287 1 pfleura2
      DO I=1,MPLUS
288 1 pfleura2
         DO J=1,MPLUS
289 1 pfleura2
            IJ = MPLUS * (I-1) + J
290 1 pfleura2
            B(IJ) = B(IJ) * GSAVE(IGeom,I) * GSAVE(IGeom,J)
291 1 pfleura2
         END DO
292 1 pfleura2
      END DO
293 1 pfleura2
294 1 pfleura2
     ! INVERT THE GEDIIS MATRIX B:
295 1 pfleura2
    DO I=1, MPLUS
296 1 pfleura2
       !WRITE(*,'(10(1X,F20.4))') B((I-1)*MPLUS+1:I*(MPLUS))
297 1 pfleura2
    END DO
298 1 pfleura2
299 1 pfleura2
      CALL MINV(B,MPLUS,DET) ! matrix inversion.
300 1 pfleura2
301 1 pfleura2
    DO I=1, MPLUS
302 1 pfleura2
       !WRITE(*,'(10(1X,F20.16))') B((I-1)*MPLUS+1:I*(MPLUS))
303 1 pfleura2
    END DO
304 1 pfleura2
305 1 pfleura2
      DO I=1,MPLUS
306 1 pfleura2
         DO J=1,MPLUS
307 1 pfleura2
            IJ = MPLUS * (I-1) + J
308 1 pfleura2
            B(IJ) = B(IJ) * GSAVE(IGeom,I) * GSAVE(IGeom,J)
309 1 pfleura2
         END DO
310 1 pfleura2
      END DO
311 1 pfleura2
312 1 pfleura2
      ! COMPUTE THE NEW INTERPOLATED PARAMETER VECTOR (Geometry):
313 1 pfleura2
    ci=0.d0
314 1 pfleura2
    ci_tmp=0.d0
315 1 pfleura2
316 1 pfleura2
    ci_lt_zero= .FALSE.
317 1 pfleura2
    DO I=1, MSET(IGeom)
318 1 pfleura2
     DO J=1, MSET(IGeom) ! B matrix is read column-wise
319 1 pfleura2
        ci(I)=ci(I)+B((J-1)*(MPLUS)+I)*ESET(IGeom,J) !ESET is energy set.
320 1 pfleura2
     END DO
321 1 pfleura2
     ci(I)=ci(I)+B((MPLUS-1)*(MPLUS)+I)
322 1 pfleura2
     !Print *, 'NO ci < 0 yet, c(',I,')=', ci(I)
323 1 pfleura2
     IF((ci(I) .LT. 0.0D0) .OR. (ci(I) .GT. 1.0D0)) THEN
324 1 pfleura2
       ci_lt_zero=.TRUE.
325 1 pfleura2
       EXIT
326 1 pfleura2
     END IF
327 1 pfleura2
      END DO !matches DO I=1, MSET(IGeom)
328 1 pfleura2
329 1 pfleura2
    IF (ci_lt_zero) Then
330 1 pfleura2
       cis_zero = 0
331 1 pfleura2
         ER_star = 0.D0
332 1 pfleura2
     ER_star_tmp = 1e32
333 1 pfleura2
334 1 pfleura2
     ! B_ij calculations from <B_ij=(g_i-g_j)(R_i-R_j)>, Full B matrix created first and then rows and columns are removed.
335 1 pfleura2
         JJ=0
336 1 pfleura2
         INV=-NFree
337 1 pfleura2
         DO IX=1,MSET(IGeom)
338 1 pfleura2
            INV=INV+NFree
339 1 pfleura2
            JNV=-NFree
340 1 pfleura2
            DO JX=1,MSET(IGeom)
341 1 pfleura2
               JNV=JNV+NFree
342 1 pfleura2
               JJ = JJ + 1
343 1 pfleura2
               BST(JJ)=0.D0
344 1 pfleura2
           DO KX=1, NFree
345 1 pfleura2
              BST(JJ) = BST(JJ) + (((GradSet_free(IGeom,INV+KX)-GradSet_free(IGeom,JNV+KX))* &
346 1 pfleura2
                    (GeomSet_free(IGeom,INV+KX)-GeomSet_free(IGeom,JNV+KX)))/2.D0)
347 1 pfleura2
           END DO
348 1 pfleura2
            END DO
349 1 pfleura2
       END DO
350 1 pfleura2
351 1 pfleura2
     DO I=1, (2**MSET(IGeom))-2 ! all (2**MSET(IGeom))-2 combinations of cis, except the one where all cis are .GT. 0 and .LT. 1
352 1 pfleura2
         !Print *, 'Entering into DO I=1, (2**MSET(IGeom))-2  loop, MSET(IGeom)=', MSET(IGeom), ', I=', I
353 1 pfleura2
         ci(:)=1.D0
354 1 pfleura2
         ! find out which cis are zero in I:
355 1 pfleura2
       DO IX=1, MSET(IGeom)
356 1 pfleura2
          JJ=IAND(I, 2**(IX-1))
357 1 pfleura2
        IF(JJ .EQ. 0) Then
358 1 pfleura2
          ci(IX)=0.D0
359 1 pfleura2
          END IF
360 1 pfleura2
       END DO
361 1 pfleura2
362 1 pfleura2
       ci_lt_zero = .FALSE.
363 1 pfleura2
       ! B_ij calculations from <B_ij=(g_i-g_j)(R_i-R_j)>, Full B matrix created first and then rows and columns are removed.
364 1 pfleura2
       DO IX=1, MSET(IGeom)*MSET(IGeom)
365 1 pfleura2
                B(IX) = BST(IX) !just a copy of the original B (GEDIIS) matrix
366 1 pfleura2
             END DO
367 1 pfleura2
368 1 pfleura2
             ! Removal of KXth row and KXth column in order to accomodate cI to be zero:
369 1 pfleura2
       current_size_B_mat=MSET(IGeom)
370 1 pfleura2
       cis_zero = 0
371 1 pfleura2
       ! The bits of I (index of the upper loop 'DO I=1, (2**MSET(IGeom))-2'), gives which cis are zero.
372 1 pfleura2
       DO KX=1, MSET(IGeom) ! searching for each bit of I (index of the upper loop 'DO I=1, (2**MSET(IGeom))-2')
373 1 pfleura2
          IF (ci(KX) .EQ. 0.D0) Then !remove KXth row and KXth column
374 1 pfleura2
           cis_zero = cis_zero + 1
375 1 pfleura2
376 1 pfleura2
             ! First row removal: (B matrix is read column-wise)
377 1 pfleura2
             JJ=0
378 1 pfleura2
                   DO IX=1,current_size_B_mat ! columns reading
379 1 pfleura2
                      DO JX=1,current_size_B_mat ! rows reading
380 1 pfleura2
                 IF (JX .NE. KX) Then
381 1 pfleura2
                     JJ = JJ + 1
382 1 pfleura2
                     B_tmp(JJ) = B((IX-1)*current_size_B_mat+JX)
383 1 pfleura2
                 END IF
384 1 pfleura2
              END DO
385 1 pfleura2
             END DO
386 1 pfleura2
387 1 pfleura2
             DO IX=1,current_size_B_mat*(current_size_B_mat-1)
388 1 pfleura2
                B(IX) = B_tmp(IX)
389 1 pfleura2
             END DO
390 1 pfleura2
391 1 pfleura2
             ! Now column removal:
392 1 pfleura2
             JJ=0
393 1 pfleura2
                   DO IX=1,KX-1 ! columns reading
394 1 pfleura2
                      DO JX=1,current_size_B_mat-1 ! rows reading
395 1 pfleura2
                 JJ = JJ + 1
396 1 pfleura2
                 B_tmp(JJ) = B(JJ)
397 1 pfleura2
              END DO
398 1 pfleura2
             END DO
399 1 pfleura2
400 1 pfleura2
                   DO IX=KX+1,current_size_B_mat
401 1 pfleura2
                      DO JX=1,current_size_B_mat-1
402 1 pfleura2
                 JJ = JJ + 1
403 1 pfleura2
                    B_tmp(JJ) = B(JJ+current_size_B_mat-1)
404 1 pfleura2
              END DO
405 1 pfleura2
             END DO
406 1 pfleura2
407 1 pfleura2
             DO IX=1,(current_size_B_mat-1)*(current_size_B_mat-1)
408 1 pfleura2
                B(IX) = B_tmp(IX)
409 1 pfleura2
             END DO
410 1 pfleura2
           current_size_B_mat = current_size_B_mat - 1
411 1 pfleura2
        END IF ! matches IF (ci(KX) .EQ. 0.D0) Then !remove
412 1 pfleura2
           END DO ! matches DO KX=1, MSET(IGeom)
413 1 pfleura2
414 1 pfleura2
       ! The following shifting is required to correct indices of B_ij elements in the GEDIIS matrix.
415 1 pfleura2
       ! The correction is needed because the last coloumn and row of the matrix contains all 1 and one zero.
416 1 pfleura2
       DO IX=MSET(IGeom)-cis_zero-1,1,-1
417 1 pfleura2
        DO JX=MSET(IGeom)-cis_zero,1,-1
418 1 pfleura2
           B(IX*(MSET(IGeom)-cis_zero)+JX+IX) = B(IX*(MSET(IGeom)-cis_zero)+JX)
419 1 pfleura2
        END DO
420 1 pfleura2
       END DO
421 1 pfleura2
422 1 pfleura2
       ! for last row and last column of GEDIIS matrix
423 1 pfleura2
       DO IX=1,MPLUS-cis_zero
424 1 pfleura2
        B((MPLUS-cis_zero)*IX) = 1.D0
425 1 pfleura2
        B((MPLUS-cis_zero)*(MSET(IGeom)-cis_zero)+IX) = 1.D0
426 1 pfleura2
       END DO
427 1 pfleura2
       B((MPLUS-cis_zero) * (MPLUS-cis_zero)) = 0.D0
428 1 pfleura2
429 1 pfleura2
           DO IX=1, MPLUS
430 1 pfleura2
              !WRITE(*,'(10(1X,F20.4))') B((IX-1)*MPLUS+1:IX*(MPLUS))
431 1 pfleura2
           END DO
432 1 pfleura2
433 1 pfleura2
       ! ELIMINATE ERROR VECTORS WITH THE LARGEST NORM:
434 1 pfleura2
             IF (NGEDIIS .NE. MSET(IGeom)) THEN
435 1 pfleura2
          JX = min(MSET(IGeom)-NGEDIIS,MSET(IGeom)-cis_zero-1)
436 1 pfleura2
                DO II=1,JX
437 1 pfleura2
                   XMAX = -1.D10
438 1 pfleura2
                   ITERA = 0
439 1 pfleura2
                   DO IX=1,MSET(IGeom)-cis_zero
440 1 pfleura2
                      XNORM = 0.D0
441 1 pfleura2
                      INV = (IX-1) * (MPLUS-cis_zero)
442 1 pfleura2
                      DO J=1,MSET(IGeom)-cis_zero
443 1 pfleura2
                         XNORM = XNORM + ABS(B(INV + J))
444 1 pfleura2
                      END DO
445 1 pfleura2
                      IF (XMAX.LT.XNORM .AND. XNORM.NE.1.0D0) THEN
446 1 pfleura2
                         XMAX = XNORM
447 1 pfleura2
                         ITERA = IX
448 1 pfleura2
                         IONE = INV + IX
449 1 pfleura2
                      ENDIF
450 1 pfleura2
           END DO
451 1 pfleura2
452 1 pfleura2
                   DO IX=1,MPLUS-cis_zero
453 1 pfleura2
                      INV = (IX-1) * (MPLUS-cis_zero)
454 1 pfleura2
                      DO J=1,MPLUS-cis_zero
455 1 pfleura2
                         JNV = (J-1) * (MPLUS-cis_zero)
456 1 pfleura2
                         IF (J.EQ.ITERA) B(INV + J) = 0.D0
457 1 pfleura2
                         B(JNV + IX) = B(INV + J)
458 1 pfleura2
            END DO
459 1 pfleura2
               END DO
460 1 pfleura2
                   B(IONE) = 1.0D0
461 1 pfleura2
          END DO
462 1 pfleura2
           END IF ! matches IF (NGEDIIS .NE. MSET(IGeom)) THEN
463 1 pfleura2
464 1 pfleura2
       ! SCALE GEDIIS MATRIX BEFORE INVERSION:
465 1 pfleura2
       DO IX=1,MPLUS-cis_zero
466 1 pfleura2
        II = (MPLUS-cis_zero) * (IX-1) + IX ! B(II)=diagonal elements of B matrix
467 1 pfleura2
        GSAVE(IGeom,IX) = 1.D0 / DSQRT(1.D-20+DABS(B(II)))
468 1 pfleura2
       END DO
469 1 pfleura2
       GSAVE(IGeom,MPLUS-cis_zero) = 1.D0
470 1 pfleura2
       DO IX=1,MPLUS-cis_zero
471 1 pfleura2
        DO JX=1,MPLUS-cis_zero
472 1 pfleura2
           IJ = (MPLUS-cis_zero) * (IX-1) + JX
473 1 pfleura2
           B(IJ) = B(IJ) * GSAVE(IGeom,IX) * GSAVE(IGeom,JX)
474 1 pfleura2
        END DO
475 1 pfleura2
       END DO
476 1 pfleura2
477 1 pfleura2
       ! INVERT THE GEDIIS MATRIX B:
478 1 pfleura2
       CALL MINV(B,MPLUS-cis_zero,DET) ! matrix inversion.
479 1 pfleura2
480 1 pfleura2
       DO IX=1,MPLUS-cis_zero
481 1 pfleura2
        DO JX=1,MPLUS-cis_zero
482 1 pfleura2
           IJ = (MPLUS-cis_zero) * (IX-1) + JX
483 1 pfleura2
           B(IJ) = B(IJ) * GSAVE(IGeom,IX) * GSAVE(IGeom,JX)
484 1 pfleura2
        END DO
485 1 pfleura2
       END DO
486 1 pfleura2
487 1 pfleura2
           DO IX=1, MPLUS
488 1 pfleura2
              !WRITE(*,'(10(1X,F20.4))') B((IX-1)*MPLUS+1:IX*(MPLUS))
489 1 pfleura2
           END DO
490 1 pfleura2
491 1 pfleura2
             ! ESET is rearranged to handle zero cis and stored in ESET_tmp:
492 1 pfleura2
       JJ=0
493 1 pfleura2
       DO IX=1, MSET(IGeom)
494 1 pfleura2
        IF (ci(IX) .NE. 0) Then
495 1 pfleura2
           JJ=JJ+1
496 1 pfleura2
           ESET_tmp(JJ) = ESET(IGeom,IX)
497 1 pfleura2
        END IF
498 1 pfleura2
       END DO
499 1 pfleura2
500 1 pfleura2
       ! DETERMINATION OF nonzero cis:
501 1 pfleura2
       MyPointer=1
502 1 pfleura2
            DO IX=1, MSET(IGeom)-cis_zero
503 1 pfleura2
          tmp = 0.D0
504 1 pfleura2
            DO J=1, MSET(IGeom)-cis_zero ! B matrix is read column-wise
505 1 pfleura2
           tmp=tmp+B((J-1)*(MPLUS-cis_zero)+IX)*ESET_tmp(J)
506 1 pfleura2
        END DO
507 1 pfleura2
            tmp=tmp+B((MPLUS-cis_zero-1)*(MPLUS-cis_zero)+IX)
508 1 pfleura2
            IF((tmp .LT. 0.0D0) .OR. (tmp .GT. 1.0D0)) THEN
509 1 pfleura2
               ci_lt_zero=.TRUE.
510 1 pfleura2
               EXIT
511 1 pfleura2
        ELSE
512 1 pfleura2
           DO JX=MyPointer,MSET(IGeom)
513 1 pfleura2
              IF (ci(JX) .NE. 0) Then
514 1 pfleura2
                 ci(JX) = tmp
515 1 pfleura2
               MyPointer=JX+1
516 1 pfleura2
               EXIT
517 1 pfleura2
            END IF
518 1 pfleura2
           END DO
519 1 pfleura2
            END IF
520 1 pfleura2
             END DO !matches DO I=1, MSET(IGeom)-cis_zero
521 1 pfleura2
           !Print *, 'Local set of cis, first 10:, MSET(IGeom)=', MSET(IGeom), ', I of (2**MSET(IGeom))-2=', I
522 1 pfleura2
       !WRITE(*,'(10(1X,F20.4))') ci(1:MSET(IGeom))
523 1 pfleura2
           !Print *, 'Local set of cis ends:****************************************'
524 1 pfleura2
525 1 pfleura2
             ! new set of cis determined based on the lower energy (ER_star):
526 1 pfleura2
       IF (.NOT. ci_lt_zero) Then
527 1 pfleura2
                Call Energy_GEDIIS(MRESET,MSET(IGeom),ci,GeomSet_free(IGeom,:),GradSet_free(IGeom,:),ESET(IGeom,:),NFree,ER_star)
528 1 pfleura2
        IF (ER_star .LT. ER_star_tmp) Then
529 1 pfleura2
           ci_tmp=ci
530 1 pfleura2
           ER_star_tmp = ER_star
531 1 pfleura2
        END IF
532 1 pfleura2
             END IF ! matches IF (.NOT. ci_lt_zero) Then
533 1 pfleura2
          END DO !matches DO I=1, (2**K)-2 ! all (2**K)-2 combinations of cis, except the one where all cis are .GT. 0 and .LT. 1
534 1 pfleura2
      ci = ci_tmp
535 1 pfleura2
    END IF! matches IF (ci_lt_zero) Then
536 1 pfleura2
537 1 pfleura2
    Print *, 'Final set of cis, first 10:***********************************'
538 1 pfleura2
    WRITE(*,'(10(1X,F20.4))') ci(1:MSET(IGeom))
539 1 pfleura2
    Print *, 'Final set of cis ends:****************************************'
540 1 pfleura2
    Geom_new_free(:) = 0.D0
541 1 pfleura2
    DO I=1, MSET(IGeom)
542 1 pfleura2
         Geom_new_free(:) = Geom_new_free(:) + (ci(I)*GeomSet_free(IGeom,(I-1)*NFree+1:I*NFree)) !MPLUS=MSET(IGeom)+1
543 1 pfleura2
     ! R_(N+1)=R*+DeltaR:
544 1 pfleura2
     DO J=1, NFree
545 1 pfleura2
      tmp=0.D0
546 1 pfleura2
      DO K=1,NFree
547 1 pfleura2
         ! this can be commented:
548 1 pfleura2
         !tmp=tmp+HFree((J-1)*NFree+K)*GradSet_free(IGeom,(I-1)*NFree+K) ! If Hinv=.False., then we need to invert Hess
549 1 pfleura2
      END DO
550 1 pfleura2
      Geom_new_free(J) = Geom_new_free(J) - (ci(I)*tmp)
551 1 pfleura2
     END DO
552 1 pfleura2
    END DO
553 1 pfleura2
554 1 pfleura2
    Step_free(:) =  Geom_new_free(:) - Geom_free(:)
555 1 pfleura2
556 1 pfleura2
      XNORM = SQRT(DOT_PRODUCT(Step_free,Step_free))
557 1 pfleura2
      IF (PRINT) THEN
558 1 pfleura2
         WRITE (6,'(/10X,''DEVIATION IN X '',F10.4,8X,''DETERMINANT '',G9.3)') XNORM, DET
559 1 pfleura2
         !WRITE(*,'(10X,''GEDIIS COEFFICIENTS'')')
560 1 pfleura2
         !WRITE(*,'(10X,5F12.5)') (B(MPLUS*MSET(IGeom)+I),I=1,MSET(IGeom))
561 1 pfleura2
      ENDIF
562 1 pfleura2
563 1 pfleura2
      ! THE FOLLOWING TOLERENCES FOR XNORM AND DET ARE SOMEWHAT ARBITRARY!
564 1 pfleura2
      THRES = MAX(10.D0**(-NFree), 1.D-25)
565 1 pfleura2
      IF (XNORM.GT.2.D0 .OR. DABS(DET) .LT. THRES) THEN
566 1 pfleura2
         IF (PRINT)THEN
567 1 pfleura2
            WRITE(*,*) "THE GEDIIS MATRIX IS ILL CONDITIONED"
568 1 pfleura2
            WRITE(*,*) " - PROBABLY, VECTORS ARE LINEARLY DEPENDENT - "
569 1 pfleura2
            WRITE(*,*) "THE GEDIIS STEP WILL BE REPEATED WITH A SMALLER SPACE"
570 1 pfleura2
         END IF
571 1 pfleura2
         DO K=1,MM
572 1 pfleura2
      B(K) = BS(K) ! why this is reverted? Because "IF (NGEDIIS .GT. 0) GO TO 80", see below
573 1 pfleura2
         END DO
574 1 pfleura2
         NGEDIIS = NGEDIIS - 1
575 1 pfleura2
         IF (NGEDIIS .GT. 0) GO TO 80
576 1 pfleura2
         IF (PRINT) WRITE(*,'(10X,''NEWTON-RAPHSON STEP TAKEN'')')
577 1 pfleura2
         Geom_new_free(:) = Geom_free(:) ! Geom_new is set to original Geom, thus Step = Geom(:) - Geom_new(:)=zero, the whole
578 1 pfleura2
                           ! new update to Geom_new is discarded, since XNORM.GT.2.D0 .OR. DABS(DET) .LT. THRES
579 1 pfleura2
      END IF ! matches IF (XNORM.GT.2.D0 .OR. DABS(DET).LT. THRES) THEN
580 1 pfleura2
581 1 pfleura2
    !******************************************************************************************************************
582 1 pfleura2
    Geom_new_free(:) = 0.D0
583 1 pfleura2
    DO I=1, MSET(IGeom)
584 1 pfleura2
         Geom_new_free(:) = Geom_new_free(:) + (ci(I)*GeomSet_free(IGeom,(I-1)*NFree+1:I*NFree)) !MPLUS=MSET(IGeom)+1
585 1 pfleura2
     ! R_(N+1)=R*+DeltaR:
586 1 pfleura2
     DO J=1, NFree
587 1 pfleura2
      tmp=0.D0
588 1 pfleura2
      DO K=1,NFree
589 1 pfleura2
         tmp=tmp+HFree((J-1)*NFree+K)*GradSet_free(IGeom,(I-1)*NFree+K) ! If Hinv=.False., then we need to invert Hess
590 1 pfleura2
      END DO
591 1 pfleura2
      Geom_new_free(J) = Geom_new_free(J) - (ci(I)*tmp)
592 1 pfleura2
     END DO
593 1 pfleura2
    END DO
594 1 pfleura2
595 1 pfleura2
    Step_free(:) =  Geom_new_free(:) - Geom_free(:)
596 1 pfleura2
   !******************************************************************************************************************
597 1 pfleura2
    Step = 0.d0
598 1 pfleura2
      DO I=1,NFree
599 1 pfleura2
     Step = Step + Step_free(I)*Vfree(:,I)
600 1 pfleura2
      END DO
601 1 pfleura2
602 1 pfleura2
      DEALLOCATE(Hfree,Htmp,Grad_free,Step_free,Geom_free,Geom_new_free)
603 1 pfleura2
604 1 pfleura2
      IF (PRINT)  WRITE(*,'(/,''       END Step_GEDIIS_ALL   '',/)')
605 1 pfleura2
606 1 pfleura2
      END SUBROUTINE Step_GEDIIS_All