root / src / egrad_vasp.f90
Historique | Voir | Annoter | Télécharger (5,19 ko)
1 |
subroutine egrad_vasp(e,geomcart,gradcart) |
---|---|
2 |
|
3 |
! This routines calculates the energy and the gradient of |
4 |
! a molecule, using VASP |
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, Coord, & |
38 |
ProgExe,Latr,FrozAtoms |
39 |
use Io_module |
40 |
|
41 |
! |
42 |
IMPLICIT NONE |
43 |
|
44 |
! Energy (calculated if F300K=.F., else estimated) |
45 |
REAL(KREAL), INTENT (OUT) :: e |
46 |
! Nb degree of freedom |
47 |
REAL(KREAL), INTENT (IN) :: geomcart(Nat,3) |
48 |
! Gradient calculated at Geom geometry |
49 |
REAL(KREAL), INTENT (OUT) :: gradcart(3*NAt) |
50 |
|
51 |
! ====================================================================== |
52 |
|
53 |
character(LCHARS) :: LINE |
54 |
|
55 |
REAL(KREAL) :: Pi |
56 |
|
57 |
INTEGER(KINT) :: iat, i, n3at |
58 |
REAL(KREAL) :: RTmp1, RTmp2, RTmp3 |
59 |
|
60 |
! |
61 |
CHARACTER(132), PARAMETER :: FileIn='POSCAR', FileOut='OUTCAR' |
62 |
|
63 |
CHARACTER(VLCHARS), SAVE :: RunCommand |
64 |
|
65 |
! REAL(KREAL), ALLOCATABLE :: GeomTmpC(:), GeomTmp(:) |
66 |
|
67 |
LOGICAL :: Debug |
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 |
Pi=dacos(-1.0d0) |
77 |
debug=valid('EGRAD') |
78 |
if (debug) Call header('Entering Egrad_VASP') |
79 |
|
80 |
n3at=3*nat |
81 |
|
82 |
! ALLOCATE(GeomTmpC(n3at)) |
83 |
|
84 |
|
85 |
gradcart=0. |
86 |
E=0. |
87 |
RunCommand=Trim(ProgExe) |
88 |
|
89 |
IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(RunCommand) |
90 |
|
91 |
! we create the input file |
92 |
|
93 |
! write(*,*) "Coucou 2" |
94 |
|
95 |
OPEN(IOTMP,File=FileIn) |
96 |
|
97 |
! we convert the coordinates into Cartesian coordinates |
98 |
! GeomTmpC=Reshape(GeomCart,SHAPE=(/3*Nat/)) |
99 |
|
100 |
Call PrintGeomVasp(Vasp_Title,geomcart,Latr,FrozAtoms,IoTmp) |
101 |
|
102 |
CLOSE(IOTMP) |
103 |
|
104 |
|
105 |
IF (debug) THEN |
106 |
WRITE(*,*) 'DBG EGRAD ' // Trim(FileIn),' COORD=',TRIM(COORD) |
107 |
call system('cat ' // Trim(FileIn)) |
108 |
END IF |
109 |
|
110 |
call system(RunCommand) |
111 |
|
112 |
if (debug) WRITE(*,*) 'DBG EGRAD_VASP, back from calculation' |
113 |
! if (debug) WRITE(*,*) "DBG EGRAD_VASP Order(I)",Order(1:Nat) |
114 |
|
115 |
! Whatever the coordinate system, we read the cartesian gradient |
116 |
! and we convert it to Zmat or other if needed |
117 |
|
118 |
OPEN(IOTMP,FILE=FileOut, STATUS='old') |
119 |
! We first search for the forces |
120 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
121 |
! We search for the last part of the OUTCAR file, after wavefunction is converged |
122 |
DO WHILE (INDEX(LINE,'EDIFF is reached')==0) |
123 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
124 |
END DO |
125 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
126 |
DO WHILE (INDEX(LINE,'TOTEN')==0) |
127 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
128 |
END DO |
129 |
Line=Line(26:) |
130 |
READ(LINE,*) E |
131 |
if (debug) WRITE(*,'(1X,A,F20.6)') "DBG Egrad_VASP E=",E |
132 |
|
133 |
! We search for the forces |
134 |
DO WHILE (INDEX(LINE,'TOTAL-FORCE')==0) |
135 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
136 |
END DO |
137 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
138 |
DO I=1,Nat |
139 |
Iat=I |
140 |
IF (renum) Iat=Order(I) |
141 |
READ(IOTMP,*,END=999,ERR=999) Rtmp1,RTmp2,RTmp3,GradCart(3*Iat-2:3*Iat) |
142 |
END DO |
143 |
|
144 |
|
145 |
if (debug) THEN |
146 |
WRITE(*,*) "DBG EGRAD_VASP - GradCart=FORCES - original order" |
147 |
DO I=1,Nat |
148 |
Iat=Order(I) |
149 |
WRITE(*,'(1X,I5,3(1X,F20.6))') I,GradCart(3*Iat-2:3*Iat) |
150 |
END DO |
151 |
WRITE(*,*) "DBG EGRAD_VASP - GradCart=FORCES - internal order" |
152 |
DO I=1,Nat |
153 |
Iat=Order(I) |
154 |
WRITE(*,'(3(1X,I5),3(1X,F20.6))') I,Iat,OrderInv(I),GradCart(3*I-2:3*I) |
155 |
END DO |
156 |
END IF |
157 |
|
158 |
! VASP gives the Forces in eV/Angstrom, so we convert it into the |
159 |
! gradient in ua/Angstrom |
160 |
gradCart=-1.0d0*ev2Au*gradCart |
161 |
|
162 |
! if (debug.AND.(COORD.EQ.'ZMAT')) THEN |
163 |
! DO WHILE (INDEX(LINE,'Internal Coordinate Forces')==0) |
164 |
! READ(IOTMP,'(A)') LINE |
165 |
!! WRITE(*,*) 'Pas bon:' // TRIM(LINE) |
166 |
! END DO |
167 |
! WRITE(*,*) TRIM(LINE) |
168 |
! DO I=1,Nat+2 |
169 |
! READ(IOTMP,'(A)') LINE |
170 |
! WRITE(*,*) TRIM(LINE) |
171 |
! END DO |
172 |
! END IF |
173 |
|
174 |
CLOSE(IOTMP) |
175 |
|
176 |
! DeALLOCATE(GeomTmpC) |
177 |
|
178 |
|
179 |
if (debug) Call header('Egrad_VASP Over') |
180 |
|
181 |
RETURN |
182 |
|
183 |
999 CONTINUE |
184 |
WRITE(*,*) 'We should not be here !!!!' |
185 |
STOP |
186 |
! ====================================================================== |
187 |
end subroutine egrad_vasp |
188 |
|
189 |
|
190 |
|
191 |
|