root / src / egrad_vasp.f90 @ 3
Historique | Voir | Annoter | Télécharger (4,11 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 |
use Path_module, only : Nat, renum,Order,OrderInv,AtName, Coord, dzdc, & |
8 |
indzmat,Nom,Atome, massat, unit,Nfroz,Frozen,ProgExe,Latr,FrozAtoms |
9 |
use Io_module |
10 |
|
11 |
! |
12 |
IMPLICIT NONE |
13 |
|
14 |
! Energy (calculated if F300K=.F., else estimated) |
15 |
REAL(KREAL), INTENT (OUT) :: e |
16 |
! Nb degree of freedom |
17 |
REAL(KREAL), INTENT (IN) :: geomcart(Nat,3) |
18 |
! Gradient calculated at Geom geometry |
19 |
REAL(KREAL), INTENT (OUT) :: gradcart(3*NAt) |
20 |
|
21 |
! ====================================================================== |
22 |
|
23 |
character(LCHARS) :: LINE |
24 |
|
25 |
REAL(KREAL) :: d, a_val, Pi |
26 |
|
27 |
INTEGER(KINT) :: iat, jat, kat, i, j, n3at, Nfr |
28 |
REAL(KREAL) :: RTmp1, RTmp2, RTmp3 |
29 |
|
30 |
! |
31 |
CHARACTER(132), PARAMETER :: FileIn='POSCAR', FileOut='OUTCAR' |
32 |
|
33 |
CHARACTER(VLCHARS), SAVE :: RunCommand |
34 |
INTEGER(kint) :: NbLists,LastIt,Nat4,NStep, Firstit |
35 |
|
36 |
INTEGER(KINT) :: IGeom, Idx |
37 |
! REAL(KREAL), ALLOCATABLE :: GeomTmpC(:), GeomTmp(:) |
38 |
REAL(KREAL) :: MRot(3,3), Rmsd |
39 |
|
40 |
REAL(KREAL) :: X,Y,Z |
41 |
LOGICAL :: Debug |
42 |
|
43 |
INTERFACE |
44 |
function valid(string) result (isValid) |
45 |
CHARACTER(*), intent(in) :: string |
46 |
logical :: isValid |
47 |
END function VALID |
48 |
END INTERFACE |
49 |
|
50 |
Pi=dacos(-1.0d0) |
51 |
debug=valid('EGRAD') |
52 |
if (debug) Call header('Entering Egrad_VASP') |
53 |
|
54 |
n3at=3*nat |
55 |
|
56 |
! ALLOCATE(GeomTmpC(n3at)) |
57 |
|
58 |
|
59 |
gradcart=0. |
60 |
E=0. |
61 |
RunCommand=Trim(ProgExe) |
62 |
|
63 |
IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(RunCommand) |
64 |
|
65 |
! we create the input file |
66 |
|
67 |
! write(*,*) "Coucou 2" |
68 |
|
69 |
OPEN(IOTMP,File=FileIn) |
70 |
|
71 |
! we convert the coordinates into Cartesian coordinates |
72 |
! GeomTmpC=Reshape(GeomCart,SHAPE=(/3*Nat/)) |
73 |
|
74 |
Call PrintGeomVasp(Vasp_Title,geomcart,Latr,FrozAtoms,IoTmp) |
75 |
|
76 |
CLOSE(IOTMP) |
77 |
|
78 |
|
79 |
IF (debug) THEN |
80 |
WRITE(*,*) 'DBG EGRAD ' // Trim(FileIn),' COORD=',TRIM(COORD) |
81 |
call system('cat ' // Trim(FileIn)) |
82 |
END IF |
83 |
|
84 |
call system(RunCommand) |
85 |
|
86 |
if (debug) WRITE(*,*) 'DBG EGRAD_VASP, back from calculation' |
87 |
! if (debug) WRITE(*,*) "DBG EGRAD_VASP Order(I)",Order(1:Nat) |
88 |
|
89 |
! Whatever the coordinate system, we read the cartesian gradient |
90 |
! and we convert it to Zmat or other if needed |
91 |
|
92 |
OPEN(IOTMP,FILE=FileOut, STATUS='old') |
93 |
! We first search for the forces |
94 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
95 |
! We search for the last part of the OUTCAR file, after wavefunction is converged |
96 |
DO WHILE (INDEX(LINE,'EDIFF is reached')==0) |
97 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
98 |
END DO |
99 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
100 |
DO WHILE (INDEX(LINE,'TOTEN')==0) |
101 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
102 |
END DO |
103 |
Line=Line(26:) |
104 |
READ(LINE,*) E |
105 |
if (debug) WRITE(*,'(1X,A,F20.6)') "DBG Egrad_VASP E=",E |
106 |
|
107 |
! We search for the forces |
108 |
DO WHILE (INDEX(LINE,'TOTAL-FORCE')==0) |
109 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
110 |
END DO |
111 |
READ(IOTMP,'(A)',END=999,ERR=999) LINE |
112 |
DO I=1,Nat |
113 |
Iat=I |
114 |
IF (renum) Iat=Order(I) |
115 |
READ(IOTMP,*,END=999,ERR=999) Rtmp1,RTmp2,RTmp3,GradCart(3*Iat-2:3*Iat) |
116 |
END DO |
117 |
|
118 |
|
119 |
if (debug) THEN |
120 |
WRITE(*,*) "DBG EGRAD_VASP - GradCart=FORCES - original order" |
121 |
DO I=1,Nat |
122 |
Iat=Order(I) |
123 |
WRITE(*,'(1X,I5,3(1X,F20.6))') I,GradCart(3*Iat-2:3*Iat) |
124 |
END DO |
125 |
WRITE(*,*) "DBG EGRAD_VASP - GradCart=FORCES - internal order" |
126 |
DO I=1,Nat |
127 |
Iat=Order(I) |
128 |
WRITE(*,'(3(1X,I5),3(1X,F20.6))') I,Iat,OrderInv(I),GradCart(3*I-2:3*I) |
129 |
END DO |
130 |
END IF |
131 |
|
132 |
! VASP gives the Forces in eV/Angstrom, so we convert it into the |
133 |
! gradient in ua/Angstrom |
134 |
gradCart=-1.0d0*ev2Au*gradCart |
135 |
|
136 |
! if (debug.AND.(COORD.EQ.'ZMAT')) THEN |
137 |
! DO WHILE (INDEX(LINE,'Internal Coordinate Forces')==0) |
138 |
! READ(IOTMP,'(A)') LINE |
139 |
!! WRITE(*,*) 'Pas bon:' // TRIM(LINE) |
140 |
! END DO |
141 |
! WRITE(*,*) TRIM(LINE) |
142 |
! DO I=1,Nat+2 |
143 |
! READ(IOTMP,'(A)') LINE |
144 |
! WRITE(*,*) TRIM(LINE) |
145 |
! END DO |
146 |
! END IF |
147 |
|
148 |
CLOSE(IOTMP) |
149 |
|
150 |
! DeALLOCATE(GeomTmpC) |
151 |
|
152 |
|
153 |
if (debug) Call header('Egrad_VASP Over') |
154 |
|
155 |
RETURN |
156 |
|
157 |
999 CONTINUE |
158 |
WRITE(*,*) 'We should not be here !!!!' |
159 |
STOP |
160 |
! ====================================================================== |
161 |
end subroutine egrad_vasp |
162 |
|
163 |
|
164 |
|
165 |
|