root / src / egrad_mopac.f90
Historique | Voir | Annoter | Télécharger (6,23 ko)
1 |
subroutine egrad_mopac(e,geomcart,gradcart) |
---|---|
2 |
|
3 |
! This routines calculates the energy and the gradient of |
4 |
! a molecule, using mopac |
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, ProgExe, & |
38 |
FPBC, Lat_a, Lat_b, Lat_c, IPer |
39 |
use Io_module |
40 |
|
41 |
IMPLICIT NONE |
42 |
|
43 |
! Energy (calculated if F300K=.F., else estimated) |
44 |
REAL(KREAL), INTENT (OUT) :: e |
45 |
! Nb degree of freedom |
46 |
REAL(KREAL), INTENT (IN) :: geomcart(Nat,3) |
47 |
! Gradient calculated at Geom geometry |
48 |
REAL(KREAL), INTENT (OUT) :: gradcart(3*NAt) |
49 |
|
50 |
! ====================================================================== |
51 |
|
52 |
character(LCHARS) :: LINE |
53 |
|
54 |
logical :: debug |
55 |
LOGICAL, SAVE :: first=.true. |
56 |
|
57 |
REAL(KREAL) :: Pi |
58 |
|
59 |
INTEGER(KINT) :: iat, i, j, k, n3at |
60 |
INTEGER(KINT) :: ITmp1,ITmp2, Idx |
61 |
|
62 |
CHARACTER(132) :: FileIn, FileOut, FileChk |
63 |
|
64 |
CHARACTER(3) :: SepIn=' < ', SepOut=' > ' |
65 |
CHARACTER(VLCHARS), SAVE :: RunCommand |
66 |
|
67 |
! ====================================================================== |
68 |
|
69 |
INTERFACE |
70 |
function valid(string) result (isValid) |
71 |
CHARACTER(*), intent(in) :: string |
72 |
logical :: isValid |
73 |
END function VALID |
74 |
END INTERFACE |
75 |
|
76 |
! ====================================================================== |
77 |
|
78 |
Pi=dacos(-1.0d0) |
79 |
n3at=3*nat |
80 |
|
81 |
debug=valid('EGRAD') |
82 |
if (debug) WRITE(*,*) '================ Entering Egrad_Mopac ====================' |
83 |
|
84 |
if (first) THEN |
85 |
First=.FALSE. |
86 |
SepIn=' ' |
87 |
SepOut=' ' |
88 |
END IF |
89 |
|
90 |
gradcart=0. |
91 |
E=0. |
92 |
|
93 |
FileIn=Trim(CalcName) // "." // Trim(ISuffix) |
94 |
FileOut=Trim(CalcName) // "." // Trim(OSuffix) |
95 |
FileChk=Trim(CalcName) // ".arc" |
96 |
|
97 |
RunCommand=Trim(ProgExe) // SepIn // Trim(FileIn) |
98 |
|
99 |
IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(RunCommand) |
100 |
|
101 |
! we create the input file |
102 |
|
103 |
OPEN(IOTMP,File=FileIn) |
104 |
WRITE(IOTMP,'(A)') '* All comments are put at the begining of the file' |
105 |
Current => Mopac_comment |
106 |
DO WHILE (ASSOCIATED(Current%next)) |
107 |
WRITE(IOTMP,'(A)') Trim(current%line) |
108 |
Current => current%next |
109 |
END DO |
110 |
|
111 |
Current => Mopac_root |
112 |
DO WHILE (ASSOCIATED(Current%next)) |
113 |
WRITE(IOTMP,'(1X,A)') Trim(current%line) |
114 |
WRITE(*,'(1X,A)') Trim(current%line) |
115 |
Current => current%next |
116 |
END DO |
117 |
|
118 |
DO I=1,Nat |
119 |
If (renum) THEN |
120 |
Iat=Order(I) |
121 |
WRITE(IOTMP,'(1X,A10,3(1X,F15.8))') Trim(AtName(I)),GeomCart(Iat,:) |
122 |
ELSE |
123 |
Iat=OrderInv(I) |
124 |
WRITE(IOTMP,'(1X,A10,3(1X,F15.8))') Trim(AtName(Iat)),GeomCart(I,:) |
125 |
END IF |
126 |
END DO |
127 |
! If we have periodic boundaries, we have to print them |
128 |
IF (FPBC) THEN |
129 |
IF (IPER>=1) THEN |
130 |
WRITE(IOTMP,'(1X," Tv ",3(1X,F15.8))') Lat_a(1:3) |
131 |
END IF |
132 |
IF (IPER>=2) THEN |
133 |
WRITE(IOTMP,'(1X," Tv ",3(1X,F15.8))') Lat_b(1:3) |
134 |
END IF |
135 |
IF (IPER>=3) THEN |
136 |
WRITE(IOTMP,'(1X," Tv ",3(1X,F15.8))') Lat_c(1:3) |
137 |
END IF |
138 |
END IF |
139 |
|
140 |
WRITE(IOTMP,'(1X,A)') Trim(Mopac_EndGeom) |
141 |
|
142 |
Current => Mopac_End |
143 |
DO WHILE (ASSOCIATED(Current%next)) |
144 |
WRITE(IOTMP,'(1X,A)') Trim(current%line) |
145 |
Current => current%next |
146 |
END DO |
147 |
|
148 |
WRITE(IOTMP,*) |
149 |
CLOSE(IOTMP) |
150 |
|
151 |
IF (debug) THEN |
152 |
WRITE(*,*) 'DBG EGRAD ' // Trim(FileIn),' COORD=',TRIM(COORD) |
153 |
call system('cat ' // Trim(FileIn)) |
154 |
END IF |
155 |
|
156 |
|
157 |
call system(RunCommand) |
158 |
|
159 |
if (debug) WRITE(*,*) 'DBG EGRAD, back from calculation' |
160 |
|
161 |
! Whatever the coordinate system, we read the cartesian gradient |
162 |
! and we convert it to Zmat or other if needed |
163 |
|
164 |
OPEN(IOTMP,FILE=FileOut, STATUS='old') |
165 |
! We first search for the forces |
166 |
READ(IOTMP,'(A)') LINE |
167 |
|
168 |
DO WHILE (INDEX(LINE,'FINAL POINT AND DERIVATIVES')==0) |
169 |
IF (INDEX(Line,'FINAL HEAT OF FORMATION').NE.0) THEN |
170 |
Itmp1=Index(Line,"=")+1 |
171 |
ITmp2=Index(Line,"KCAL")-1 |
172 |
Read(Line(Itmp1:Itmp2),*) E |
173 |
E=E/au2kcal |
174 |
END IF |
175 |
READ(IOTMP,'(A)',END=999) LINE |
176 |
END DO |
177 |
|
178 |
if (debug) WRITE(*,*) "GradCart reading - renum=",Renum |
179 |
READ(IOTMP,'(A)') LINE |
180 |
READ(IOTMP,'(A)') LINE |
181 |
DO I=1,Nat |
182 |
Iat=I |
183 |
IF (renum) Iat=Order(I) |
184 |
DO J=1,3 |
185 |
READ(IOTMP,'(A)') Line |
186 |
! Gradient is on the 7th column |
187 |
DO K=1,6 |
188 |
Line=AdjustL(Line) |
189 |
Idx=Index(Line," ") |
190 |
Line=AdjustL(Line(Idx+1:)) |
191 |
END DO |
192 |
READ(Line,*) GradCart(3*Iat-3+j) |
193 |
END DO |
194 |
if (debug) WRITE(*,*) "I,Iat,GradCart:",I,Iat,GradCart(3*Iat-2:3*Iat) |
195 |
END DO |
196 |
|
197 |
! Mopac gives the Forces in kcal/angstrom, so we convert it into the |
198 |
! gradient in ua/Angstrom |
199 |
gradCart=gradCart/au2kcal |
200 |
|
201 |
! if (debug.AND.(COORD.EQ.'ZMAT')) THEN |
202 |
! DO WHILE (INDEX(LINE,'Internal Coordinate Forces')==0) |
203 |
! READ(IOTMP,'(A)') LINE |
204 |
!! WRITE(*,*) 'Pas bon:' // TRIM(LINE) |
205 |
! END DO |
206 |
! WRITE(*,*) TRIM(LINE) |
207 |
! DO I=1,Nat+2 |
208 |
! READ(IOTMP,'(A)') LINE |
209 |
! WRITE(*,*) TRIM(LINE) |
210 |
! END DO |
211 |
! END IF |
212 |
|
213 |
CLOSE(IOTMP) |
214 |
|
215 |
if (debug) WRITE(*,*) '================ Egrad_Mopac Over ====================' |
216 |
|
217 |
RETURN |
218 |
|
219 |
999 CONTINUE |
220 |
WRITE(*,*) 'We should not be here !!!!' |
221 |
STOP |
222 |
! ====================================================================== |
223 |
end subroutine egrad_mopac |