root / src / Egrad.f90 @ 7
Historique | Voir | Annoter | Télécharger (7,4 ko)
1 |
subroutine Egrad(E,Geom,Grad,NCoord,IGeom,IOpt,GeomCart) |
---|---|
2 |
|
3 |
! This routines calculates the energy E and the gradient Grad of |
4 |
! a molecule with Geometry Geom (may be in internal coordinates), |
5 |
! using for now, either gaussian or Ext, more general later. |
6 |
|
7 |
use Path_module, only : Nat,AtName,Coord,dzdc,indzmat,Nom,Atome,massat,unit, & |
8 |
prog,NCart,XyzGeomF,IntCoordF,XyzGeomI, renum, OrderInv |
9 |
! IntCoordF(NGeomF,NCoord),GeomOld_all(NGeomF,NCoord) |
10 |
! allocated in Path.f90 |
11 |
use Io_module |
12 |
IMPLICIT NONE |
13 |
|
14 |
! Energy (calculated if F300K=.F., else estimated) |
15 |
REAL(KREAL), INTENT (OUT) :: E |
16 |
! NCoord: Number of the degrees of freedom |
17 |
! IGeom: index of the geometry. |
18 |
INTEGER(KINT), INTENT (IN) :: NCoord, IGeom, IOpt |
19 |
! Geometry at which gradient is calculated (cf Factual also): |
20 |
REAL(KREAL), INTENT (INOUT) :: Geom(NCoord) |
21 |
! Gradient calculated at Geom geometry: |
22 |
REAL(KREAL), INTENT (OUT) :: Grad(NCoord) |
23 |
! Cartesian geometry corresponding to (Internal Geometry) Geom: |
24 |
REAL(KREAL), INTENT (OUT) :: GeomCart(Nat,3) |
25 |
|
26 |
! ====================================================================== |
27 |
|
28 |
CHARACTER(LCHARS) :: LINE |
29 |
LOGICAL :: debug |
30 |
LOGICAL, SAVE :: first=.true. |
31 |
REAL(KREAL), SAVE :: Eold=1.e6 |
32 |
REAL(KREAL), ALLOCATABLE :: GradTmp(:) |
33 |
REAL(KREAL), ALLOCATABLE :: GradCart(:) |
34 |
REAL(KREAL), ALLOCATABLE :: x(:), y(:), z(:) |
35 |
REAL(KREAL) :: d, a_val, Pi |
36 |
|
37 |
INTEGER(KINT) :: iat,jat,kat,i,j,IBeg,Idx |
38 |
INTEGER(KINT), PARAMETER :: IOLOG=12, IOCOM=11 |
39 |
|
40 |
CHARACTER(LCHARS) :: ListName, CH32SVAR1 |
41 |
CHARACTER(VLCHARS), SAVE :: RstrtCopy, RunCommand |
42 |
LOGICAL, SAVE :: FCopyRstrt=.False., FOrderChecked=.False. |
43 |
LOGICAL :: FTmp |
44 |
INTEGER(kint) :: NbLists,NStep |
45 |
INTEGER(KINT), PARAMETER :: NbExtName=4 |
46 |
|
47 |
|
48 |
! ====================================================================== |
49 |
INTERFACE |
50 |
function valid(string) result (isValid) |
51 |
CHARACTER(*), intent(in) :: string |
52 |
logical :: isValid |
53 |
END function VALID |
54 |
END INTERFACE |
55 |
! ====================================================================== |
56 |
|
57 |
INTEGER(KINT) :: OptGeom |
58 |
Namelist/path/OptGeom |
59 |
|
60 |
Pi=dacos(-1.0d0) |
61 |
|
62 |
debug=valid('EGRAD') |
63 |
if (debug) Call header("Entering Egrad") |
64 |
|
65 |
Print *, 'EGrad.f90, L73, IOpt=', IOpt |
66 |
|
67 |
Grad=0. |
68 |
E=0. |
69 |
|
70 |
ALLOCATE(GradCart(3*Nat)) |
71 |
ALLOCATE(x(Nat),y(Nat),z(Nat)) |
72 |
SELECT CASE (COORD) |
73 |
CASE ('ZMAT') |
74 |
! In order to avoid problem with numbering and linear angles/diehedral, we convert the |
75 |
! zmat into cartesian coordinates |
76 |
! A remplacer par Int2Cart :-) |
77 |
Call Int2Cart(nat,indzmat,Geom,GeomCart) |
78 |
CASE ('BAKER') |
79 |
WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad_baker.f90, instead of Egrad.f90. Stop" |
80 |
CASE ('MIXED') |
81 |
! write(*,*) "Coucou 4-mixed" |
82 |
Call Mixed2Cart(Nat,indzmat,geom,GeomCart) |
83 |
CASE ('CART') |
84 |
! write(*,*) "Coucou 4-cart" |
85 |
if (debug) WRITE(*,*) "Coord=CART... in egrad" |
86 |
GeomCart=reshape(Geom,(/Nat,3/)) |
87 |
DO Iat=1,Nat |
88 |
WRITE(IOTMP,'(1X,A10,3(1X,F15.8),A)') Trim(AtName(Iat)),Geom(3*Iat-2:3*Iat) |
89 |
END DO |
90 |
CASE ('HYBRID') |
91 |
! write(*,*) "Coucou 4-hybrid" |
92 |
GeomCart=reshape(Geom,(/Nat,3/)) |
93 |
CASE DEFAULT |
94 |
WRITE (*,*) "Coord=",TRIM(COORD)," not recognized. Stop" |
95 |
STOP |
96 |
END SELECT |
97 |
|
98 |
SELECT CASE (Prog) |
99 |
CASE ('GAUSSIAN') |
100 |
! we call the Gaussian routine. |
101 |
! it will return the energy and the gradient in cartesian coordinates. |
102 |
! Gradient calculated at GeomCart geometry and stored in GradCart, which has INTENT(OUT) |
103 |
Call egrad_gaussian(E,GeomCart,GradCart) |
104 |
CASE ('MOPAC') |
105 |
! we call the Mopac routine. |
106 |
! it will return the energy and the gradient in cartesian coordinates. |
107 |
! Gradient calculated at GeomCart geometry and stored in GradCart, which has INTENT(OUT) |
108 |
Call egrad_mopac(E,GeomCart,GradCart) |
109 |
CASE ('VASP') |
110 |
Call egrad_vasp(E,Geomcart,GradCart) |
111 |
CASE ('TURBOMOLE') |
112 |
Call egrad_turbomole(E,Geomcart,GradCart) |
113 |
CASE ('EXT') |
114 |
Call egrad_ext(E,Geomcart,GradCart) |
115 |
CASE ('TEST') |
116 |
Call egrad_test(Nat,E,Geomcart,GradCart) |
117 |
CASE ('CHAMFRE') |
118 |
Call egrad_chamfre(Nat,E,Geomcart,GradCart) |
119 |
END SELECT |
120 |
if (debug) THEN |
121 |
WRITE(*,*) 'DBG EGRAD, GradCart read' |
122 |
DO I=1,Nat |
123 |
Iat=I |
124 |
if (renum) Iat=orderInv(I) |
125 |
WRITE(*,'(1X,A10,3(1X,F15.8),A)') Trim(AtName(Iat)),GradCart(3*I-2:3*I) |
126 |
END DO |
127 |
END IF |
128 |
|
129 |
! If (PROG=="VASP") STOP |
130 |
|
131 |
! We have the gradient in Cartesian coordinates... which is ok for COORD=Cart, Hybrid |
132 |
! but we have to convert it into Zmat if COORD=Zmat |
133 |
SELECT CASE (COORD) |
134 |
CASE ("ZMAT") |
135 |
ALLOCATE(GradTmp(3*Nat)) |
136 |
if (debug) WRITE(*,*) "DBG EGRAD Calling CalcBMat_int" |
137 |
CALL CalcBMat_int (nat, reshape(GeomCart,(/3,Nat/),ORDER=(/2,1/)), indzmat, dzdc, massat,atome) |
138 |
|
139 |
if (debug) WRITE(*,*) "DBG EGRAD Calling ConvGrad_Cart2Int" |
140 |
CALL ConvGrad_Cart2Int(Nat,dzdc,IndZmat,GradCart,GradTmp) |
141 |
|
142 |
if (debug) WRITE(*,*) "DBG EGRAD Storing Grad" |
143 |
Grad(1)=GradTmp(4) |
144 |
Grad(2)=GradTmp(7) |
145 |
! We might have troubles whith angles that are not in the [0:pi] range because, |
146 |
! when converted to catesian coord, they do appear in the [0:Pi] range to gaussian... |
147 |
! so that the derivative is not good, and a multiplicative factor should be added... |
148 |
! |
149 |
! This in fact should be taken care of in B mat calculation... |
150 |
! |
151 |
IF ((Geom(3).LE.0).OR.(Geom(3).GE.Pi)) THEN |
152 |
Grad(3)=-1.0d0*GradTmp(8) |
153 |
ELSE |
154 |
Grad(3)=GradTmp(8) |
155 |
END IF |
156 |
Idx=4 |
157 |
DO I=4,Nat |
158 |
Grad(Idx)=GradTmp(3*i-2) |
159 |
Grad(Idx+1)=GradTmp(3*i-1) |
160 |
Grad(Idx+2)=GradTmp(3*i) |
161 |
Idx=Idx+3 |
162 |
END DO |
163 |
DEALLOCATE(GradTmp) |
164 |
CASE ('BAKER') |
165 |
WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad_baker.f90, instead of Egrad.f90. Stop" |
166 |
CASE ('MIXED') |
167 |
ALLOCATE(GradTmp(3*Nat)) |
168 |
if (debug) WRITE(*,*) "DBG EGRAD Calling CalcBMat_mixed" |
169 |
CALL CalcBMat_mixed (nat, reshape(GeomCart,(/3,Nat/),ORDER=(/2,1/)), indzmat, dzdc, massat,atome) |
170 |
|
171 |
if (Debug) THEN |
172 |
WRITE(*,*) "DzDc" |
173 |
DO I=1,Nat |
174 |
DO J=1,3 |
175 |
WRITE(*,'(50(1X,F12.6))') dzdc(J,I,:,:) |
176 |
END DO |
177 |
END DO |
178 |
END IF |
179 |
|
180 |
|
181 |
CALL ConvGrad_Cart2Int(Nat,dzdc,IndZmat,GradCart,GradTmp) |
182 |
DO I=1,Nat |
183 |
WRITE(*,'(1X,3(1X,F10.5))') GradTmp(3*I-2:3*I) |
184 |
END DO |
185 |
SELECT CASE (NCART) |
186 |
CASE (1) |
187 |
Grad(1:3)=GradTmp(1:3) |
188 |
Grad(4)=GradTmp(4) |
189 |
Grad(5)=GradTmp(7) |
190 |
Grad(6)=GradTmp(8) |
191 |
Idx=7 |
192 |
IBeg=4 |
193 |
CASE(2) |
194 |
Grad(1:3)=GradTmp(1:3) |
195 |
Grad(4:6)=GradTmp(4:6) |
196 |
Grad(7)=GradTmp(7) |
197 |
Grad(8)=GradTmp(8) |
198 |
Idx=9 |
199 |
IBeg=4 |
200 |
CASE DEFAULT |
201 |
Idx=1 |
202 |
IBeg=1 |
203 |
END SELECT |
204 |
DO I=IBeg,Nat |
205 |
Grad(Idx)=GradTmp(3*i-2) |
206 |
Grad(Idx+1)=GradTmp(3*i-1) |
207 |
Grad(Idx+2)=GradTmp(3*i) |
208 |
Idx=Idx+3 |
209 |
END DO |
210 |
DEALLOCATE(GradTmp) |
211 |
CASE ("CART","HYBRID") |
212 |
Grad=GradCart |
213 |
CASE DEFAULT |
214 |
WRITE(*,*) "Coord=",AdjustL(Trim(COORD))," not yet implemented in Egradients. STOP" |
215 |
STOP |
216 |
END SELECT |
217 |
|
218 |
DEALLOCATE(GradCart) |
219 |
DEALLOCATE(x,y,z) |
220 |
|
221 |
IF (Debug) WRITE(*,*) 'DBG EGRAD grad',grad(1:NCoord) |
222 |
if (debug) Call Header("Egrad Over") |
223 |
|
224 |
!999 CONTINUE |
225 |
!if (.NOT.Ftmp) WRITE(*,*) 'We should not be here !!!!' |
226 |
!STOP |
227 |
! ====================================================================== |
228 |
end subroutine Egrad |
229 |
|