Statistiques
| Révision :

root / src / egrad_vasp.f90 @ 12

Historique | Voir | Annoter | Télécharger (5,19 ko)

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