root / src / EgradPath.f90
Historique | Voir | Annoter | Télécharger (30,58 ko)
1 | 2 | pfleura2 | SUBROUTINE EGradPath(Iopt) |
---|---|---|---|
2 | 1 | pfleura2 | |
3 | 1 | pfleura2 | ! This subroutine calculates the energy+gradient for all points of |
4 | 1 | pfleura2 | ! the path |
5 | 1 | pfleura2 | |
6 | 12 | pfleura2 | !---------------------------------------------------------------------- |
7 | 12 | pfleura2 | ! Copyright 2003-2014 Ecole Normale Supérieure de Lyon, |
8 | 12 | pfleura2 | ! Centre National de la Recherche Scientifique, |
9 | 12 | pfleura2 | ! Université Claude Bernard Lyon 1. All rights reserved. |
10 | 12 | pfleura2 | ! |
11 | 12 | pfleura2 | ! This work is registered with the Agency for the Protection of Programs |
12 | 12 | pfleura2 | ! as IDDN.FR.001.100009.000.S.P.2014.000.30625 |
13 | 12 | pfleura2 | ! |
14 | 12 | pfleura2 | ! Authors: P. Fleurat-Lessard, P. Dayal |
15 | 12 | pfleura2 | ! Contact: optnpath@gmail.com |
16 | 12 | pfleura2 | ! |
17 | 12 | pfleura2 | ! This file is part of "Opt'n Path". |
18 | 12 | pfleura2 | ! |
19 | 12 | pfleura2 | ! "Opt'n Path" is free software: you can redistribute it and/or modify |
20 | 12 | pfleura2 | ! it under the terms of the GNU Affero General Public License as |
21 | 12 | pfleura2 | ! published by the Free Software Foundation, either version 3 of the License, |
22 | 12 | pfleura2 | ! or (at your option) any later version. |
23 | 12 | pfleura2 | ! |
24 | 12 | pfleura2 | ! "Opt'n Path" is distributed in the hope that it will be useful, |
25 | 12 | pfleura2 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of |
26 | 12 | pfleura2 | ! |
27 | 12 | pfleura2 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
28 | 12 | pfleura2 | ! GNU Affero General Public License for more details. |
29 | 12 | pfleura2 | ! |
30 | 12 | pfleura2 | ! You should have received a copy of the GNU Affero General Public License |
31 | 12 | pfleura2 | ! along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>. |
32 | 12 | pfleura2 | ! |
33 | 12 | pfleura2 | ! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr, |
34 | 12 | pfleura2 | ! for commercial licensing opportunities. |
35 | 12 | pfleura2 | !---------------------------------------------------------------------- |
36 | 12 | pfleura2 | |
37 | 1 | pfleura2 | use Path_module, only : ProgExe, NCoord, PROG, Nat, NGeomF, Coord, IntCoordF, & |
38 | 1 | pfleura2 | IndZmat, XyzGeomF, Atome, Order, OrderInv, Latr, FrozAtoms, & |
39 | 1 | pfleura2 | Grad, Energies, renum, dzdc, massat, Pi, NCART, OptReac, & |
40 | 8 | pfleura2 | OptProd, AtName, Unit, CalcEReac, CalcEProd |
41 | 1 | pfleura2 | ! GeomOld_all(NGeomF,NCoord), BTransInv (3*Nat,NCoord): used for Baker case |
42 | 1 | pfleura2 | use Io_module |
43 | 1 | pfleura2 | |
44 | 1 | pfleura2 | IMPLICIT NONE |
45 | 1 | pfleura2 | |
46 | 1 | pfleura2 | INTEGER(KINT), INTENT(IN) :: Iopt |
47 | 1 | pfleura2 | |
48 | 1 | pfleura2 | INTEGER(KINT) :: IGeom,IGeom0,IGeomF |
49 | 1 | pfleura2 | CHARACTER(LCHARS) :: Line,FileInp,RunCommand |
50 | 1 | pfleura2 | CHARACTER(SCHARS) :: Line2 |
51 | 1 | pfleura2 | |
52 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: GeomTmp(:),GradTmp(:) ! NCoord |
53 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: GradCart(:) ! 3*Nat |
54 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: PathCart(:,:) ! (3*Nat,NGeomF) |
55 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: x(:), y(:), z(:) |
56 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: x_k(:), y_k(:), z_k(:) |
57 | 1 | pfleura2 | REAL(KREAL), ALLOCATABLE :: GeomCart(:,:) !(Nat,3) |
58 | 1 | pfleura2 | REAL(KREAL) :: E |
59 | 12 | pfleura2 | ! Flag to see if we should print debug info |
60 | 1 | pfleura2 | LOGICAL :: Debug |
61 | 12 | pfleura2 | ! Flags to check if we have found the information we need |
62 | 12 | pfleura2 | LOGICAL, SAVE :: FReadE, FReadGrad |
63 | 1 | pfleura2 | LOGICAL, SAVE :: First=.TRUE. |
64 | 12 | pfleura2 | ! Flag to check if a file/directory exists |
65 | 12 | pfleura2 | LOGICAL :: FExists |
66 | 1 | pfleura2 | |
67 | 1 | pfleura2 | INTEGER(KINT) :: I, Iat, IBeg, Idx |
68 | 1 | pfleura2 | REAL(KREAL) :: RTmp1, RTmp2, RTmp3 |
69 | 1 | pfleura2 | INTEGER(KINT) :: Itmp1, ITmp2 |
70 | 1 | pfleura2 | |
71 | 1 | pfleura2 | INTERFACE |
72 | 1 | pfleura2 | function valid(string) result (isValid) |
73 | 1 | pfleura2 | CHARACTER(*), intent(in) :: string |
74 | 1 | pfleura2 | logical :: isValid |
75 | 1 | pfleura2 | END function VALID |
76 | 1 | pfleura2 | |
77 | 1 | pfleura2 | subroutine Egrad(E,Geom,Grad,NCoord,IGeom,IOpt,GeomCart,FOptGeom,GeomOld,GeomCart_old) |
78 | 1 | pfleura2 | |
79 | 1 | pfleura2 | ! This routines calculates the energy E and the gradient Grad of |
80 | 1 | pfleura2 | ! a molecule with Geometry Geom (may be in internal coordinates), |
81 | 1 | pfleura2 | ! using for now, either Gaussian or Ext, more general later. |
82 | 1 | pfleura2 | |
83 | 1 | pfleura2 | use Path_module, only : Nat,AtName,Coord,dzdc,indzmat,Nom,Atome,massat,unit, & |
84 | 1 | pfleura2 | prog,NCart,XyzGeomF,IntCoordF,BTransInv, XyzGeomI, & |
85 | 1 | pfleura2 | GeomOld_all,BTransInvF,BTransInv_local,UMatF,UMat_local & |
86 | 1 | pfleura2 | , BprimT,a0,ScanCoord, Coordinate,NPrim,BMat_BakerT,BBT,BBT_inv & |
87 | 1 | pfleura2 | , Order,OrderInv, XPrimitiveF |
88 | 1 | pfleura2 | |
89 | 1 | pfleura2 | ! IntCoordF(NGeomF,NCoord),GeomOld_all(NGeomF,NCoord) |
90 | 1 | pfleura2 | ! allocated in Path.f90 |
91 | 1 | pfleura2 | |
92 | 1 | pfleura2 | use Io_module |
93 | 1 | pfleura2 | |
94 | 1 | pfleura2 | ! Energy (calculated if F300K=.F., else estimated) |
95 | 1 | pfleura2 | REAL(KREAL), INTENT (OUT) :: E |
96 | 1 | pfleura2 | ! NCoord: Number of the degrees of freedom |
97 | 1 | pfleura2 | ! IGeom: index of the geometry. |
98 | 1 | pfleura2 | INTEGER(KINT), INTENT (IN) :: NCoord, IGeom, IOpt |
99 | 1 | pfleura2 | ! Geometry at which gradient is calculated (cf Factual also): |
100 | 1 | pfleura2 | REAL(KREAL), INTENT (INOUT) :: Geom(NCoord) |
101 | 1 | pfleura2 | ! Gradient calculated at Geom geometry: |
102 | 1 | pfleura2 | REAL(KREAL), INTENT (OUT) :: Grad(NCoord) |
103 | 1 | pfleura2 | ! Cartesian geometry corresponding to (Internal Geometry) Geom: |
104 | 1 | pfleura2 | REAL(KREAL), INTENT (OUT) :: GeomCart(Nat,3) |
105 | 1 | pfleura2 | !!! Optional, just for geometry optimization with Baker coordinates |
106 | 1 | pfleura2 | REAL(KREAL), INTENT (IN), OPTIONAL :: GeomCart_old(Nat,3) |
107 | 1 | pfleura2 | REAL(KREAL), INTENT (INOUT), OPTIONAL :: GeomOld(NCoord) |
108 | 1 | pfleura2 | ! FOptGeom is a flag indicating if we are doing a geom optimization |
109 | 1 | pfleura2 | ! it can be omitted so that we use a local flag: Flag_Opt_Geom |
110 | 1 | pfleura2 | LOGICAL,INTENT (IN), OPTIONAL :: FOptGeom |
111 | 1 | pfleura2 | ! if FOptGeom is given Flag_Opt_Geom=FOptGeom, else Flag_Opt_Geom=F |
112 | 1 | pfleura2 | LOGICAL :: Flag_Opt_Geom |
113 | 1 | pfleura2 | |
114 | 1 | pfleura2 | END subroutine Egrad |
115 | 1 | pfleura2 | |
116 | 1 | pfleura2 | END INTERFACE |
117 | 1 | pfleura2 | |
118 | 1 | pfleura2 | debug=valid('EGradPath') |
119 | 1 | pfleura2 | |
120 | 1 | pfleura2 | if (debug) Call header("Entering EgradPath") |
121 | 1 | pfleura2 | |
122 | 1 | pfleura2 | ALLOCATE(GeomTmp(NCoord)) |
123 | 1 | pfleura2 | ALLOCATE(x(Nat),y(Nat),z(Nat)) |
124 | 1 | pfleura2 | ALLOCATE(x_k(Nat),y_k(Nat),z_k(Nat)) |
125 | 1 | pfleura2 | |
126 | 1 | pfleura2 | |
127 | 1 | pfleura2 | IF (RunMode=="PARA") THEN ! matches at L315 |
128 | 1 | pfleura2 | |
129 | 1 | pfleura2 | ALLOCATE(GradCart(3*Nat),PathCart(3*Nat,NGeomF)) |
130 | 1 | pfleura2 | |
131 | 1 | pfleura2 | |
132 | 1 | pfleura2 | SELECT CASE(Prog) |
133 | 1 | pfleura2 | CASE ("VASP") |
134 | 1 | pfleura2 | ! We will use the NEB routine of VASP to get all forces at once on Para8 queue |
135 | 1 | pfleura2 | |
136 | 1 | pfleura2 | ! First, we create all the POSCAR into the 0x directories |
137 | 1 | pfleura2 | |
138 | 12 | pfleura2 | WRITE(IOOUT,*) " Creating the POSCAR files " |
139 | 12 | pfleura2 | |
140 | 12 | pfleura2 | |
141 | 1 | pfleura2 | ! With RunMode="Para" one cannot optimize reactants or products (due to Vasp constraints) |
142 | 1 | pfleura2 | ! so this has been done beforehand :-) |
143 | 12 | pfleura2 | |
144 | 12 | pfleura2 | |
145 | 1 | pfleura2 | IGeom0=2 |
146 | 1 | pfleura2 | IGeomF=NGeomF-1 |
147 | 1 | pfleura2 | PathCart=0.d0 |
148 | 1 | pfleura2 | |
149 | 1 | pfleura2 | IF (First) THEN |
150 | 1 | pfleura2 | ! For the first iteration, we do write a POSCAR into 00 and 0Final |
151 | 1 | pfleura2 | ! Vasp needs it ! |
152 | 1 | pfleura2 | First=.FALSE. |
153 | 1 | pfleura2 | IGeom0=1 |
154 | 1 | pfleura2 | IGeomF=NGeomF |
155 | 12 | pfleura2 | ! PFL 2013 OCt -> |
156 | 12 | pfleura2 | ! We might have to create the directories |
157 | 12 | pfleura2 | Line2="00" |
158 | 12 | pfleura2 | DO IGeom=IGeom0,IGeomF |
159 | 12 | pfleura2 | WRITE(Line,'(I5)') IGeom-1 |
160 | 12 | pfleura2 | Idx=2-Len_TRIM(ADJUSTL(Line)) |
161 | 12 | pfleura2 | FileInp=Line2(1:Idx) // TRIM(AdjustL(Line)) |
162 | 12 | pfleura2 | INQUIRE(File=FileInp,Exist=FExists) |
163 | 12 | pfleura2 | If (.NOT.Fexists) THEN |
164 | 12 | pfleura2 | FileInp="mkdir " // TRIM(FileInp) |
165 | 12 | pfleura2 | if (debug) WRITE(*,'(1X,A)') "Creating dir: " // TRIM(FileInp) |
166 | 12 | pfleura2 | CALL SYSTEM(FileInp) |
167 | 12 | pfleura2 | ELSE IF (DEBUG) THEN |
168 | 12 | pfleura2 | WRITE(*,'(1X,A)') "Directory " // TRIM(FileInp) // " already exists" |
169 | 12 | pfleura2 | END IF |
170 | 12 | pfleura2 | END DO |
171 | 12 | pfleura2 | END IF ! If (First) |
172 | 1 | pfleura2 | |
173 | 1 | pfleura2 | Line2="00" |
174 | 1 | pfleura2 | DO IGeom=IGeom0,IGeomF |
175 | 1 | pfleura2 | |
176 | 1 | pfleura2 | ! NOTE THAT HERE IGeom0 (IGeomF) IS NOT 1 (NGeomF) FOR NON-FIRST OPTIMIZATION. |
177 | 1 | pfleura2 | ! THIS WILL HAVE EFFECT ON XyzGeomF UPDATE IN BAKER CASE. |
178 | 1 | pfleura2 | WRITE(Line,'(I5)') IGeom-1 |
179 | 1 | pfleura2 | Idx=2-Len_TRIM(ADJUSTL(Line)) |
180 | 12 | pfleura2 | FileInp="." // FileDelim // Line2(1:Idx) // TRIM(AdjustL(Line)) // FileDelim // "POSCAR" |
181 | 1 | pfleura2 | if (debug) WRITE(*,*) "Creating ",TRIM(FileInp) |
182 | 1 | pfleura2 | OPEN(IOTMP,File=TRIM(FileInp)) |
183 | 1 | pfleura2 | |
184 | 1 | pfleura2 | WRITE(Line,"('Geometry ',I3,'/',I3,' for iteration ',I3)") Igeom,NgeomF,Iopt |
185 | 1 | pfleura2 | if (debug) WRITE(*,*) Line |
186 | 1 | pfleura2 | SELECT CASE(Coord) |
187 | 1 | pfleura2 | CASE ('ZMAT') |
188 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
189 | 1 | pfleura2 | Call Int2cart(nat,indzmat,geomTmp,PathCart(1,IGeom)) |
190 | 1 | pfleura2 | CASE ('BAKER') |
191 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
192 | 1 | pfleura2 | ! |
193 | 1 | pfleura2 | IF (IOpt .EQ. 0) THEN |
194 | 1 | pfleura2 | ! EgradPath(...) is called only after the call of PathCreate(...) |
195 | 1 | pfleura2 | ! and in PathCreate, IF (MOD(IOpt,IReparam).EQ.0) (which is always |
196 | 1 | pfleura2 | ! true if IOpt==0), THEN ExtraPol_baker(...) is called. |
197 | 1 | pfleura2 | ! In ExtraPol_baker(...), XyzGeomF(...) is calculated upon converting |
198 | 1 | pfleura2 | ! the interpolated internal coordinates into cartesian coordinates. |
199 | 1 | pfleura2 | ! Thus here we don't need to reconvert the internal coordinates again |
200 | 1 | pfleura2 | ! to the cartesian coordinates. |
201 | 1 | pfleura2 | DO I=1,Nat |
202 | 1 | pfleura2 | PathCart(3*I-2,IGeom) = XyzGeomF(IGeom,1,I) |
203 | 1 | pfleura2 | PathCart(3*I-1,IGeom) = XyzGeomF(IGeom,2,I) |
204 | 1 | pfleura2 | PathCart(3*I,IGeom) = XyzGeomF(IGeom,3,I) |
205 | 1 | pfleura2 | END DO |
206 | 1 | pfleura2 | ELSE |
207 | 1 | pfleura2 | ! XyzGeomF(NGeomF,3,Nat). |
208 | 1 | pfleura2 | x_k(:) = XyzGeomF(IGeom,1,:) |
209 | 1 | pfleura2 | y_k(:) = XyzGeomF(IGeom,2,:) |
210 | 1 | pfleura2 | z_k(:) = XyzGeomF(IGeom,3,:) |
211 | 1 | pfleura2 | |
212 | 1 | pfleura2 | !Call ConvertBakerInternal_cart(GeomOld_all(IGeom,:),IntCoordI(IGeom,:), & |
213 | 1 | pfleura2 | ! x_k,y_k,z_k,x,y,z) |
214 | 1 | pfleura2 | ! PathCart(3*Nat,NGeomF) |
215 | 1 | pfleura2 | DO I=1,Nat |
216 | 1 | pfleura2 | PathCart(3*I-2,IGeom) = x(I) |
217 | 1 | pfleura2 | PathCart(3*I-1,IGeom) = y(I) |
218 | 1 | pfleura2 | PathCart(3*I,IGeom) = z(I) |
219 | 1 | pfleura2 | END DO |
220 | 1 | pfleura2 | |
221 | 1 | pfleura2 | XyzGeomF(IGeom,1,:)=x(:) |
222 | 1 | pfleura2 | XyzGeomF(IGeom,2,:)=y(:) |
223 | 1 | pfleura2 | XyzGeomF(IGeom,3,:)=z(:) |
224 | 1 | pfleura2 | END IF |
225 | 1 | pfleura2 | CASE ('MIXED') |
226 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
227 | 1 | pfleura2 | Call Mixed2cart(nat,indzmat,geomTmp,PathCart(1,IGeom)) |
228 | 1 | pfleura2 | CASE ('CART','HYBRID') |
229 | 1 | pfleura2 | !!! CAUTION : PFL 29.AUG.2008 -> |
230 | 1 | pfleura2 | ! XyzGeomI/F stored in (NGeomI/F,3,Nat) and NOT (NGeomI/F,Nat,3) |
231 | 1 | pfleura2 | ! This should be made consistent !!!!!!!!!!!!!!! |
232 | 1 | pfleura2 | GeomTmp=Reshape(reshape(XyzGeomF(IGeom,:,:),(/Nat,3/),ORDER=(/2,1/)),(/3*Nat/)) |
233 | 1 | pfleura2 | !!! <- CAUTION : PFL 29.AUG.2008 |
234 | 1 | pfleura2 | PathCart(:,IGeom)=GeomTmp |
235 | 1 | pfleura2 | CASE DEFAULT |
236 | 1 | pfleura2 | WRITE(*,*) "Coord=",TRIM(Coord)," not recognized. EgradPath L134. STOP" |
237 | 1 | pfleura2 | STOP |
238 | 1 | pfleura2 | END SELECT |
239 | 1 | pfleura2 | |
240 | 1 | pfleura2 | IF (debug) THEN |
241 | 1 | pfleura2 | WRITE(*,*) "Calling PrintGeomVasp for geom",Igeom |
242 | 1 | pfleura2 | Call PrintGeom(Line,Nat,NCoord,PathCart(1,Igeom),Coord,IoOut,Atome,Order,OrderInv,IndZmat) |
243 | 1 | pfleura2 | END IF |
244 | 1 | pfleura2 | Call PrintGeomVasp(Line,PathCart(1,IGeom),Latr,FrozAtoms,IoTmp) |
245 | 1 | pfleura2 | |
246 | 1 | pfleura2 | ! IF (debug) THEN |
247 | 1 | pfleura2 | ! WRITE(*,*) 'DBG MAin, COORD=',TRIM(COORD), ' Igeom=',IGeom |
248 | 1 | pfleura2 | ! WRITE(*,'(3(1X,F9.4))') Reshape(Reshape(XyzGeomF(IGeom,:,:),(/Nat,3/),ORDER=(/2,1/)),(/3*Nat/)) |
249 | 1 | pfleura2 | ! WRITE(*,'(3(1X,F9.4))') Reshape(XyzGeomF(IGeom,:,:),(/3*Nat/)) |
250 | 1 | pfleura2 | ! DO Iat=1,Nat |
251 | 1 | pfleura2 | ! WRITE(*,*) XyzGeomF(IGeom,1:3,Iat) |
252 | 1 | pfleura2 | ! END DO |
253 | 1 | pfleura2 | ! END IF |
254 | 1 | pfleura2 | |
255 | 1 | pfleura2 | IF (debug) THEN |
256 | 1 | pfleura2 | WRITE(*,*) 'DBG MAin, COORD=',TRIM(COORD),' Igeom=',IGeom,' Iopt=',Iopt |
257 | 1 | pfleura2 | WRITE(*,*) GeomTmp(1:min(NCoord,9)) |
258 | 1 | pfleura2 | END IF |
259 | 1 | pfleura2 | |
260 | 1 | pfleura2 | CLOSE(IOTMP) |
261 | 1 | pfleura2 | END DO ! matches DO IGeom=IGeom0,IGeomF |
262 | 1 | pfleura2 | |
263 | 1 | pfleura2 | !!!!!!!!!!!!!!!!!!! |
264 | 1 | pfleura2 | ! |
265 | 1 | pfleura2 | ! We calculate the energies and gradients |
266 | 1 | pfleura2 | ! |
267 | 1 | pfleura2 | |
268 | 1 | pfleura2 | RunCommand=Trim(ProgExe) |
269 | 1 | pfleura2 | if (debug) WRITE(*,*) "Launching:", TRIM(RunCommand) |
270 | 12 | pfleura2 | WRITE(IOOUT,*) " Launching VASP to get the Forces " |
271 | 1 | pfleura2 | call system(RunCommand) |
272 | 1 | pfleura2 | if (debug) WRITE(*,*) "Back from:", TRIM(RunCommand) |
273 | 1 | pfleura2 | |
274 | 1 | pfleura2 | !!!!!!!!!!!!! |
275 | 1 | pfleura2 | ! |
276 | 1 | pfleura2 | ! We read the gradients and energies |
277 | 1 | pfleura2 | ! |
278 | 12 | pfleura2 | WRITE(IOOUT,*) " Reading Forces " |
279 | 12 | pfleura2 | |
280 | 1 | pfleura2 | IGeom0=2 |
281 | 1 | pfleura2 | IGeomF=NGeomF-1 |
282 | 1 | pfleura2 | Grad=0.d0 |
283 | 1 | pfleura2 | |
284 | 1 | pfleura2 | Line2="00" |
285 | 1 | pfleura2 | DO IGeom=IGeom0,IGeomF |
286 | 1 | pfleura2 | WRITE(Line,'(I5)') IGeom-1 |
287 | 1 | pfleura2 | Idx=2-Len_TRIM(ADJUSTL(Line)) |
288 | 12 | pfleura2 | FileInp="." // FileDelim // Line2(1:Idx) // TRIM(AdjustL(Line)) // FileDelim // "OUTCAR" |
289 | 1 | pfleura2 | if (debug) WRITE(*,*) "Reading E and forces from ", TRIM(FileInp) |
290 | 1 | pfleura2 | OPEN(IOTMP,File=TRIM(FileInp)) |
291 | 1 | pfleura2 | |
292 | 1 | pfleura2 | ! We first search for the forces |
293 | 1 | pfleura2 | READ(IOTMP,'(A)',END=999,ERR=999) LINE |
294 | 1 | pfleura2 | ! We search for the last part of the OUTCAR file, after wavefunction is converged |
295 | 1 | pfleura2 | DO WHILE (INDEX(LINE,'EDIFF is reached')==0) |
296 | 12 | pfleura2 | READ(IOTMP,'(A)',END=999,ERR=9991) LINE |
297 | 1 | pfleura2 | END DO |
298 | 1 | pfleura2 | |
299 | 12 | pfleura2 | ! PFL oct 2013 -> |
300 | 12 | pfleura2 | ! The order of the things that we search might change |
301 | 12 | pfleura2 | ! so that we look for TOTEN and TOTAL-FORCE in parallel :) |
302 | 12 | pfleura2 | FReadE=.FALSE. |
303 | 12 | pfleura2 | FReadGrad=.FALSE. |
304 | 12 | pfleura2 | DO WHILE (.NOT.(FReadE.AND.FReadGrad)) |
305 | 12 | pfleura2 | READ(IOTMP,'(A)',END=999,ERR=9992) LINE |
306 | 12 | pfleura2 | |
307 | 12 | pfleura2 | IF (INDEX(LINE,'TOTEN')>0) THEN |
308 | 12 | pfleura2 | Line=Line(26:) |
309 | 12 | pfleura2 | READ(LINE,*) E |
310 | 12 | pfleura2 | if (debug) WRITE(*,'(1X,A,I5,A,F20.6)') "DBG Egrad_VASP E(",Igeom,")=",E |
311 | 12 | pfleura2 | Energies(Igeom)=E |
312 | 12 | pfleura2 | FReadE=.TRUE. |
313 | 12 | pfleura2 | END IF |
314 | 12 | pfleura2 | |
315 | 12 | pfleura2 | IF (INDEX(LINE,'TOTAL-FORCE')>0) THEN |
316 | 12 | pfleura2 | READ(IOTMP,'(A)',END=999,ERR=9993) LINE |
317 | 12 | pfleura2 | DO I=1,Nat |
318 | 12 | pfleura2 | Iat=I |
319 | 12 | pfleura2 | IF (renum) Iat=Order(I) |
320 | 12 | pfleura2 | READ(IOTMP,*,END=999,ERR=9994) Rtmp1,RTmp2,RTmp3,GradCart(3*Iat-2:3*Iat) |
321 | 12 | pfleura2 | END DO |
322 | 12 | pfleura2 | FReadGrad=.TRUE. |
323 | 12 | pfleura2 | END IF |
324 | 1 | pfleura2 | END DO |
325 | 1 | pfleura2 | |
326 | 1 | pfleura2 | Close(IoTmp) |
327 | 1 | pfleura2 | |
328 | 12 | pfleura2 | WRITE(IOOUT,*) " Forces read for geometry ",IGeom |
329 | 12 | pfleura2 | |
330 | 1 | pfleura2 | IF (debug) THEN |
331 | 1 | pfleura2 | WRITE(*,*) "DBG EGRAD_VASP - GradCart=FORCES - original order" |
332 | 1 | pfleura2 | DO I=1,Nat |
333 | 1 | pfleura2 | Iat=Order(I) |
334 | 1 | pfleura2 | WRITE(*,'(1X,I5,3(1X,F20.6))') I,GradCart(3*Iat-2:3*Iat) |
335 | 1 | pfleura2 | END DO |
336 | 1 | pfleura2 | WRITE(*,*) "DBG EGRAD_VASP - GradCart=FORCES - internal order" |
337 | 1 | pfleura2 | DO I=1,Nat |
338 | 1 | pfleura2 | Iat=Order(I) |
339 | 1 | pfleura2 | WRITE(*,'(3(1X,I5),3(1X,F20.6))') I,Iat,OrderInv(I),GradCart(3*I-2:3*I) |
340 | 1 | pfleura2 | END DO |
341 | 1 | pfleura2 | END IF |
342 | 1 | pfleura2 | |
343 | 1 | pfleura2 | ! VASP gives the Forces in eV/Angstrom, so we convert it into the |
344 | 1 | pfleura2 | ! gradient in ua/Angstrom |
345 | 1 | pfleura2 | GradCart=-1.0d0*ev2Au*GradCart |
346 | 1 | pfleura2 | |
347 | 12 | pfleura2 | |
348 | 1 | pfleura2 | ! We convert it into the right coordinate system |
349 | 1 | pfleura2 | ! This should be put into a subroutine (because it is also in egradient) |
350 | 1 | pfleura2 | SELECT CASE (COORD) |
351 | 1 | pfleura2 | CASE ("ZMAT") |
352 | 1 | pfleura2 | ALLOCATE(GradTmp(3*Nat)) |
353 | 1 | pfleura2 | if (debug) WRITE(*,*) "DBG EGRAD Calling CalcBMat_int" |
354 | 1 | pfleura2 | CALL CalcBmat_int(nat, PathCart(1:3*Nat,Igeom), indzmat, dzdc, massat,atome) |
355 | 1 | pfleura2 | CALL ConvGrad_Cart2Int(Nat,dzdc,IndZmat,GradCart,GradTmp) |
356 | 1 | pfleura2 | |
357 | 1 | pfleura2 | Grad(IGeom,1)=GradTmp(4) |
358 | 1 | pfleura2 | Grad(IGeom,2)=GradTmp(7) |
359 | 1 | pfleura2 | ! We might have troubles whith angles that are not in the [0:pi] range because, |
360 | 1 | pfleura2 | ! when converted to catesian coord, they do appear in the [0:Pi] range to gaussian... |
361 | 1 | pfleura2 | ! so that the derivative is not good, and a multiplicative factor should be added... |
362 | 1 | pfleura2 | ! |
363 | 1 | pfleura2 | ! This in fact should be taken care of in calc Bmat ... |
364 | 1 | pfleura2 | ! |
365 | 1 | pfleura2 | IF ((IntCoordF(Igeom,3).LE.0).OR.(IntCoordF(Igeom,3).GE.Pi)) THEN |
366 | 1 | pfleura2 | Grad(IGeom,3)=-1.0d0*GradTmp(8) |
367 | 1 | pfleura2 | ELSE |
368 | 1 | pfleura2 | Grad(IGeom,3)=GradTmp(8) |
369 | 1 | pfleura2 | END IF |
370 | 1 | pfleura2 | Idx=4 |
371 | 1 | pfleura2 | DO I=4,Nat |
372 | 1 | pfleura2 | Grad(IGeom,Idx)=GradTmp(3*i-2) |
373 | 1 | pfleura2 | Grad(IGeom,Idx+1)=GradTmp(3*i-1) |
374 | 1 | pfleura2 | Grad(IGeom,Idx+2)=GradTmp(3*i) |
375 | 1 | pfleura2 | Idx=Idx+3 |
376 | 1 | pfleura2 | END DO |
377 | 1 | pfleura2 | DEALLOCATE(GradTmp) |
378 | 1 | pfleura2 | ! We have the gradient in Cartesian coordinates... which is ok for COORD=Cart, Hybrid |
379 | 1 | pfleura2 | ! but we have to convert it into internal coordinates if COORD=Baker |
380 | 1 | pfleura2 | CASE ('BAKER') |
381 | 1 | pfleura2 | ALLOCATE(GradTmp(3*Nat)) |
382 | 1 | pfleura2 | GradTmp=0.d0 |
383 | 1 | pfleura2 | ! Below is to be corrected. |
384 | 1 | pfleura2 | ! This is inside 'IF ((PROG=="VASP").AND.(RunMode=="PARA")) THEN' |
385 | 1 | pfleura2 | !DO I=1, 3*Nat |
386 | 1 | pfleura2 | ! here GradTmp and Grad are gradients in Baker coordinates |
387 | 1 | pfleura2 | ! Not implemented IF ((PROG=="VASP").AND.(RunMode=="PARA")) |
388 | 1 | pfleura2 | ! GradTmp(:) = GradTmp(:) + BTransInv(:,I)*GradCart(I) |
389 | 1 | pfleura2 | !END DO |
390 | 1 | pfleura2 | Grad(IGeom,:) = GradTmp |
391 | 1 | pfleura2 | DEALLOCATE(GradTmp) |
392 | 1 | pfleura2 | CASE ('MIXED') |
393 | 1 | pfleura2 | ALLOCATE(GradTmp(3*Nat)) |
394 | 1 | pfleura2 | if (debug) WRITE(*,*) "DBG EGRAD Calling CalcBMat_mixed" |
395 | 1 | pfleura2 | CALL CalcBMat_mixed (nat,PathCart(1:3*Nat,IGeom), indzmat, dzdc, massat,atome) |
396 | 1 | pfleura2 | CALL ConvGrad_Cart2Int(Nat,dzdc,IndZmat,GradCart,GradTmp) |
397 | 1 | pfleura2 | DO I=1,Nat |
398 | 1 | pfleura2 | WRITE(*,'(1X,3(1X,F10.5))') GradTmp(3*I-2:3*I) |
399 | 1 | pfleura2 | END DO |
400 | 1 | pfleura2 | SELECT CASE (NCART) |
401 | 1 | pfleura2 | CASE (1) |
402 | 1 | pfleura2 | Grad(IGeom,1:3)=GradTmp(1:3) |
403 | 1 | pfleura2 | Grad(IGeom,4)=GradTmp(4) |
404 | 1 | pfleura2 | Grad(IGeom,5)=GradTmp(7) |
405 | 1 | pfleura2 | Grad(IGeom,6)=GradTmp(8) |
406 | 1 | pfleura2 | Idx=7 |
407 | 1 | pfleura2 | IBeg=4 |
408 | 1 | pfleura2 | CASE(2) |
409 | 1 | pfleura2 | Grad(IGeom,1:3)=GradTmp(1:3) |
410 | 1 | pfleura2 | Grad(IGeom,4:6)=GradTmp(4:6) |
411 | 1 | pfleura2 | Grad(IGeom,7)=GradTmp(7) |
412 | 1 | pfleura2 | Grad(IGeom,8)=GradTmp(8) |
413 | 1 | pfleura2 | Idx=9 |
414 | 1 | pfleura2 | IBeg=4 |
415 | 1 | pfleura2 | CASE DEFAULT |
416 | 1 | pfleura2 | Idx=1 |
417 | 1 | pfleura2 | IBeg=1 |
418 | 1 | pfleura2 | END SELECT |
419 | 1 | pfleura2 | DO I=IBeg,Nat |
420 | 1 | pfleura2 | Grad(IGeom,Idx)=GradTmp(3*i-2) |
421 | 1 | pfleura2 | Grad(IGeom,Idx+1)=GradTmp(3*i-1) |
422 | 1 | pfleura2 | Grad(IGeom,Idx+2)=GradTmp(3*i) |
423 | 1 | pfleura2 | Idx=Idx+3 |
424 | 1 | pfleura2 | END DO |
425 | 1 | pfleura2 | DEALLOCATE(GradTmp) |
426 | 1 | pfleura2 | CASE ("CART","HYBRID") |
427 | 1 | pfleura2 | Grad(IGeom,:)=GradCart |
428 | 1 | pfleura2 | CASE DEFAULT |
429 | 1 | pfleura2 | WRITE(*,*) "Coord=",AdjustL(Trim(COORD))," not yet implemented in EgradPath.L257.STOP" |
430 | 1 | pfleura2 | STOP |
431 | 1 | pfleura2 | END SELECT |
432 | 1 | pfleura2 | |
433 | 1 | pfleura2 | IF (debug) THEN |
434 | 1 | pfleura2 | Line='DBG Path, GradTmp' |
435 | 1 | pfleura2 | Call PrintGeom(Line,Nat,NCoord,reshape(GRad(Igeom,1:NCoord),(/Ncoord/)),Coord,6,Atome,Order,OrderInv,IndZmat) |
436 | 1 | pfleura2 | END IF |
437 | 1 | pfleura2 | END DO ! matches DO IGeom=IGeom0,IGeomF |
438 | 12 | pfleura2 | |
439 | 12 | pfleura2 | CASE ("GAUSSIAN") |
440 | 1 | pfleura2 | |
441 | 12 | pfleura2 | WRITE(IOOUT,*) " Creating the Input files and ListJob" |
442 | 1 | pfleura2 | ! We create the output files. |
443 | 1 | pfleura2 | IGeom0=2 |
444 | 1 | pfleura2 | IGeomF=NGeomF-1 |
445 | 1 | pfleura2 | PathCart=0.d0 |
446 | 1 | pfleura2 | |
447 | 1 | pfleura2 | If (OptReac) IGeom0=1 |
448 | 1 | pfleura2 | If (OptProd) IGeomF=NGeomF |
449 | 1 | pfleura2 | |
450 | 1 | pfleura2 | OPEN(IOTMP2,File="ListJob",STATUS="REPLACE") |
451 | 1 | pfleura2 | |
452 | 1 | pfleura2 | DO Igeom=IGeom0,IGeomF |
453 | 1 | pfleura2 | ! We first convert the geometry into a Cartesian one |
454 | 1 | pfleura2 | SELECT CASE(Coord) |
455 | 1 | pfleura2 | CASE ('ZMAT') |
456 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
457 | 1 | pfleura2 | Call Int2cart(nat,indzmat,geomTmp,PathCart(1,IGeom)) |
458 | 1 | pfleura2 | CASE ('BAKER') |
459 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
460 | 1 | pfleura2 | ! |
461 | 1 | pfleura2 | IF (IOpt .EQ. 0) THEN |
462 | 1 | pfleura2 | ! EgradPath(...) is called only after the call of PathCreate(...) |
463 | 1 | pfleura2 | ! and in PathCreate, IF (MOD(IOpt,IReparam).EQ.0) (which is always |
464 | 1 | pfleura2 | ! true if IOpt==0), THEN ExtraPol_baker(...) is called. |
465 | 1 | pfleura2 | ! In ExtraPol_baker(...), XyzGeomF(...) is calculated upon converting |
466 | 1 | pfleura2 | ! the interpolated internal coordinates into cartesian coordinates. |
467 | 1 | pfleura2 | ! Thus here we don't need to reconvert the internal coordinates again |
468 | 1 | pfleura2 | ! to the cartesian coordinates. |
469 | 1 | pfleura2 | DO I=1,Nat |
470 | 1 | pfleura2 | PathCart(3*I-2,IGeom) = XyzGeomF(IGeom,1,I) |
471 | 1 | pfleura2 | PathCart(3*I-1,IGeom) = XyzGeomF(IGeom,2,I) |
472 | 1 | pfleura2 | PathCart(3*I,IGeom) = XyzGeomF(IGeom,3,I) |
473 | 1 | pfleura2 | END DO |
474 | 1 | pfleura2 | ELSE |
475 | 1 | pfleura2 | |
476 | 1 | pfleura2 | ! XyzGeomF(NGeomF,3,Nat). |
477 | 1 | pfleura2 | x_k(:) = XyzGeomF(IGeom,1,:) |
478 | 1 | pfleura2 | y_k(:) = XyzGeomF(IGeom,2,:) |
479 | 1 | pfleura2 | z_k(:) = XyzGeomF(IGeom,3,:) |
480 | 1 | pfleura2 | |
481 | 1 | pfleura2 | !Call ConvertBakerInternal_cart(GeomOld_all(IGeom,:),IntCoordI(IGeom,:), & |
482 | 1 | pfleura2 | ! x_k,y_k,z_k,x,y,z) |
483 | 1 | pfleura2 | ! PathCart(3*Nat,NGeomF) |
484 | 1 | pfleura2 | DO I=1,Nat |
485 | 1 | pfleura2 | PathCart(3*I-2,IGeom) = x(I) |
486 | 1 | pfleura2 | PathCart(3*I-1,IGeom) = y(I) |
487 | 1 | pfleura2 | PathCart(3*I,IGeom) = z(I) |
488 | 1 | pfleura2 | END DO |
489 | 1 | pfleura2 | |
490 | 1 | pfleura2 | XyzGeomF(IGeom,1,:)=x(:) |
491 | 1 | pfleura2 | XyzGeomF(IGeom,2,:)=y(:) |
492 | 1 | pfleura2 | XyzGeomF(IGeom,3,:)=z(:) |
493 | 1 | pfleura2 | END IF |
494 | 1 | pfleura2 | CASE ('MIXED') |
495 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
496 | 1 | pfleura2 | Call Mixed2cart(nat,indzmat,geomTmp,PathCart(1,IGeom)) |
497 | 1 | pfleura2 | CASE ('CART','HYBRID') |
498 | 1 | pfleura2 | !!! CAUTION : PFL 29.AUG.2008 -> |
499 | 1 | pfleura2 | ! XyzGeomI/F stored in (NGeomI/F,3,Nat) and NOT (NGeomI/F,Nat,3) |
500 | 1 | pfleura2 | ! This should be made consistent !!!!!!!!!!!!!!! |
501 | 1 | pfleura2 | GeomTmp=Reshape(reshape(XyzGeomF(IGeom,:,:),(/Nat,3/),ORDER=(/2,1/)),(/3*Nat/)) |
502 | 1 | pfleura2 | !!! <- CAUTION : PFL 29.AUG.2008 |
503 | 1 | pfleura2 | PathCart(:,IGeom)=GeomTmp |
504 | 1 | pfleura2 | CASE DEFAULT |
505 | 1 | pfleura2 | WRITE(*,*) "Coord=",TRIM(Coord)," not recognized. EgradPath L75.STOP" |
506 | 1 | pfleura2 | STOP |
507 | 1 | pfleura2 | END SELECT |
508 | 1 | pfleura2 | |
509 | 1 | pfleura2 | WRITE(Line,'(I5)') IGeom |
510 | 1 | pfleura2 | FileInp="PathPar_" // TRIM(ADJUSTL(Line)) // ".com" |
511 | 1 | pfleura2 | WRITE(IOTMP2,'(1X,A)') TRIM(FileInp) |
512 | 1 | pfleura2 | |
513 | 1 | pfleura2 | OPEN(IOTMP,File=TRIM(FileInp)) |
514 | 1 | pfleura2 | Current => Gauss_root |
515 | 1 | pfleura2 | DO WHILE (ASSOCIATED(Current%next)) |
516 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A)') Trim(current%line) |
517 | 1 | pfleura2 | WRITE(*,'(1X,A)') Trim(current%line) |
518 | 1 | pfleura2 | Current => current%next |
519 | 1 | pfleura2 | END DO |
520 | 1 | pfleura2 | |
521 | 1 | pfleura2 | WRITE(IOTMP,*) |
522 | 1 | pfleura2 | |
523 | 1 | pfleura2 | Current => Gauss_comment |
524 | 1 | pfleura2 | DO WHILE (ASSOCIATED(Current%next)) |
525 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A)') Trim(current%line) |
526 | 1 | pfleura2 | Current => current%next |
527 | 1 | pfleura2 | END DO |
528 | 1 | pfleura2 | ! WRITE(IOOUT,"('Geometry ',I3,'/',I3,' for iteration ',I3)"),Igeom,NgeomF,Iopt |
529 | 1 | pfleura2 | |
530 | 1 | pfleura2 | WRITE(IOTMP,*) |
531 | 1 | pfleura2 | WRITE(IOTMP,*) Trim(Gauss_charge) |
532 | 1 | pfleura2 | |
533 | 1 | pfleura2 | DO I=1,Nat |
534 | 1 | pfleura2 | If (renum) THEN |
535 | 1 | pfleura2 | Iat=Order(I) |
536 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A10,3(1X,F15.8),A)') Trim(AtName(I)) & |
537 | 1 | pfleura2 | ,PathCart(Iat,IGeom),PathCart(Iat+Nat,IGeom),PathCart(Iat+2*Nat,IGeom) & |
538 | 1 | pfleura2 | ,TRIM(Gauss_Paste(I)) |
539 | 1 | pfleura2 | ELSE |
540 | 1 | pfleura2 | Iat=OrderInv(I) |
541 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A10,3(1X,F15.8),A)') Trim(AtName(Iat)) & |
542 | 1 | pfleura2 | ,PathCart(Iat,IGeom),PathCart(Iat+Nat,IGeom),PathCart(Iat+2*Nat,IGeom) & |
543 | 1 | pfleura2 | , TRIM(Gauss_Paste(Iat)) |
544 | 1 | pfleura2 | END IF |
545 | 1 | pfleura2 | END DO |
546 | 1 | pfleura2 | |
547 | 1 | pfleura2 | WRITE(IOTMP,*) |
548 | 1 | pfleura2 | |
549 | 1 | pfleura2 | Current => Gauss_End |
550 | 1 | pfleura2 | DO WHILE (ASSOCIATED(Current%next)) |
551 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A)') Trim(current%line) |
552 | 1 | pfleura2 | Current => current%next |
553 | 1 | pfleura2 | END DO |
554 | 1 | pfleura2 | |
555 | 1 | pfleura2 | WRITE(IOTMP,*) |
556 | 1 | pfleura2 | CLOSE(IOTMP) |
557 | 1 | pfleura2 | |
558 | 1 | pfleura2 | END DO |
559 | 1 | pfleura2 | |
560 | 1 | pfleura2 | Close(IOTMP2) |
561 | 1 | pfleura2 | |
562 | 1 | pfleura2 | ! We launch the parallel calculations |
563 | 1 | pfleura2 | |
564 | 1 | pfleura2 | Line=Trim(ProgExe) |
565 | 1 | pfleura2 | IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(Line) |
566 | 12 | pfleura2 | WRITE(IOOUT,*) " Launching Gaussian to compute the forces" |
567 | 1 | pfleura2 | call system(Line) |
568 | 1 | pfleura2 | IF (DEBUG) WRITE(*,*)'Back from RunCommand:' |
569 | 1 | pfleura2 | |
570 | 1 | pfleura2 | ! We read energies and gradients from the log files |
571 | 1 | pfleura2 | |
572 | 12 | pfleura2 | WRITE(*,*) "Gaussian over. Reading the Forces" |
573 | 12 | pfleura2 | |
574 | 1 | pfleura2 | DO IGeom=IGeom0,IGeomF |
575 | 1 | pfleura2 | WRITE(Line,'(I5)') IGeom |
576 | 1 | pfleura2 | FileInp="PathPar_" // TRIM(ADJUSTL(Line)) // ".log" |
577 | 1 | pfleura2 | OPEN(IOTMP,File=TRIM(FileInp)) |
578 | 1 | pfleura2 | |
579 | 1 | pfleura2 | ! We first search for the forces |
580 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
581 | 1 | pfleura2 | DO WHILE (INDEX(LINE,'Forces (Hartrees/Bohr)')==0) |
582 | 1 | pfleura2 | ! here, we have a problem, because there are so many ways to write E for Gaussian :-) |
583 | 1 | pfleura2 | ! but we do not really need E except if we want to weight the different points |
584 | 1 | pfleura2 | ! on the path... that will be done later :-) |
585 | 1 | pfleura2 | ! And I should use ConvGaussXmol ... |
586 | 1 | pfleura2 | ! PFL 3rd Jan 2008 -> |
587 | 1 | pfleura2 | ! I have finally upgraded the energy reading phase so that it should be able to read |
588 | 1 | pfleura2 | ! many formats corresponding to HF, DFT, MP2, AM1, PM3, ONIOM with or w/o PCM |
589 | 1 | pfleura2 | ! |
590 | 1 | pfleura2 | |
591 | 1 | pfleura2 | ! For MP2, energy is after the last = |
592 | 1 | pfleura2 | IF ((Line(2:9)=="E2 = ")) THEN |
593 | 1 | pfleura2 | Itmp1=Index(LINE,"=",BACK=.TRUE.)+1 |
594 | 1 | pfleura2 | READ(LINE(Itmp1:),*) e |
595 | 1 | pfleura2 | END IF |
596 | 1 | pfleura2 | ! For other methods, it is after the first = |
597 | 1 | pfleura2 | IF ((LINE(2:9)=='Energy= ').OR.(Line(2:9)=="SCF Done").OR.(Line(2:9)=="ONIOM: e") & |
598 | 1 | pfleura2 | .OR.(Line(2:9)==" with al") & |
599 | 1 | pfleura2 | .OR.(Line(2:31)=="Counterpoise: corrected energy")) THEN |
600 | 1 | pfleura2 | Itmp1=Index(LINE,"=")+1 |
601 | 1 | pfleura2 | READ(LINE(Itmp1:),*) e |
602 | 1 | pfleura2 | END IF |
603 | 1 | pfleura2 | ! <- PFL 3 Jan 2008 |
604 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
605 | 1 | pfleura2 | END DO |
606 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
607 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
608 | 1 | pfleura2 | DO I=1,Nat |
609 | 1 | pfleura2 | Iat=I |
610 | 1 | pfleura2 | IF (renum) Iat=Order(I) |
611 | 1 | pfleura2 | READ(IOTMP,*) Itmp1,ITmp2, GradCart(3*Iat-2:3*Iat) |
612 | 1 | pfleura2 | END DO |
613 | 1 | pfleura2 | |
614 | 1 | pfleura2 | ! Gaussian gives the Forces in ua/a0, so we convert it into the |
615 | 1 | pfleura2 | ! gradient in ua/Angstrom |
616 | 1 | pfleura2 | gradCart=-1.0d0*Unit*gradCart |
617 | 1 | pfleura2 | CLOSE(IOTMP) |
618 | 1 | pfleura2 | Energies(IGeom)=E |
619 | 12 | pfleura2 | |
620 | 12 | pfleura2 | WRITE(*,*) "Forces and energy read for geometry ",IGeom |
621 | 1 | pfleura2 | ! We convert it into the right coordinate system |
622 | 1 | pfleura2 | ! This should be put into a subroutine (because it is also in egradient) |
623 | 1 | pfleura2 | SELECT CASE (COORD) |
624 | 1 | pfleura2 | CASE ("ZMAT") |
625 | 1 | pfleura2 | ALLOCATE(GradTmp(3*Nat)) |
626 | 1 | pfleura2 | if (debug) WRITE(*,*) "DBG EGRAD Calling CalcBMat_int" |
627 | 1 | pfleura2 | CALL CalcBMat_int (nat, PathCart(1:3*Nat,Igeom), indzmat, dzdc, massat,atome) |
628 | 1 | pfleura2 | CALL ConvGrad_Cart2Int(Nat,dzdc,IndZmat,GradCart,GradTmp) |
629 | 1 | pfleura2 | |
630 | 1 | pfleura2 | Grad(IGeom,1)=GradTmp(4) |
631 | 1 | pfleura2 | Grad(IGeom,2)=GradTmp(7) |
632 | 1 | pfleura2 | ! We might have troubles whith angles that are not in the [0:pi] range because, |
633 | 1 | pfleura2 | ! when converted to catesian coord, they do appear in the [0:Pi] range to gaussian... |
634 | 1 | pfleura2 | ! so that the derivative is not good, and a multiplicative factor should be added... |
635 | 1 | pfleura2 | ! |
636 | 1 | pfleura2 | ! This in fact should be taken care of in B-mat calculation... |
637 | 1 | pfleura2 | ! |
638 | 1 | pfleura2 | IF ((IntCoordF(Igeom,3).LE.0).OR.(IntCoordF(Igeom,3).GE.Pi)) THEN |
639 | 1 | pfleura2 | Grad(IGeom,3)=-1.0d0*GradTmp(8) |
640 | 1 | pfleura2 | ELSE |
641 | 1 | pfleura2 | Grad(IGeom,3)=GradTmp(8) |
642 | 1 | pfleura2 | END IF |
643 | 1 | pfleura2 | Idx=4 |
644 | 1 | pfleura2 | DO I=4,Nat |
645 | 1 | pfleura2 | Grad(IGeom,Idx)=GradTmp(3*i-2) |
646 | 1 | pfleura2 | Grad(IGeom,Idx+1)=GradTmp(3*i-1) |
647 | 1 | pfleura2 | Grad(IGeom,Idx+2)=GradTmp(3*i) |
648 | 1 | pfleura2 | Idx=Idx+3 |
649 | 1 | pfleura2 | END DO |
650 | 1 | pfleura2 | DEALLOCATE(GradTmp) |
651 | 1 | pfleura2 | ! We have the gradient in Cartesian coordinates... which is ok for COORD=Cart, Hybrid |
652 | 1 | pfleura2 | ! but we have to convert it into internal coordinates if COORD=Baker |
653 | 1 | pfleura2 | CASE ('BAKER') |
654 | 1 | pfleura2 | ALLOCATE(GradTmp(3*Nat)) |
655 | 1 | pfleura2 | GradTmp=0.d0 |
656 | 1 | pfleura2 | ! Below is to be corrected. |
657 | 1 | pfleura2 | ! This is inside 'IF ((PROG=="VASP").AND.(RunMode=="PARA")) THEN' |
658 | 1 | pfleura2 | !DO I=1, 3*Nat |
659 | 1 | pfleura2 | ! here GradTmp and Grad are gradients in Baker coordinates |
660 | 1 | pfleura2 | ! Not implemented IF ((PROG=="VASP").AND.(RunMode=="PARA")) |
661 | 1 | pfleura2 | ! GradTmp(:) = GradTmp(:) + BTransInv(:,I)*GradCart(I) |
662 | 1 | pfleura2 | !END DO |
663 | 1 | pfleura2 | Grad(IGeom,:) = GradTmp |
664 | 1 | pfleura2 | DEALLOCATE(GradTmp) |
665 | 1 | pfleura2 | CASE ('MIXED') |
666 | 1 | pfleura2 | ALLOCATE(GradTmp(3*Nat)) |
667 | 1 | pfleura2 | if (debug) WRITE(*,*) "DBG EGRAD Calling CalcBMat_mixed" |
668 | 1 | pfleura2 | CALL CalcBMat_mixed(nat,PathCart(1:3*Nat,IGeom), indzmat, dzdc, massat,atome) |
669 | 1 | pfleura2 | CALL ConvGrad_Cart2Int(Nat,dzdc,IndZmat,GradCart,GradTmp) |
670 | 1 | pfleura2 | DO I=1,Nat |
671 | 1 | pfleura2 | WRITE(*,'(1X,3(1X,F10.5))') GradTmp(3*I-2:3*I) |
672 | 1 | pfleura2 | END DO |
673 | 1 | pfleura2 | SELECT CASE (NCART) |
674 | 1 | pfleura2 | CASE (1) |
675 | 1 | pfleura2 | Grad(IGeom,1:3)=GradTmp(1:3) |
676 | 1 | pfleura2 | Grad(IGeom,4)=GradTmp(4) |
677 | 1 | pfleura2 | Grad(IGeom,5)=GradTmp(7) |
678 | 1 | pfleura2 | Grad(IGeom,6)=GradTmp(8) |
679 | 1 | pfleura2 | Idx=7 |
680 | 1 | pfleura2 | IBeg=4 |
681 | 1 | pfleura2 | CASE(2) |
682 | 1 | pfleura2 | Grad(IGeom,1:3)=GradTmp(1:3) |
683 | 1 | pfleura2 | Grad(IGeom,4:6)=GradTmp(4:6) |
684 | 1 | pfleura2 | Grad(IGeom,7)=GradTmp(7) |
685 | 1 | pfleura2 | Grad(IGeom,8)=GradTmp(8) |
686 | 1 | pfleura2 | Idx=9 |
687 | 1 | pfleura2 | IBeg=4 |
688 | 1 | pfleura2 | CASE DEFAULT |
689 | 1 | pfleura2 | Idx=1 |
690 | 1 | pfleura2 | IBeg=1 |
691 | 1 | pfleura2 | END SELECT |
692 | 1 | pfleura2 | DO I=IBeg,Nat |
693 | 1 | pfleura2 | Grad(IGeom,Idx)=GradTmp(3*i-2) |
694 | 1 | pfleura2 | Grad(IGeom,Idx+1)=GradTmp(3*i-1) |
695 | 1 | pfleura2 | Grad(IGeom,Idx+2)=GradTmp(3*i) |
696 | 1 | pfleura2 | Idx=Idx+3 |
697 | 1 | pfleura2 | END DO |
698 | 1 | pfleura2 | DEALLOCATE(GradTmp) |
699 | 1 | pfleura2 | CASE ("CART","HYBRID") |
700 | 1 | pfleura2 | Grad(IGeom,:)=GradCart |
701 | 1 | pfleura2 | CASE DEFAULT |
702 | 1 | pfleura2 | WRITE(*,*) "Coord=",AdjustL(Trim(COORD))," not yet implemented in EgradPath.L257.STOP" |
703 | 1 | pfleura2 | STOP |
704 | 1 | pfleura2 | END SELECT |
705 | 1 | pfleura2 | |
706 | 1 | pfleura2 | IF (debug) THEN |
707 | 1 | pfleura2 | Line='DBG Path, GradTmp' |
708 | 1 | pfleura2 | Call PrintGeom(Line,Nat,NCoord,reshape(GRad(Igeom,1:NCoord),(/Ncoord/)),Coord,6,Atome,Order,OrderInv,IndZmat) |
709 | 1 | pfleura2 | END IF |
710 | 1 | pfleura2 | |
711 | 1 | pfleura2 | END DO |
712 | 1 | pfleura2 | |
713 | 1 | pfleura2 | CASE DEFAULT |
714 | 1 | pfleura2 | WRITE(*,*) "Prog=",TRIM(Prog)," not yet supported for Para mode." |
715 | 1 | pfleura2 | STOP |
716 | 1 | pfleura2 | END SELECT |
717 | 1 | pfleura2 | |
718 | 1 | pfleura2 | DEALLOCATE(GradCart,PathCart) |
719 | 1 | pfleura2 | |
720 | 12 | pfleura2 | |
721 | 4 | pfleura2 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
722 | 4 | pfleura2 | ! |
723 | 4 | pfleura2 | ! Serial executions |
724 | 4 | pfleura2 | ! |
725 | 4 | pfleura2 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
726 | 1 | pfleura2 | ELSE ! matches IF (RunMode=="PARA") THEN |
727 | 1 | pfleura2 | ! We will launch all calculations sequentially |
728 | 1 | pfleura2 | ALLOCATE(GradTmp(NCoord)) |
729 | 1 | pfleura2 | ! We have the new path, we calculate its energy and gradient |
730 | 1 | pfleura2 | IGeom0=2 |
731 | 1 | pfleura2 | IGeomF=NGeomF-1 |
732 | 1 | pfleura2 | IF (OptReac.OR.CalcEReac) IGeom0=1 |
733 | 1 | pfleura2 | IF (OptProd.OR.CalcEProd) IGeomF=NGeomF |
734 | 1 | pfleura2 | CalcEReac=.FALSE. |
735 | 1 | pfleura2 | CalcEProd=.FALSE. |
736 | 1 | pfleura2 | ALLOCATE(GeomCart(Nat,3)) |
737 | 1 | pfleura2 | |
738 | 1 | pfleura2 | DO IGeom=IGeom0,IGeomF |
739 | 1 | pfleura2 | WRITE(Line,"('Geometry ',I3,'/',I3,' for iteration ',I3)") Igeom,NgeomF,Iopt |
740 | 12 | pfleura2 | WRITE(IOOUT,*) "Creating Input for ", TRIM(Line) |
741 | 12 | pfleura2 | ! if (debug) WRITE(*,*) Line |
742 | 1 | pfleura2 | SELECT CASE(Coord) |
743 | 1 | pfleura2 | CASE ('ZMAT','MIXED') |
744 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
745 | 1 | pfleura2 | CASE ('BAKER') |
746 | 1 | pfleura2 | GeomTmp=IntCoordF(IGeom,:) |
747 | 1 | pfleura2 | CASE ('CART','HYBRID') |
748 | 1 | pfleura2 | !!! CAUTION : PFL 29.AUG.2008 -> |
749 | 1 | pfleura2 | ! XyzGeomI/F stored in (NGeomI/F,3,Nat) and NOT (NGeomI/F,Nat,3) |
750 | 1 | pfleura2 | ! This should be made consistent !!!!!!!!!!!!!!! |
751 | 1 | pfleura2 | GeomTmp=Reshape(reshape(XyzGeomF(IGeom,:,:),(/Nat,3/),ORDER=(/2,1/)),(/3*Nat/)) |
752 | 1 | pfleura2 | ! GeomTmp=Reshape(XyzGeomF(IGeom,:,:),(/3*Nat/)) |
753 | 1 | pfleura2 | !!! <- CAUTION : PFL 29.AUG.2008 |
754 | 1 | pfleura2 | CASE DEFAULT |
755 | 1 | pfleura2 | WRITE(*,*) "Coord=",TRIM(Coord)," not recognized in EgradPath. L299" |
756 | 1 | pfleura2 | STOP |
757 | 1 | pfleura2 | END SELECT |
758 | 1 | pfleura2 | |
759 | 1 | pfleura2 | IF (debug) THEN |
760 | 1 | pfleura2 | WRITE(*,*) 'DBG Main, COORD=',TRIM(COORD),' IGeom=',IGeom,' Iopt=',Iopt |
761 | 1 | pfleura2 | WRITE(*,*) GeomTmp(1:min(NCoord,12)) |
762 | 1 | pfleura2 | END IF |
763 | 1 | pfleura2 | |
764 | 1 | pfleura2 | GeomCart=0.d0 |
765 | 1 | pfleura2 | |
766 | 1 | pfleura2 | ! GradTmp is gradient and calculated in Egrad.F, INTENT(OUT). |
767 | 1 | pfleura2 | ! GeomCart has INTENT(OUT) |
768 | 12 | pfleura2 | |
769 | 12 | pfleura2 | WRITE(IOOUT,*) "Computing and reading forces for ", TRIM(Line) |
770 | 1 | pfleura2 | Call EGrad(E,GeomTmp,GradTmp,NCoord,IGeom,IOpt,GeomCart) !GeomTmp=IntCoordF(IGeom,:) |
771 | 1 | pfleura2 | |
772 | 1 | pfleura2 | ! Egrad calls ConvertBakerInternal_cart to convert GeomTmp=IntCoordF(IGeom,:) into |
773 | 1 | pfleura2 | ! GeomCart (Cartesian Coordinates) so that energy and gradients can be calculated |
774 | 1 | pfleura2 | ! for each geometry. |
775 | 1 | pfleura2 | XyzGeomF(IGeom,1,:)=GeomCart(:,1) |
776 | 1 | pfleura2 | XyzGeomF(IGeom,2,:)=GeomCart(:,2) |
777 | 1 | pfleura2 | XyzGeomF(IGeom,3,:)=GeomCart(:,3) |
778 | 1 | pfleura2 | |
779 | 1 | pfleura2 | Energies(IGeom)=E |
780 | 1 | pfleura2 | |
781 | 1 | pfleura2 | if (debug) THEN |
782 | 1 | pfleura2 | Line='DBG Path, GradTmp' |
783 | 1 | pfleura2 | Call PrintGeom(Line,Nat,NCoord,GRadTmp,Coord,6,Atome,Order,OrderInv,IndZmat) |
784 | 1 | pfleura2 | END IF |
785 | 1 | pfleura2 | Grad(IGeom,:)=GradTmp |
786 | 1 | pfleura2 | END DO ! matches DO IGeom=IGeom0,IGeomF |
787 | 1 | pfleura2 | DEALLOCATE(GradTmp,GeomCart) |
788 | 12 | pfleura2 | |
789 | 1 | pfleura2 | END IF ! matches IF (RunMode=="PARA") THEN AFTER ELSE |
790 | 1 | pfleura2 | |
791 | 12 | pfleura2 | WRITE(IOOUT,*) " Forces all read. Moving on" |
792 | 12 | pfleura2 | |
793 | 12 | pfleura2 | |
794 | 1 | pfleura2 | DEALLOCATE(GeomTmp) |
795 | 1 | pfleura2 | DEALLOCATE(x,y,z) |
796 | 1 | pfleura2 | DEALLOCATE(x_k,y_k,z_k) |
797 | 1 | pfleura2 | if (debug) Call header('Exiting EgradPath') |
798 | 1 | pfleura2 | RETURN |
799 | 1 | pfleura2 | 999 WRITE(*,*) "EgradPath : We should not be here !" |
800 | 12 | pfleura2 | WRITE(*,*) "End of file happened too soon" |
801 | 1 | pfleura2 | STOP |
802 | 12 | pfleura2 | 9991 WRITE(*,*) "EgradPath : We should not be here !" |
803 | 12 | pfleura2 | WRITE(*,*) "EgradPath was Looking for 'EDIFF is reached'" |
804 | 12 | pfleura2 | STOP |
805 | 12 | pfleura2 | 9992 WRITE(*,*) "EgradPath : We should not be here !" |
806 | 12 | pfleura2 | WRITE(*,*) "EgradPath was Looking for 'TOTEN'" |
807 | 12 | pfleura2 | STOP |
808 | 12 | pfleura2 | 9993 WRITE(*,*) "EgradPath : We should not be here !" |
809 | 12 | pfleura2 | WRITE(*,*) "EgradPath was Looking for 'TOTAL-FORCE'" |
810 | 12 | pfleura2 | STOP |
811 | 12 | pfleura2 | 9994 WRITE(*,*) "EgradPath : We should not be here !" |
812 | 12 | pfleura2 | WRITE(*,*) "EgradPath was reading the forces " |
813 | 12 | pfleura2 | STOP |
814 | 1 | pfleura2 | |
815 | 1 | pfleura2 | END SUBROUTINE EGRADPATH |