Statistiques
| Révision :

root / src / egrad_mopac.f90 @ 1

Historique | Voir | Annoter | Télécharger (4,71 ko)

1 1 equemene
subroutine egrad_mopac(e,geomcart,gradcart)
2 1 equemene
3 1 equemene
  ! This routines calculates the energy and the gradient of
4 1 equemene
  ! a molecule, using mopac
5 1 equemene
6 1 equemene
  use Path_module, only : Nat, renum,Order,OrderInv,AtName, Coord, dzdc, indzmat,Nom,Atome, massat, unit,ProgExe, au2kcal
7 1 equemene
  use Io_module
8 1 equemene
9 1 equemene
  IMPLICIT NONE
10 1 equemene
11 1 equemene
  ! Energy (calculated if F300K=.F., else estimated)
12 1 equemene
  REAL(KREAL), INTENT (OUT) :: e
13 1 equemene
  ! Nb degree of freedom
14 1 equemene
  REAL(KREAL), INTENT (IN) :: geomcart(Nat,3)
15 1 equemene
  ! Gradient calculated at Geom geometry
16 1 equemene
  REAL(KREAL), INTENT (OUT) :: gradcart(3*NAt)
17 1 equemene
18 1 equemene
  ! ======================================================================
19 1 equemene
20 1 equemene
  character(LCHARS) :: LINE
21 1 equemene
22 1 equemene
  logical           :: debug
23 1 equemene
  LOGICAL, SAVE :: first=.true.
24 1 equemene
25 1 equemene
  REAL(KREAL), ALLOCATABLE :: GradTmp(:)
26 1 equemene
27 1 equemene
  REAL(KREAL) :: d, a_val, Pi
28 1 equemene
29 1 equemene
  REAL(KREAL) :: coef,x
30 1 equemene
  INTEGER(KINT) :: iat, jat, kat, i, j, k,n3at
31 1 equemene
  INTEGER(KINT) :: ITmp1,ITmp2, Idx
32 1 equemene
  INTEGER(KINT), PARAMETER :: IOLOG=12, IOCOM=11
33 1 equemene
34 1 equemene
  CHARACTER(132) :: FileIn, FileOut, FileChk
35 1 equemene
36 1 equemene
  CHARACTER(3) :: SepIn=' < ', SepOut=' > '
37 1 equemene
  CHARACTER(VLCHARS), SAVE ::  RunCommand
38 1 equemene
39 1 equemene
  ! ======================================================================
40 1 equemene
41 1 equemene
  INTERFACE
42 1 equemene
     function valid(string) result (isValid)
43 1 equemene
       CHARACTER(*), intent(in) :: string
44 1 equemene
       logical                  :: isValid
45 1 equemene
     END function VALID
46 1 equemene
  END INTERFACE
47 1 equemene
48 1 equemene
  ! ======================================================================
49 1 equemene
50 1 equemene
  Pi=dacos(-1.0d0)
51 1 equemene
  n3at=3*nat
52 1 equemene
53 1 equemene
  debug=valid('EGRAD')
54 1 equemene
  if (debug) WRITE(*,*) '================ Entering Egrad_Mopac ===================='
55 1 equemene
56 1 equemene
  if (first) THEN
57 1 equemene
     First=.FALSE.
58 1 equemene
     SepIn='   '
59 1 equemene
     SepOut='   '
60 1 equemene
  END IF
61 1 equemene
62 1 equemene
  gradcart=0.
63 1 equemene
  E=0.
64 1 equemene
65 1 equemene
  FileIn=Trim(CalcName) // "." // Trim(ISuffix)
66 1 equemene
  FileOut=Trim(CalcName) // "." // Trim(OSuffix)
67 1 equemene
  FileChk=Trim(CalcName) // ".arc"
68 1 equemene
69 1 equemene
  RunCommand=Trim(ProgExe) // SepIn // Trim(FileIn)
70 1 equemene
71 1 equemene
  IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(RunCommand)
72 1 equemene
73 1 equemene
  ! we create the input file
74 1 equemene
75 1 equemene
  OPEN(IOTMP,File=FileIn)
76 1 equemene
  WRITE(IOTMP,'(A)') '* All comments are put at the begining of the file'
77 1 equemene
  Current => Mopac_comment
78 1 equemene
  DO WHILE (ASSOCIATED(Current%next))
79 1 equemene
     WRITE(IOTMP,'(A)') Trim(current%line)
80 1 equemene
     Current => current%next
81 1 equemene
  END DO
82 1 equemene
83 1 equemene
  Current => Mopac_root
84 1 equemene
  DO WHILE (ASSOCIATED(Current%next))
85 1 equemene
     WRITE(IOTMP,'(1X,A)') Trim(current%line)
86 1 equemene
     WRITE(*,'(1X,A)') Trim(current%line)
87 1 equemene
     Current => current%next
88 1 equemene
  END DO
89 1 equemene
90 1 equemene
  DO I=1,Nat
91 1 equemene
     If (renum) THEN
92 1 equemene
        Iat=Order(I)
93 1 equemene
        WRITE(IOTMP,'(1X,A10,3(1X,F15.8))') Trim(AtName(I)),GeomCart(Iat,:)
94 1 equemene
     ELSE
95 1 equemene
        Iat=OrderInv(I)
96 1 equemene
        WRITE(IOTMP,'(1X,A10,3(1X,F15.8))') Trim(AtName(Iat)),GeomCart(I,:)
97 1 equemene
     END IF
98 1 equemene
  END DO
99 1 equemene
100 1 equemene
     WRITE(IOTMP,'(1X,A)') Trim(Mopac_EndGeom)
101 1 equemene
102 1 equemene
  Current => Mopac_End
103 1 equemene
  DO WHILE (ASSOCIATED(Current%next))
104 1 equemene
     WRITE(IOTMP,'(1X,A)') Trim(current%line)
105 1 equemene
     Current => current%next
106 1 equemene
  END DO
107 1 equemene
108 1 equemene
  WRITE(IOTMP,*)
109 1 equemene
  CLOSE(IOTMP)
110 1 equemene
111 1 equemene
  IF (debug) THEN
112 1 equemene
     WRITE(*,*) 'DBG EGRAD ' //  Trim(FileIn),' COORD=',TRIM(COORD)
113 1 equemene
     call system('cat ' // Trim(FileIn))
114 1 equemene
  END IF
115 1 equemene
116 1 equemene
117 1 equemene
  call system(RunCommand)
118 1 equemene
119 1 equemene
  if (debug) WRITE(*,*) 'DBG EGRAD, back from calculation'
120 1 equemene
121 1 equemene
  ! Whatever the coordinate system, we read the cartesian gradient
122 1 equemene
  ! and we convert it to Zmat or other if needed
123 1 equemene
124 1 equemene
  OPEN(IOTMP,FILE=FileOut, STATUS='old')
125 1 equemene
  ! We first search for the forces
126 1 equemene
  READ(IOTMP,'(A)') LINE
127 1 equemene
128 1 equemene
  DO WHILE (INDEX(LINE,'FINAL  POINT  AND  DERIVATIVES')==0)
129 1 equemene
     IF (INDEX(Line,'FINAL HEAT OF FORMATION').NE.0)  THEN
130 1 equemene
        Itmp1=Index(Line,"=")+1
131 1 equemene
        ITmp2=Index(Line,"KCAL")-1
132 1 equemene
        Read(Line(Itmp1:Itmp2),*) E
133 1 equemene
        E=E/au2kcal
134 1 equemene
     END IF
135 1 equemene
     READ(IOTMP,'(A)',END=999) LINE
136 1 equemene
  END DO
137 1 equemene
138 1 equemene
  if (debug) WRITE(*,*) "GradCart reading - renum=",Renum
139 1 equemene
  READ(IOTMP,'(A)') LINE
140 1 equemene
  READ(IOTMP,'(A)') LINE
141 1 equemene
  DO I=1,Nat
142 1 equemene
     Iat=I
143 1 equemene
     IF (renum) Iat=Order(I)
144 1 equemene
     DO J=1,3
145 1 equemene
        READ(IOTMP,'(A)') Line
146 1 equemene
! Gradient is on the 7th column
147 1 equemene
        DO K=1,6
148 1 equemene
           Line=AdjustL(Line)
149 1 equemene
           Idx=Index(Line," ")
150 1 equemene
           Line=AdjustL(Line(Idx+1:))
151 1 equemene
        END DO
152 1 equemene
        READ(Line,*) GradCart(3*Iat-3+j)
153 1 equemene
     END DO
154 1 equemene
     if (debug) WRITE(*,*) "I,Iat,GradCart:",I,Iat,GradCart(3*Iat-2:3*Iat)
155 1 equemene
  END DO
156 1 equemene
157 1 equemene
! Mopac gives the Forces in kcal/angstrom, so we convert it into the
158 1 equemene
! gradient in ua/Angstrom
159 1 equemene
  gradCart=gradCart/au2kcal
160 1 equemene
161 1 equemene
!  if (debug.AND.(COORD.EQ.'ZMAT')) THEN
162 1 equemene
!     DO WHILE (INDEX(LINE,'Internal Coordinate Forces')==0)
163 1 equemene
!        READ(IOTMP,'(A)') LINE
164 1 equemene
!!        WRITE(*,*) 'Pas bon:' // TRIM(LINE)
165 1 equemene
!     END DO
166 1 equemene
!     WRITE(*,*) TRIM(LINE)
167 1 equemene
!     DO I=1,Nat+2
168 1 equemene
!        READ(IOTMP,'(A)') LINE
169 1 equemene
!        WRITE(*,*) TRIM(LINE)
170 1 equemene
!     END DO
171 1 equemene
!  END IF
172 1 equemene
173 1 equemene
  CLOSE(IOTMP)
174 1 equemene
175 1 equemene
  if (debug) WRITE(*,*) '================  Egrad_Mopac Over ===================='
176 1 equemene
177 1 equemene
  RETURN
178 1 equemene
179 1 equemene
999 CONTINUE
180 1 equemene
  WRITE(*,*) 'We should not be here !!!!'
181 1 equemene
  STOP
182 1 equemene
  ! ======================================================================
183 1 equemene
end subroutine egrad_mopac