root / src / egrad_gaussian.f90
Historique | Voir | Annoter | Télécharger (7,16 ko)
1 |
subroutine egrad_gaussian(e,geomcart,gradcart) |
---|---|
2 |
|
3 |
! This routines calculates the energy and the gradient of |
4 |
! a molecule, using Gaussian |
5 |
|
6 |
!---------------------------------------------------------------------- |
7 |
! Copyright 2003-2014 Ecole Normale Supérieure de Lyon, |
8 |
! Centre National de la Recherche Scientifique, |
9 |
! Université Claude Bernard Lyon 1. All rights reserved. |
10 |
! |
11 |
! This work is registered with the Agency for the Protection of Programs |
12 |
! as IDDN.FR.001.100009.000.S.P.2014.000.30625 |
13 |
! |
14 |
! Authors: P. Fleurat-Lessard, P. Dayal |
15 |
! Contact: optnpath@gmail.com |
16 |
! |
17 |
! This file is part of "Opt'n Path". |
18 |
! |
19 |
! "Opt'n Path" is free software: you can redistribute it and/or modify |
20 |
! it under the terms of the GNU Affero General Public License as |
21 |
! published by the Free Software Foundation, either version 3 of the License, |
22 |
! or (at your option) any later version. |
23 |
! |
24 |
! "Opt'n Path" is distributed in the hope that it will be useful, |
25 |
! but WITHOUT ANY WARRANTY; without even the implied warranty of |
26 |
! |
27 |
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
28 |
! GNU Affero General Public License for more details. |
29 |
! |
30 |
! You should have received a copy of the GNU Affero General Public License |
31 |
! along with "Opt'n Path". If not, see <http://www.gnu.org/licenses/>. |
32 |
! |
33 |
! Contact The Office of Technology Licensing, valorisation@ens-lyon.fr, |
34 |
! for commercial licensing opportunities. |
35 |
!---------------------------------------------------------------------- |
36 |
|
37 |
use Path_module, only : Nat, renum,Order,OrderInv,AtName, Coord,unit,ProgExe |
38 |
use Io_module |
39 |
|
40 |
IMPLICIT NONE |
41 |
|
42 |
! Energy (calculated if F300K=.F., else estimated) |
43 |
REAL(KREAL), INTENT (OUT) :: e |
44 |
! Nb degree of freedom |
45 |
REAL(KREAL), INTENT (IN) :: geomcart(Nat,3) |
46 |
! Gradient calculated at Geom geometry |
47 |
REAL(KREAL), INTENT (OUT) :: gradcart(3*NAt) |
48 |
|
49 |
! ====================================================================== |
50 |
|
51 |
character(LCHARS) :: LINE |
52 |
|
53 |
logical :: debug, G03Ok |
54 |
LOGICAL, SAVE :: first=.true. |
55 |
|
56 |
REAL(KREAL) :: Pi |
57 |
|
58 |
INTEGER(KINT) :: ItG03 |
59 |
INTEGER(KINT) :: iat, i, n3at |
60 |
INTEGER(KINT) :: ITmp1, ITmp2 |
61 |
|
62 |
CHARACTER(132) :: FileIn, FileOut, FileChk |
63 |
|
64 |
CHARACTER(3) :: SepIn=' < ', SepOut=' > ' |
65 |
CHARACTER(VLCHARS), SAVE :: RunCommand |
66 |
! LOGICAL, SAVE :: FCopyRstrt=.False., FOrderChecked=.False. |
67 |
|
68 |
! ====================================================================== |
69 |
|
70 |
LOGICAL, EXTERNAL :: valid |
71 |
|
72 |
! ====================================================================== |
73 |
|
74 |
Pi=dacos(-1.0d0) |
75 |
n3at=3*nat |
76 |
|
77 |
debug=valid('EGRAD') |
78 |
if (debug) Call Header("Entering Egrad_gaussian") |
79 |
|
80 |
if (first) THEN |
81 |
First=.FALSE. |
82 |
call system("uname > Path.tmp") |
83 |
OPEN(IOTMP,FILE='Path.tmp') |
84 |
READ(IOTMP,'(A)') Line |
85 |
CLOSE(IOTMP) |
86 |
IF (index(Line,'CYGWIN').NE.0) THEN |
87 |
SepIn=' ' |
88 |
SepOut=' ' |
89 |
END IF |
90 |
END IF |
91 |
|
92 |
gradcart=0. |
93 |
E=0. |
94 |
FileIn=Trim(CalcName) // "." // Trim(ISuffix) |
95 |
FileOut=Trim(CalcName) // "." // Trim(OSuffix) |
96 |
FileChk=Trim(CalcName) // ".chk" |
97 |
RunCommand=Trim(ProgExe) // SepIn // Trim(FileIn) // SepOut // Trim(FileOut) |
98 |
|
99 |
IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(RunCommand) |
100 |
|
101 |
! we create the input file |
102 |
|
103 |
OPEN(IOTMP,File=FileIn) |
104 |
|
105 |
Current => Gauss_root |
106 |
DO WHILE (ASSOCIATED(Current%next)) |
107 |
WRITE(IOTMP,'(1X,A)') Trim(current%line) |
108 |
WRITE(*,'(1X,A)') Trim(current%line) |
109 |
Current => current%next |
110 |
END DO |
111 |
|
112 |
WRITE(IOTMP,*) |
113 |
|
114 |
Current => Gauss_comment |
115 |
DO WHILE (ASSOCIATED(Current%next)) |
116 |
WRITE(IOTMP,'(1X,A)') Trim(current%line) |
117 |
Current => current%next |
118 |
END DO |
119 |
! WRITE(IOOUT,"('Geometry ',I3,'/',I3,' for iteration ',I3)"),Igeom,NgeomF,Iopt |
120 |
|
121 |
WRITE(IOTMP,*) |
122 |
WRITE(IOTMP,*) Trim(Gauss_charge) |
123 |
|
124 |
DO I=1,Nat |
125 |
If (renum) THEN |
126 |
Iat=Order(I) |
127 |
WRITE(IOTMP,'(1X,A10,3(1X,F15.8),1X,A)') Trim(AtName(I)),GeomCart(Iat,:), TRIM(Gauss_Paste(I)) |
128 |
ELSE |
129 |
Iat=OrderInv(I) |
130 |
WRITE(IOTMP,'(1X,A10,3(1X,F15.8),1X,A)') Trim(AtName(Iat)),GeomCart(I,:), TRIM(Gauss_Paste(Iat)) |
131 |
END IF |
132 |
END DO |
133 |
|
134 |
WRITE(IOTMP,*) |
135 |
|
136 |
Current => Gauss_End |
137 |
DO WHILE (ASSOCIATED(Current%next)) |
138 |
WRITE(IOTMP,'(1X,A)') Trim(current%line) |
139 |
Current => current%next |
140 |
END DO |
141 |
|
142 |
WRITE(IOTMP,*) |
143 |
CLOSE(IOTMP) |
144 |
|
145 |
IF (debug) THEN |
146 |
WRITE(*,*) 'DBG EGRAD ' // Trim(FileIn),' COORD=',TRIM(COORD) |
147 |
call system('cat ' // Trim(FileIn)) |
148 |
END IF |
149 |
|
150 |
!!! PFL 15th March 2008 -> |
151 |
! Sometimes for unknown reasons G03 fails. |
152 |
! We thus add one test to check that Gaussian has finisehd properly. |
153 |
! We also check that we do not run an infinite loop... |
154 |
|
155 |
G03Ok=.FALSE. |
156 |
ItG03=0 |
157 |
|
158 |
DO WHILE ((.NOT.G03Ok).AND.(ItG03<3)) |
159 |
ItG03=ItG03+1 |
160 |
call system(RunCommand) |
161 |
OPEN(IOTMP,FILE=FileOut, STATUS='old', POSITION='APPEND') |
162 |
BACKSPACE(IOTMP) |
163 |
READ(IOTMP,'(A)') Line |
164 |
IF (LINE(2:19).EQ.'Normal termination') G03Ok=.TRUE. |
165 |
END DO |
166 |
CLOSE(IOTMP) |
167 |
IF (.NOT.G03Ok) THEN |
168 |
WRITE(*,*) "Egrad_gaussian Failed..." |
169 |
WRITE(*,*) "Wrong termination of Gaussian after 3 tries. Check Input and geometry" |
170 |
STOP |
171 |
END IF |
172 |
|
173 |
if (debug) WRITE(*,*) 'DBG EGRAD, back from calculation' |
174 |
|
175 |
! Whatever the coordinate system, we read the cartesian gradient |
176 |
! and we convert it to Zmat or other if needed |
177 |
|
178 |
OPEN(IOTMP,FILE=FileOut, STATUS='old') |
179 |
! We first search for the forces |
180 |
READ(IOTMP,'(A)') LINE |
181 |
DO WHILE (INDEX(LINE,'Forces (Hartrees/Bohr)')==0) |
182 |
! here, we have a problem, because there are so many ways to write E for Gaussian :-) |
183 |
! but we do not really need E except if we want to weight the different points |
184 |
! on the path... that will be done later :-) |
185 |
! And I should use ConvGaussXmol ... |
186 |
! PFL 3rd Jan 2008 -> |
187 |
! I have finally upgraded the energy reading phase so that it should be able to read |
188 |
! many formats corresponding to HF, DFT, MP2, AM1, PM3, ONIOM with or w/o PCM |
189 |
|
190 |
! For MP2, energy is after the last = |
191 |
IF ((Line(2:9)=="E2 = ")) THEN |
192 |
Itmp1=Index(LINE,"=",BACK=.TRUE.)+1 |
193 |
READ(LINE(Itmp1:),*) e |
194 |
END IF |
195 |
! For other methods, it is after the first = |
196 |
IF ((LINE(2:9)=='Energy= ').OR.(Line(2:9)=="SCF Done").OR.(Line(2:9)=="ONIOM: e") & |
197 |
.OR.(Line(2:9)==" with al") & |
198 |
.OR.(Line(2:31)=="Counterpoise: corrected energy")) THEN |
199 |
Itmp1=Index(LINE,"=")+1 |
200 |
READ(LINE(Itmp1:),*) e |
201 |
END IF |
202 |
! <- PFL 3 Jan 2008 |
203 |
READ(IOTMP,'(A)') LINE |
204 |
END DO |
205 |
if (debug) WRITE(*,*) "GradCart reading - renum=",Renum |
206 |
READ(IOTMP,'(A)') LINE |
207 |
READ(IOTMP,'(A)') LINE |
208 |
DO I=1,Nat |
209 |
Iat=I |
210 |
IF (renum) Iat=Order(I) |
211 |
READ(IOTMP,*) Itmp1,ITmp2, GradCart(3*Iat-2:3*Iat) |
212 |
if (debug) WRITE(*,*) "I,Iat,GradCart:",I,Iat,GradCart(3*Iat-2:3*Iat) |
213 |
END DO |
214 |
|
215 |
! Gaussian gives the Forces in ua/a0, so we convert it into the |
216 |
! gradient in ua/Angstrom |
217 |
gradCart=-1.0d0*gradCart*Unit |
218 |
|
219 |
! if (debug.AND.(COORD.EQ.'ZMAT')) THEN |
220 |
! DO WHILE (INDEX(LINE,'Internal Coordinate Forces')==0) |
221 |
! READ(IOTMP,'(A)') LINE |
222 |
!! WRITE(*,*) 'Pas bon:' // TRIM(LINE) |
223 |
! END DO |
224 |
! WRITE(*,*) TRIM(LINE) |
225 |
! DO I=1,Nat+2 |
226 |
! READ(IOTMP,'(A)') LINE |
227 |
! WRITE(*,*) TRIM(LINE) |
228 |
! END DO |
229 |
! END IF |
230 |
|
231 |
CLOSE(IOTMP) |
232 |
|
233 |
if (debug) Call header("Egrad_gaussian Over") |
234 |
|
235 |
RETURN |
236 |
|
237 |
!999 CONTINUE |
238 |
! if (.NOT.Ftmp) WRITE(*,*) 'We should not be here !!!!' |
239 |
! STOP |
240 |
! ====================================================================== |
241 |
end subroutine egrad_gaussian |