Statistiques
| Révision :

root / src / egrad_turbomole.f90 @ 12

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

1 1 pfleura2
subroutine egrad_turbomole(e,geomcart,gradcart)
2 1 pfleura2
3 1 pfleura2
  ! This routines calculates the energy and the gradient of
4 1 pfleura2
  ! a molecule, using TurboMole
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 12 pfleura2
37 1 pfleura2
  use Path_module, only : Nat, AtName, unit,ProgExe,renum,order,OrderInv
38 1 pfleura2
  use Io_module
39 1 pfleura2
40 1 pfleura2
  IMPLICIT NONE
41 1 pfleura2
42 1 pfleura2
  ! Energy (calculated if F300K=.F., else estimated)
43 1 pfleura2
  REAL(KREAL), INTENT (OUT) :: e
44 1 pfleura2
  ! Nb degree of freedom
45 1 pfleura2
  REAL(KREAL), INTENT (IN) :: geomcart(Nat,3)
46 1 pfleura2
  ! Gradient calculated at Geom geometry
47 1 pfleura2
  REAL(KREAL), INTENT (OUT) :: gradcart(3*NAt)
48 1 pfleura2
49 1 pfleura2
  ! ======================================================================
50 1 pfleura2
51 1 pfleura2
  character(LCHARS) :: LINE
52 1 pfleura2
53 1 pfleura2
  logical           :: debug
54 1 pfleura2
55 1 pfleura2
56 2 pfleura2
  REAL(KREAL) :: Pi
57 1 pfleura2
58 2 pfleura2
  INTEGER(KINT) :: iat, i, n3at
59 2 pfleura2
  INTEGER(KINT) :: Idx
60 1 pfleura2
61 2 pfleura2
  CHARACTER(132) :: FileIn, FileOut
62 1 pfleura2
63 2 pfleura2
!  CHARACTER(3) :: SepIn=' ', SepOut=' '
64 1 pfleura2
  CHARACTER(VLCHARS), SAVE ::  RunCommand
65 1 pfleura2
66 1 pfleura2
  ! ======================================================================
67 1 pfleura2
68 1 pfleura2
  INTERFACE
69 1 pfleura2
     function valid(string) result (isValid)
70 1 pfleura2
       CHARACTER(*), intent(in) :: string
71 1 pfleura2
       logical                  :: isValid
72 1 pfleura2
     END function VALID
73 1 pfleura2
  END INTERFACE
74 1 pfleura2
75 1 pfleura2
  ! ======================================================================
76 1 pfleura2
77 1 pfleura2
  Pi=dacos(-1.0d0)
78 1 pfleura2
  n3at=3*nat
79 1 pfleura2
80 1 pfleura2
  debug=valid('EGRAD')
81 1 pfleura2
  if (debug) WRITE(*,*) '================ Entering Egrad_TurboMole ===================='
82 1 pfleura2
83 1 pfleura2
84 1 pfleura2
  gradcart=0.
85 1 pfleura2
  E=0.
86 1 pfleura2
87 1 pfleura2
  FileIn="coord"
88 1 pfleura2
  FileOut="gradient"
89 1 pfleura2
90 1 pfleura2
  RunCommand=Trim(ProgExe)
91 1 pfleura2
92 1 pfleura2
  IF (DEBUG) WRITE(*,*)'RunCommand:',TRIM(RunCommand)
93 1 pfleura2
94 1 pfleura2
  ! we create the input file
95 1 pfleura2
96 1 pfleura2
  OPEN(IOTMP,File=FileIn)
97 1 pfleura2
  WRITE(IOTMP,'(1X,A)') '$coord'
98 1 pfleura2
  DO I=1,Nat
99 1 pfleura2
     If (renum) THEN
100 1 pfleura2
        Iat=Order(I)
101 1 pfleura2
        WRITE(IOTMP,'(1X,3(1X,F15.8),A5)') GeomCart(Iat,:)*Unit,Trim(AtName(Iat))
102 1 pfleura2
     ELSE
103 1 pfleura2
        Iat=OrderInv(I)
104 1 pfleura2
        WRITE(IOTMP,'(1X,3(1X,F15.8),A5)') GeomCart(I,:)*Unit,Trim(AtName(Iat))
105 1 pfleura2
     END IF
106 1 pfleura2
  END DO
107 1 pfleura2
  WRITE(IOTMP,'(1X,A)') '$user-defined bonds'
108 1 pfleura2
  WRITE(IOTMP,'(1X,A)') '$end'
109 1 pfleura2
  CLOSE(IOTMP)
110 1 pfleura2
111 1 pfleura2
  IF (debug) THEN
112 1 pfleura2
     WRITE(*,*) 'DBG EGRAD ' //  Trim(FileIn)
113 1 pfleura2
     call system('cat ' // Trim(FileIn))
114 1 pfleura2
  END IF
115 1 pfleura2
116 1 pfleura2
117 1 pfleura2
  call system(RunCommand)
118 1 pfleura2
119 1 pfleura2
  if (debug) WRITE(*,*) 'DBG EGRAD, back from calculation'
120 1 pfleura2
121 1 pfleura2
  ! Whatever the coordinate system, we read the cartesian gradient
122 1 pfleura2
  ! and we convert it to Zmat or other if needed
123 1 pfleura2
124 1 pfleura2
  OPEN(IOTMP,FILE=FileOut, STATUS='old')
125 1 pfleura2
! We check that gradient file is ok
126 1 pfleura2
  READ(IOTMP,'(A)') LINE
127 1 pfleura2
  Line=AdjustL(Line)
128 1 pfleura2
  if (Line(1:5).NE.'$grad') THEN
129 1 pfleura2
     WRITE(*,*) "Error with Turbomole, gradient file does not start with $grad"
130 1 pfleura2
     STOP
131 1 pfleura2
  END IF
132 1 pfleura2
133 1 pfleura2
! we read the energy
134 1 pfleura2
  READ(IOTMP,'(A)') LINE
135 1 pfleura2
! energy is after the second = sign:
136 1 pfleura2
!  cycle =      1    SCF energy =    -2889.2212497430   |dE/dxyz| =  0.064137
137 1 pfleura2
  Idx=Index(Line,'=')
138 1 pfleura2
  Line=Line(Idx+1:)
139 1 pfleura2
  Idx=Index(Line,'=')
140 1 pfleura2
  Line=Line(Idx+1:)
141 1 pfleura2
 READ(Line,*) E
142 1 pfleura2
143 1 pfleura2
! We skip the coordinates
144 1 pfleura2
 DO I=1,Nat
145 1 pfleura2
    READ(IOTMP,'(A)') LINE
146 1 pfleura2
 END DO
147 1 pfleura2
148 1 pfleura2
! We read the gradient
149 1 pfleura2
  DO I=1,Nat
150 1 pfleura2
     Iat=I
151 1 pfleura2
     IF (renum) Iat=Order(I)
152 1 pfleura2
     READ(IOTMP,*)  GradCart(3*Iat-2:3*Iat)
153 1 pfleura2
     if (debug) WRITE(*,*) "I,Iat,GradCart:",I,Iat,GradCart(3*Iat-2:3*Iat)
154 1 pfleura2
  END DO
155 1 pfleura2
156 1 pfleura2
! TurboMole gives the Forces in au/a0, so we convert it into the
157 1 pfleura2
! gradient in ua/Angstrom
158 1 pfleura2
  gradCart=gradCart*unit
159 1 pfleura2
160 1 pfleura2
  CLOSE(IOTMP)
161 1 pfleura2
162 1 pfleura2
  if (debug) WRITE(*,*) '================  Egrad_TurboMole Over ===================='
163 1 pfleura2
164 1 pfleura2
  RETURN
165 1 pfleura2
166 2 pfleura2
!999 CONTINUE
167 2 pfleura2
!  WRITE(*,*) 'We should not be here !!!!'
168 2 pfleura2
!  STOP
169 1 pfleura2
  ! ======================================================================
170 1 pfleura2
end subroutine egrad_TurboMole