Statistiques
| Révision :

root / src / egrad_vasp.f90 @ 1

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