root / src / egrad_gaussian.f90 @ 5
Historique | Voir | Annoter | Télécharger (5,89 ko)
1 | 1 | pfleura2 | subroutine egrad_gaussian(e,geomcart,gradcart) |
---|---|---|---|
2 | 1 | pfleura2 | |
3 | 1 | pfleura2 | ! This routines calculates the energy and the gradient of |
4 | 1 | pfleura2 | ! a molecule, using Gaussian |
5 | 1 | pfleura2 | |
6 | 1 | pfleura2 | use Path_module, only : Nat, renum,Order,OrderInv,AtName, Coord, dzdc, indzmat,Nom,Atome, massat, unit,ProgExe |
7 | 1 | pfleura2 | use Io_module |
8 | 1 | pfleura2 | |
9 | 1 | pfleura2 | IMPLICIT NONE |
10 | 1 | pfleura2 | |
11 | 1 | pfleura2 | ! Energy (calculated if F300K=.F., else estimated) |
12 | 1 | pfleura2 | REAL(KREAL), INTENT (OUT) :: e |
13 | 1 | pfleura2 | ! Nb degree of freedom |
14 | 1 | pfleura2 | REAL(KREAL), INTENT (IN) :: geomcart(Nat,3) |
15 | 1 | pfleura2 | ! Gradient calculated at Geom geometry |
16 | 1 | pfleura2 | REAL(KREAL), INTENT (OUT) :: gradcart(3*NAt) |
17 | 1 | pfleura2 | |
18 | 1 | pfleura2 | ! ====================================================================== |
19 | 1 | pfleura2 | |
20 | 1 | pfleura2 | character(LCHARS) :: LINE |
21 | 1 | pfleura2 | |
22 | 1 | pfleura2 | logical :: debug, G03Ok |
23 | 1 | pfleura2 | LOGICAL, SAVE :: first=.true. |
24 | 1 | pfleura2 | |
25 | 2 | pfleura2 | REAL(KREAL) :: Pi |
26 | 1 | pfleura2 | |
27 | 1 | pfleura2 | INTEGER(KINT) :: ItG03 |
28 | 2 | pfleura2 | INTEGER(KINT) :: iat, i, n3at |
29 | 2 | pfleura2 | INTEGER(KINT) :: ITmp1, ITmp2 |
30 | 1 | pfleura2 | |
31 | 1 | pfleura2 | CHARACTER(132) :: FileIn, FileOut, FileChk |
32 | 1 | pfleura2 | |
33 | 1 | pfleura2 | CHARACTER(3) :: SepIn=' < ', SepOut=' > ' |
34 | 2 | pfleura2 | CHARACTER(VLCHARS), SAVE :: RunCommand |
35 | 2 | pfleura2 | ! LOGICAL, SAVE :: FCopyRstrt=.False., FOrderChecked=.False. |
36 | 1 | pfleura2 | |
37 | 1 | pfleura2 | ! ====================================================================== |
38 | 1 | pfleura2 | |
39 | 1 | pfleura2 | LOGICAL, EXTERNAL :: valid |
40 | 1 | pfleura2 | |
41 | 1 | pfleura2 | ! ====================================================================== |
42 | 1 | pfleura2 | |
43 | 1 | pfleura2 | Pi=dacos(-1.0d0) |
44 | 1 | pfleura2 | n3at=3*nat |
45 | 1 | pfleura2 | |
46 | 1 | pfleura2 | debug=valid('EGRAD') |
47 | 1 | pfleura2 | if (debug) Call Header("Entering Egrad_gaussian") |
48 | 1 | pfleura2 | |
49 | 1 | pfleura2 | if (first) THEN |
50 | 1 | pfleura2 | First=.FALSE. |
51 | 1 | pfleura2 | call system("uname > Path.tmp") |
52 | 1 | pfleura2 | OPEN(IOTMP,FILE='Path.tmp') |
53 | 1 | pfleura2 | READ(IOTMP,'(A)') Line |
54 | 1 | pfleura2 | CLOSE(IOTMP) |
55 | 1 | pfleura2 | IF (index(Line,'CYGWIN').NE.0) THEN |
56 | 1 | pfleura2 | SepIn=' ' |
57 | 1 | pfleura2 | SepOut=' ' |
58 | 1 | pfleura2 | END IF |
59 | 1 | pfleura2 | END IF |
60 | 1 | pfleura2 | |
61 | 1 | pfleura2 | gradcart=0. |
62 | 1 | pfleura2 | E=0. |
63 | 1 | pfleura2 | FileIn=Trim(CalcName) // "." // Trim(ISuffix) |
64 | 1 | pfleura2 | FileOut=Trim(CalcName) // "." // Trim(OSuffix) |
65 | 1 | pfleura2 | FileChk=Trim(CalcName) // ".chk" |
66 | 1 | pfleura2 | RunCommand=Trim(ProgExe) // SepIn // Trim(FileIn) // SepOut // Trim(FileOut) |
67 | 1 | pfleura2 | |
68 | 1 | pfleura2 | IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(RunCommand) |
69 | 1 | pfleura2 | |
70 | 1 | pfleura2 | ! we create the input file |
71 | 1 | pfleura2 | |
72 | 1 | pfleura2 | OPEN(IOTMP,File=FileIn) |
73 | 1 | pfleura2 | |
74 | 1 | pfleura2 | Current => Gauss_root |
75 | 1 | pfleura2 | DO WHILE (ASSOCIATED(Current%next)) |
76 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A)') Trim(current%line) |
77 | 1 | pfleura2 | WRITE(*,'(1X,A)') Trim(current%line) |
78 | 1 | pfleura2 | Current => current%next |
79 | 1 | pfleura2 | END DO |
80 | 1 | pfleura2 | |
81 | 1 | pfleura2 | WRITE(IOTMP,*) |
82 | 1 | pfleura2 | |
83 | 1 | pfleura2 | Current => Gauss_comment |
84 | 1 | pfleura2 | DO WHILE (ASSOCIATED(Current%next)) |
85 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A)') Trim(current%line) |
86 | 1 | pfleura2 | Current => current%next |
87 | 1 | pfleura2 | END DO |
88 | 1 | pfleura2 | ! WRITE(IOOUT,"('Geometry ',I3,'/',I3,' for iteration ',I3)"),Igeom,NgeomF,Iopt |
89 | 1 | pfleura2 | |
90 | 1 | pfleura2 | WRITE(IOTMP,*) |
91 | 1 | pfleura2 | WRITE(IOTMP,*) Trim(Gauss_charge) |
92 | 1 | pfleura2 | |
93 | 1 | pfleura2 | DO I=1,Nat |
94 | 1 | pfleura2 | If (renum) THEN |
95 | 1 | pfleura2 | Iat=Order(I) |
96 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A10,3(1X,F15.8),A)') Trim(AtName(I)),GeomCart(Iat,:), TRIM(Gauss_Paste(I)) |
97 | 1 | pfleura2 | ELSE |
98 | 1 | pfleura2 | Iat=OrderInv(I) |
99 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A10,3(1X,F15.8),A)') Trim(AtName(Iat)),GeomCart(I,:), TRIM(Gauss_Paste(Iat)) |
100 | 1 | pfleura2 | END IF |
101 | 1 | pfleura2 | END DO |
102 | 1 | pfleura2 | |
103 | 1 | pfleura2 | WRITE(IOTMP,*) |
104 | 1 | pfleura2 | |
105 | 1 | pfleura2 | Current => Gauss_End |
106 | 1 | pfleura2 | DO WHILE (ASSOCIATED(Current%next)) |
107 | 1 | pfleura2 | WRITE(IOTMP,'(1X,A)') Trim(current%line) |
108 | 1 | pfleura2 | Current => current%next |
109 | 1 | pfleura2 | END DO |
110 | 1 | pfleura2 | |
111 | 1 | pfleura2 | WRITE(IOTMP,*) |
112 | 1 | pfleura2 | CLOSE(IOTMP) |
113 | 1 | pfleura2 | |
114 | 1 | pfleura2 | IF (debug) THEN |
115 | 1 | pfleura2 | WRITE(*,*) 'DBG EGRAD ' // Trim(FileIn),' COORD=',TRIM(COORD) |
116 | 1 | pfleura2 | call system('cat ' // Trim(FileIn)) |
117 | 1 | pfleura2 | END IF |
118 | 1 | pfleura2 | |
119 | 1 | pfleura2 | !!! PFL 15th March 2008 -> |
120 | 1 | pfleura2 | ! Sometimes for unknown reasons G03 fails. |
121 | 1 | pfleura2 | ! We thus add one test to check that Gaussian has finisehd properly. |
122 | 1 | pfleura2 | ! We also check that we do not run an infinite loop... |
123 | 1 | pfleura2 | |
124 | 1 | pfleura2 | G03Ok=.FALSE. |
125 | 1 | pfleura2 | ItG03=0 |
126 | 1 | pfleura2 | |
127 | 1 | pfleura2 | DO WHILE ((.NOT.G03Ok).AND.(ItG03<3)) |
128 | 1 | pfleura2 | ItG03=ItG03+1 |
129 | 1 | pfleura2 | call system(RunCommand) |
130 | 1 | pfleura2 | OPEN(IOTMP,FILE=FileOut, STATUS='old', POSITION='APPEND') |
131 | 1 | pfleura2 | BACKSPACE(IOTMP) |
132 | 1 | pfleura2 | READ(IOTMP,'(A)') Line |
133 | 1 | pfleura2 | IF (LINE(2:19).EQ.'Normal termination') G03Ok=.TRUE. |
134 | 1 | pfleura2 | END DO |
135 | 1 | pfleura2 | CLOSE(IOTMP) |
136 | 1 | pfleura2 | IF (.NOT.G03Ok) THEN |
137 | 1 | pfleura2 | WRITE(*,*) "Egrad_gaussian Failed..." |
138 | 1 | pfleura2 | WRITE(*,*) "Wrong termination of Gaussian after 3 tries. Check Input and geometry" |
139 | 1 | pfleura2 | STOP |
140 | 1 | pfleura2 | END IF |
141 | 1 | pfleura2 | |
142 | 1 | pfleura2 | if (debug) WRITE(*,*) 'DBG EGRAD, back from calculation' |
143 | 1 | pfleura2 | |
144 | 1 | pfleura2 | ! Whatever the coordinate system, we read the cartesian gradient |
145 | 1 | pfleura2 | ! and we convert it to Zmat or other if needed |
146 | 1 | pfleura2 | |
147 | 1 | pfleura2 | OPEN(IOTMP,FILE=FileOut, STATUS='old') |
148 | 1 | pfleura2 | ! We first search for the forces |
149 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
150 | 1 | pfleura2 | DO WHILE (INDEX(LINE,'Forces (Hartrees/Bohr)')==0) |
151 | 1 | pfleura2 | ! here, we have a problem, because there are so many ways to write E for Gaussian :-) |
152 | 1 | pfleura2 | ! but we do not really need E except if we want to weight the different points |
153 | 1 | pfleura2 | ! on the path... that will be done later :-) |
154 | 1 | pfleura2 | ! And I should use ConvGaussXmol ... |
155 | 1 | pfleura2 | ! PFL 3rd Jan 2008 -> |
156 | 1 | pfleura2 | ! I have finally upgraded the energy reading phase so that it should be able to read |
157 | 1 | pfleura2 | ! many formats corresponding to HF, DFT, MP2, AM1, PM3, ONIOM with or w/o PCM |
158 | 1 | pfleura2 | |
159 | 1 | pfleura2 | ! For MP2, energy is after the last = |
160 | 1 | pfleura2 | IF ((Line(2:9)=="E2 = ")) THEN |
161 | 1 | pfleura2 | Itmp1=Index(LINE,"=",BACK=.TRUE.)+1 |
162 | 1 | pfleura2 | READ(LINE(Itmp1:),*) e |
163 | 1 | pfleura2 | END IF |
164 | 1 | pfleura2 | ! For other methods, it is after the first = |
165 | 1 | pfleura2 | IF ((LINE(2:9)=='Energy= ').OR.(Line(2:9)=="SCF Done").OR.(Line(2:9)=="ONIOM: e") & |
166 | 1 | pfleura2 | .OR.(Line(2:9)==" with al") & |
167 | 1 | pfleura2 | .OR.(Line(2:31)=="Counterpoise: corrected energy")) THEN |
168 | 1 | pfleura2 | Itmp1=Index(LINE,"=")+1 |
169 | 1 | pfleura2 | READ(LINE(Itmp1:),*) e |
170 | 1 | pfleura2 | END IF |
171 | 1 | pfleura2 | ! <- PFL 3 Jan 2008 |
172 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
173 | 1 | pfleura2 | END DO |
174 | 1 | pfleura2 | if (debug) WRITE(*,*) "GradCart reading - renum=",Renum |
175 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
176 | 1 | pfleura2 | READ(IOTMP,'(A)') LINE |
177 | 1 | pfleura2 | DO I=1,Nat |
178 | 1 | pfleura2 | Iat=I |
179 | 1 | pfleura2 | IF (renum) Iat=Order(I) |
180 | 1 | pfleura2 | READ(IOTMP,*) Itmp1,ITmp2, GradCart(3*Iat-2:3*Iat) |
181 | 1 | pfleura2 | if (debug) WRITE(*,*) "I,Iat,GradCart:",I,Iat,GradCart(3*Iat-2:3*Iat) |
182 | 1 | pfleura2 | END DO |
183 | 1 | pfleura2 | |
184 | 1 | pfleura2 | ! Gaussian gives the Forces in ua/a0, so we convert it into the |
185 | 1 | pfleura2 | ! gradient in ua/Angstrom |
186 | 1 | pfleura2 | gradCart=-1.0d0*gradCart*Unit |
187 | 1 | pfleura2 | |
188 | 1 | pfleura2 | ! if (debug.AND.(COORD.EQ.'ZMAT')) THEN |
189 | 1 | pfleura2 | ! DO WHILE (INDEX(LINE,'Internal Coordinate Forces')==0) |
190 | 1 | pfleura2 | ! READ(IOTMP,'(A)') LINE |
191 | 1 | pfleura2 | !! WRITE(*,*) 'Pas bon:' // TRIM(LINE) |
192 | 1 | pfleura2 | ! END DO |
193 | 1 | pfleura2 | ! WRITE(*,*) TRIM(LINE) |
194 | 1 | pfleura2 | ! DO I=1,Nat+2 |
195 | 1 | pfleura2 | ! READ(IOTMP,'(A)') LINE |
196 | 1 | pfleura2 | ! WRITE(*,*) TRIM(LINE) |
197 | 1 | pfleura2 | ! END DO |
198 | 1 | pfleura2 | ! END IF |
199 | 1 | pfleura2 | |
200 | 1 | pfleura2 | CLOSE(IOTMP) |
201 | 1 | pfleura2 | |
202 | 1 | pfleura2 | if (debug) Call header("Egrad_gaussian Over") |
203 | 1 | pfleura2 | |
204 | 1 | pfleura2 | RETURN |
205 | 1 | pfleura2 | |
206 | 2 | pfleura2 | !999 CONTINUE |
207 | 2 | pfleura2 | ! if (.NOT.Ftmp) WRITE(*,*) 'We should not be here !!!!' |
208 | 2 | pfleura2 | ! STOP |
209 | 1 | pfleura2 | ! ====================================================================== |
210 | 1 | pfleura2 | end subroutine egrad_gaussian |