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