root / src / Egrad_baker.f90 @ 1
Historique | Voir | Annoter | Télécharger (12,7 ko)
1 | 1 | equemene | subroutine Egrad_baker(E,Geom,Grad,NCoord,IGeom,IOpt,GeomCart,Flag_Opt_Geom,GeomOld,GeomCart_old) |
---|---|---|---|
2 | 1 | equemene | |
3 | 1 | equemene | ! This routines calculates the energy E and the gradient Grad of |
4 | 1 | equemene | ! a molecule with Geometry Geom (may be in internal coordinates), |
5 | 1 | equemene | ! using for now, either Gaussian or Ext, more general later. |
6 | 1 | equemene | |
7 | 1 | equemene | use Path_module, only : Nat,AtName,Coord,dzdc,indzmat,Nom,Atome,massat,unit, & |
8 | 1 | equemene | prog,NCart,XyzGeomF,IntCoordF,BTransInv, XyzGeomI, & |
9 | 1 | equemene | GeomOld_all,BTransInvF,BTransInv_local,UMatF,UMat_local & |
10 | 1 | equemene | , BprimT,a0,ScanCoord, Coordinate,NPrim,BMat_BakerT,BBT,BBT_inv & |
11 | 1 | equemene | , Order,OrderInv, XPrimitiveF |
12 | 1 | equemene | |
13 | 1 | equemene | ! IntCoordF(NGeomF,NCoord),GeomOld_all(NGeomF,NCoord) |
14 | 1 | equemene | ! allocated in Path.f90 |
15 | 1 | equemene | |
16 | 1 | equemene | use Io_module |
17 | 1 | equemene | IMPLICIT NONE |
18 | 1 | equemene | |
19 | 1 | equemene | ! Energy (calculated if F300K=.F., else estimated) |
20 | 1 | equemene | REAL(KREAL), INTENT (OUT) :: E |
21 | 1 | equemene | ! NCoord: Number of the degrees of freedom |
22 | 1 | equemene | ! IGeom: index of the geometry. |
23 | 1 | equemene | INTEGER(KINT), INTENT (IN) :: NCoord, IGeom, IOpt |
24 | 1 | equemene | ! Geometry at which gradient is calculated (cf Factual also): |
25 | 1 | equemene | REAL(KREAL), INTENT (INOUT) :: Geom(NCoord) |
26 | 1 | equemene | REAL(KREAL), INTENT (INOUT) :: GeomOld(NCoord) |
27 | 1 | equemene | ! Gradient calculated at Geom geometry: |
28 | 1 | equemene | REAL(KREAL), INTENT (OUT) :: Grad(NCoord) |
29 | 1 | equemene | ! Cartesian geometry corresponding to (Internal Geometry) Geom: |
30 | 1 | equemene | REAL(KREAL), INTENT (OUT) :: GeomCart(Nat,3) |
31 | 1 | equemene | REAL(KREAL), INTENT (IN) :: GeomCart_old(Nat,3) |
32 | 1 | equemene | LOGICAL,INTENT (IN) :: Flag_Opt_Geom |
33 | 1 | equemene | |
34 | 1 | equemene | ! ====================================================================== |
35 | 1 | equemene | |
36 | 1 | equemene | CHARACTER(LCHARS) :: LINE |
37 | 1 | equemene | LOGICAL :: debug |
38 | 1 | equemene | LOGICAL, SAVE :: first=.true. |
39 | 1 | equemene | LOGICAL, ALLOCATABLE :: Done(:) |
40 | 1 | equemene | REAL(KREAL), SAVE :: Eold=1.e6 |
41 | 1 | equemene | REAL(KREAL), ALLOCATABLE :: GradTmp(:) |
42 | 1 | equemene | REAL(KREAL), ALLOCATABLE :: GradCart(:) |
43 | 1 | equemene | REAL(KREAL), ALLOCATABLE :: x(:), y(:), z(:) |
44 | 1 | equemene | REAL(KREAL) :: d, a_val, Pi |
45 | 1 | equemene | |
46 | 1 | equemene | INTEGER(KINT) :: iat, jat, kat, i, j, n3at,IBeg,Idx |
47 | 1 | equemene | INTEGER(KINT), PARAMETER :: IOLOG=12, IOCOM=11 |
48 | 1 | equemene | |
49 | 1 | equemene | CHARACTER(LCHARS) :: ListName, CH32SVAR1 |
50 | 1 | equemene | CHARACTER(VLCHARS), SAVE :: RstrtCopy, RunCommand |
51 | 1 | equemene | LOGICAL, SAVE :: FCopyRstrt=.False., FOrderChecked=.False. |
52 | 1 | equemene | LOGICAL :: FTmp |
53 | 1 | equemene | INTEGER(kint) :: NbLists,LastIt,Nat4,NStep, Firstit |
54 | 1 | equemene | INTEGER(KINT), PARAMETER :: NbExtName=4 |
55 | 1 | equemene | |
56 | 1 | equemene | REAL(KREAL), ALLOCATABLE :: XPrimRef(:), XPrim(:) ! NPrim |
57 | 1 | equemene | |
58 | 1 | equemene | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
59 | 1 | equemene | ! |
60 | 1 | equemene | !!!!!!!! PFL Debug |
61 | 1 | equemene | ! |
62 | 1 | equemene | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
63 | 1 | equemene | REAL(KREAL), ALLOCATABLE :: GeomTmp(:,:) !3,Nat |
64 | 1 | equemene | REAL(KREAL), ALLOCATABLE :: V(:,:) ! (NCoord,NCoord) |
65 | 1 | equemene | REAL(KREAL), ALLOCATABLE :: XyzTmp(:,:) ! xyzatm(3,nat) |
66 | 1 | equemene | REAL(KREAL) :: tol |
67 | 1 | equemene | LOGICAL :: FAIL , FALOBBT,FALOBPrimt,FAloBBTInv |
68 | 1 | equemene | LOGICAL :: DebugPFL |
69 | 1 | equemene | |
70 | 1 | equemene | ! ====================================================================== |
71 | 1 | equemene | LOGICAL, EXTERNAL :: valid |
72 | 1 | equemene | ! ====================================================================== |
73 | 1 | equemene | |
74 | 1 | equemene | INTEGER(KINT) :: OptGeom |
75 | 1 | equemene | Namelist/path/OptGeom |
76 | 1 | equemene | |
77 | 1 | equemene | Pi=dacos(-1.0d0) |
78 | 1 | equemene | n3at=3*nat |
79 | 1 | equemene | |
80 | 1 | equemene | debug=valid('EGRAD_baker') |
81 | 1 | equemene | debugPFL=valid('bakerPFL') |
82 | 1 | equemene | if (debug) WRITE(*,*) '================ Entering Egrad_baker =================' |
83 | 1 | equemene | if (debug) Print *, 'EGrad.f90, L73, IOpt=', IOpt, ' Flag_Opt_Geom=',Flag_Opt_Geom |
84 | 1 | equemene | |
85 | 1 | equemene | Grad=0. |
86 | 1 | equemene | E=0. |
87 | 1 | equemene | |
88 | 1 | equemene | ALLOCATE(GradCart(3*Nat)) |
89 | 1 | equemene | ALLOCATE(x(Nat),y(Nat),z(Nat)) |
90 | 1 | equemene | ALLOCATE(XPrimRef(NPrim),XPrim(NPrim)) |
91 | 1 | equemene | |
92 | 1 | equemene | SELECT CASE (COORD) |
93 | 1 | equemene | CASE ('ZMAT') |
94 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad.f90, instead of Egrad_baker.f90. Stop" |
95 | 1 | equemene | CASE ('BAKER') |
96 | 1 | equemene | XPrimRef=XPrimitiveF(IGeom,:) |
97 | 1 | equemene | IF (Flag_Opt_Geom) THEN |
98 | 1 | equemene | IF (IOpt .EQ. 1) THEN |
99 | 1 | equemene | GeomCart(:,1) = XyzGeomI(IGeom,1,:) |
100 | 1 | equemene | GeomCart(:,2) = XyzGeomI(IGeom,2,:) |
101 | 1 | equemene | GeomCart(:,3) = XyzGeomI(IGeom,3,:) |
102 | 1 | equemene | ELSE |
103 | 1 | equemene | if (debug) WRITE(*,*) 'Before Call ConvertBakerInternal_cart; L82, Egrad.f90' |
104 | 1 | equemene | WRITE(*,*) 'GeomOld=' |
105 | 1 | equemene | WRITE(*,'(12(1X,F6.3))') GeomOld |
106 | 1 | equemene | WRITE(*,*) 'Geom=' |
107 | 1 | equemene | WRITE(*,'(12(1X,F6.3))') Geom |
108 | 1 | equemene | WRITE(*,*) 'GeomCart_old=' |
109 | 1 | equemene | WRITE(*,'(20(1X,F6.3))') GeomCart_old |
110 | 1 | equemene | Call ConvertBakerInternal_cart(GeomOld,Geom,GeomCart_old(:,1), & |
111 | 1 | equemene | GeomCart_old(:,2),GeomCart_old(:,3),x,y,z,XPrim,XPrimRef) |
112 | 1 | equemene | GeomCart(:,1) = x(:) |
113 | 1 | equemene | GeomCart(:,2) = y(:) |
114 | 1 | equemene | GeomCart(:,3) = z(:) |
115 | 1 | equemene | END IF |
116 | 1 | equemene | ELSE |
117 | 1 | equemene | IF (IOpt .EQ. 0) THEN |
118 | 1 | equemene | GeomCart(:,1) = XyzGeomF(IGeom,1,:) |
119 | 1 | equemene | GeomCart(:,2) = XyzGeomF(IGeom,2,:) |
120 | 1 | equemene | GeomCart(:,3) = XyzGeomF(IGeom,3,:) |
121 | 1 | equemene | ELSE |
122 | 1 | equemene | ! IntCoordF(NGeomF,NCoord), XyzGeomF(NGeomF,3,Nat) |
123 | 1 | equemene | ! Geom has to be converted into x,y,z |
124 | 1 | equemene | ! GeomOld_all and XyzGeomF are old geometry in internal and cartesian coords resp. |
125 | 1 | equemene | DO I=1,NCoord ! these variables are used in ConvertBakerInternal_cart() |
126 | 1 | equemene | BTransInv_local(I,:) = BTransInvF(IGeom,I,:) |
127 | 1 | equemene | UMat_local(:,I) = UMatF(IGeom,:,I) |
128 | 1 | equemene | END DO |
129 | 1 | equemene | if (debug) Print *, 'I am before Call ConvertBakerInternal_cart; L108, Egrad.f90' |
130 | 1 | equemene | WRITE(*,*) 'GeomOld=' |
131 | 1 | equemene | WRITE(*,'(12(1X,F6.3))') GeomOld |
132 | 1 | equemene | WRITE(*,*) 'Geom=' |
133 | 1 | equemene | WRITE(*,'(12(1X,F6.3))') Geom |
134 | 1 | equemene | WRITE(*,*) 'DBG Egrad_baker L129 GeomCart_old=' |
135 | 1 | equemene | DO I=1,Nat |
136 | 1 | equemene | WRITE(*,'(3(1X,F12.8))') XyzGeomF(IGeom,:,I) |
137 | 1 | equemene | END DO |
138 | 1 | equemene | |
139 | 1 | equemene | Call ConvertBakerInternal_cart(GeomOld_all(IGeom,:),Geom,XyzGeomF(IGeom,1,:), & |
140 | 1 | equemene | XyzGeomF(IGeom,2,:),XyzGeomF(IGeom,3,:),x,y,z,XPrim,XPrimRef) |
141 | 1 | equemene | if (debug) Print *, 'I am after Call ConvertBakerInternal_cart; L111, Egrad.f90' |
142 | 1 | equemene | DO I=1,NCoord ! these variables are used in ConvertBakerInternal_cart() |
143 | 1 | equemene | BTransInvF(IGeom,I,:) = BTransInv_local(I,:) |
144 | 1 | equemene | END DO |
145 | 1 | equemene | ! GeomCart is actually XyzCoord and has INTENT(OUT) in this subroutine. |
146 | 1 | equemene | GeomCart(:,1) = x(:) |
147 | 1 | equemene | GeomCart(:,2) = y(:) |
148 | 1 | equemene | GeomCart(:,3) = z(:) |
149 | 1 | equemene | END IF ! matches IF (IOpt .EQ. 0) THEN |
150 | 1 | equemene | END IF |
151 | 1 | equemene | CASE ('MIXED') |
152 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad.f90, instead of Egrad_baker.f90. Stop" |
153 | 1 | equemene | CASE ('CART') |
154 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad.f90, instead of Egrad_baker.f90. Stop" |
155 | 1 | equemene | CASE ('HYBRID') |
156 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad.f90, instead of Egrad_baker.f90. Stop" |
157 | 1 | equemene | CASE DEFAULT |
158 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD)," not recognized. Stop" |
159 | 1 | equemene | STOP |
160 | 1 | equemene | END SELECT |
161 | 1 | equemene | !WRITE(*,*) Nat |
162 | 1 | equemene | !WRITE(*,*) 'GeomCart in Egrad_baker' |
163 | 1 | equemene | !DO I=1,Nat |
164 | 1 | equemene | ! WRITE(*,'(1X,A,3(1X,F12.8))') AtName(I),GeomCart(I,1:3) |
165 | 1 | equemene | !END DO |
166 | 1 | equemene | |
167 | 1 | equemene | SELECT CASE (Prog) |
168 | 1 | equemene | CASE ('GAUSSIAN') |
169 | 1 | equemene | ! we call the Gaussian routine. |
170 | 1 | equemene | ! it will return the energy and the gradient in cartesian coordinates. |
171 | 1 | equemene | ! Gradient calculated at GeomCart geometry and stored in GradCart, which has INTENT(OUT) |
172 | 1 | equemene | Call egrad_gaussian(E,GeomCart,GradCart) |
173 | 1 | equemene | CASE ('MOPAC') |
174 | 1 | equemene | ! we call the Mopac routine. |
175 | 1 | equemene | ! it will return the energy and the gradient in cartesian coordinates. |
176 | 1 | equemene | ! Gradient calculated at GeomCart geometry and stored in GradCart, which has INTENT(OUT) |
177 | 1 | equemene | Call egrad_mopac(E,GeomCart,GradCart) |
178 | 1 | equemene | CASE ('VASP') |
179 | 1 | equemene | Call egrad_vasp(E,Geomcart,GradCart) |
180 | 1 | equemene | CASE ('TURBOMOLE') |
181 | 1 | equemene | Call egrad_turbomole(E,Geomcart,GradCart) |
182 | 1 | equemene | CASE ('EXT') |
183 | 1 | equemene | Call egrad_ext(E,Geomcart,GradCart) |
184 | 1 | equemene | CASE ('TEST') |
185 | 1 | equemene | Call egrad_test(Nat,E,Geomcart,GradCart) |
186 | 1 | equemene | CASE ('CHAMFRE') |
187 | 1 | equemene | Call egrad_chamfre(Nat,E,Geomcart,GradCart) |
188 | 1 | equemene | END SELECT |
189 | 1 | equemene | if (debug) THEN |
190 | 1 | equemene | WRITE(*,*) 'DBG EGRAD, GradCart read' |
191 | 1 | equemene | DO Iat=1,Nat |
192 | 1 | equemene | WRITE(*,'(1X,A10,3(1X,F15.8),A)') Trim(AtName(Iat)),GradCart(3*Iat-2:3*Iat) |
193 | 1 | equemene | END DO |
194 | 1 | equemene | END IF |
195 | 1 | equemene | |
196 | 1 | equemene | ! If (PROG=="VASP") STOP |
197 | 1 | equemene | |
198 | 1 | equemene | ! We have the gradient in Cartesian coordinates... which is ok for COORD=Cart, Hybrid |
199 | 1 | equemene | ! but we have to convert it into Zmat if COORD=Zmat |
200 | 1 | equemene | SELECT CASE (COORD) |
201 | 1 | equemene | CASE ("ZMAT") |
202 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad.f90, instead of Egrad_baker.f90. Stop" |
203 | 1 | equemene | CASE ('BAKER') |
204 | 1 | equemene | ! We have the gradient in Cartesian coordinates... which is ok for COORD=Cart, Hybrid |
205 | 1 | equemene | ! but we have to convert it into internal coordinates if COORD=Baker |
206 | 1 | equemene | !!!!!!!!!!!!!!!!!!!! |
207 | 1 | equemene | ! |
208 | 1 | equemene | ! PFL Debug |
209 | 1 | equemene | ! |
210 | 1 | equemene | !!!!!!!!!!!!!!!!!!!! |
211 | 1 | equemene | if (debugPFL) THEN |
212 | 1 | equemene | WRITE(*,*) "DBG PFL Egrad_Baker... computing dzdc" |
213 | 1 | equemene | if (.not.ALLOCATED(indzmat)) THEN |
214 | 1 | equemene | ALLOCATE(indzmat(Nat,5)) |
215 | 1 | equemene | IndZmat=0 |
216 | 1 | equemene | Indzmat(1,1)=1 |
217 | 1 | equemene | IndZmat(2,1)=2 |
218 | 1 | equemene | IndZmat(2,2)=1 |
219 | 1 | equemene | IndZmat(3,1)=3 |
220 | 1 | equemene | IndZmat(3,2)=2 |
221 | 1 | equemene | IndZmat(3,3)=1 |
222 | 1 | equemene | IndZmat(4,1)=4 |
223 | 1 | equemene | IndZmat(4,2)=3 |
224 | 1 | equemene | IndZmat(4,3)=2 |
225 | 1 | equemene | IndZmat(4,4)=1 |
226 | 1 | equemene | END IF |
227 | 1 | equemene | IF (.NOT.ALLOCATED(DzDc)) THEN |
228 | 1 | equemene | ALLOCATE(DzDc(3,nat,3,Nat)) |
229 | 1 | equemene | END IF |
230 | 1 | equemene | if (.NOT.Allocated(GeomTmp)) Allocate(GeomTmp(3,Nat)) |
231 | 1 | equemene | DO I=1,Nat |
232 | 1 | equemene | GeomTmp(1,I)=GeomCart(OrderInv(I),1) |
233 | 1 | equemene | GeomTmp(2,I)=GeomCart(OrderInv(i),2) |
234 | 1 | equemene | GeomTmp(3,I)=GeomCart(OrderInv(i),3) |
235 | 1 | equemene | END DO |
236 | 1 | equemene | |
237 | 1 | equemene | DzDc=0.d0 |
238 | 1 | equemene | CALL CalcBMat_int (nat, GeomTmp, indzmat, dzdc, massat,atome) |
239 | 1 | equemene | |
240 | 1 | equemene | END IF ! if (debugPFL) |
241 | 1 | equemene | !!!!!!!!!!!!!!!!!!!!!!!!! |
242 | 1 | equemene | ! Debugging PFL |
243 | 1 | equemene | !!!!!!!!!!!!!!!!!!!!!!!! |
244 | 1 | equemene | |
245 | 1 | equemene | WRITE(*,*) "BTransInv_local Trans that is used originally" |
246 | 1 | equemene | DO J=1,3*Nat |
247 | 1 | equemene | WRITE(*,'(50(1X,F12.8))') BTransInv_local(:,J) |
248 | 1 | equemene | END DO |
249 | 1 | equemene | |
250 | 1 | equemene | WRITE(*,*) "Calculating actual values using GeomCart" |
251 | 1 | equemene | if (.NOT.Allocated(GeomTmp)) Allocate(GeomTmp(3,Nat)) |
252 | 1 | equemene | GeomTmp(1,:)=GeomCart(1:Nat,1) |
253 | 1 | equemene | GeomTmp(2,:)=GeomCart(1:Nat,2) |
254 | 1 | equemene | GeomTmp(3,:)=GeomCart(1:Nat,3) |
255 | 1 | equemene | |
256 | 1 | equemene | ! Computing B^prim: |
257 | 1 | equemene | FAloBBT=.NOT.ALLOCATED(BBT) |
258 | 1 | equemene | FAloBBTInv=.NOT.ALLOCATED(BBT_inv) |
259 | 1 | equemene | FAloBPrimT=.NOT.ALLOCATED(BprimT) |
260 | 1 | equemene | if (FAloBPrimT) ALLOCATE(Bprimt(3*Nat,NPrim)) |
261 | 1 | equemene | if (FAloBBT) ALLOCATE(BBt(NCoord,NCoord)) |
262 | 1 | equemene | if (FAloBBTInv) ALLOCATE(BBt_inv(NCoord,NCoord)) |
263 | 1 | equemene | BprimT=0.d0 |
264 | 1 | equemene | ScanCoord=>Coordinate |
265 | 1 | equemene | I=0 |
266 | 1 | equemene | DO WHILE (Associated(ScanCoord%next)) |
267 | 1 | equemene | I=I+1 |
268 | 1 | equemene | SELECT CASE (ScanCoord%Type) |
269 | 1 | equemene | CASE ('BOND') |
270 | 1 | equemene | CALL CONSTRAINTS_BONDLENGTH_DER(Nat,ScanCoord%at1,ScanCoord%AT2,GeomTmp,BprimT(1,I)) |
271 | 1 | equemene | CASE ('ANGLE') |
272 | 1 | equemene | CALL CONSTRAINTS_BONDANGLE_DER(Nat,ScanCoord%At1,ScanCoord%AT2,ScanCoord%At3,GeomTmp,BprimT(1,I)) |
273 | 1 | equemene | CASE ('DIHEDRAL') |
274 | 1 | equemene | CALL CONSTRAINTS_TORSION_DER2(Nat,ScanCoord%At1,ScanCoord%AT2,ScanCoord%At3,ScanCoord%At4,GeomTmp,BprimT(1,I)) |
275 | 1 | equemene | END SELECT |
276 | 1 | equemene | ScanCoord => ScanCoord%next |
277 | 1 | equemene | END DO |
278 | 1 | equemene | |
279 | 1 | equemene | if (debug) THEN |
280 | 1 | equemene | WRITE(*,*) "Bprim " |
281 | 1 | equemene | DO J=1,3*Nat |
282 | 1 | equemene | WRITE(*,'(50(1X,F12.6))') BprimT(J,:) |
283 | 1 | equemene | END DO |
284 | 1 | equemene | END IF |
285 | 1 | equemene | |
286 | 1 | equemene | BMat_BakerT = 0.d0 |
287 | 1 | equemene | DO I=1,NCoord |
288 | 1 | equemene | DO J=1,NPrim |
289 | 1 | equemene | !BprimT is transpose of B^prim. |
290 | 1 | equemene | !B = UMat^T B^prim, B^T = (B^prim)^T UMat |
291 | 1 | equemene | BMat_BakerT(:,I)=BMat_BakerT(:,I)+BprimT(:,J)*UMat_local(J,I) |
292 | 1 | equemene | END DO |
293 | 1 | equemene | END DO |
294 | 1 | equemene | |
295 | 1 | equemene | BBT = 0.d0 |
296 | 1 | equemene | DO I=1,NCoord |
297 | 1 | equemene | DO J=1,3*Nat ! BBT(:,I) forms BB^T |
298 | 1 | equemene | BBT(:,I) = BBT(:,I) + BMat_BakerT(J,:)*BMat_BakerT(J,I) |
299 | 1 | equemene | END DO |
300 | 1 | equemene | END DO |
301 | 1 | equemene | |
302 | 1 | equemene | BBT_inv = 0.d0 |
303 | 1 | equemene | |
304 | 1 | equemene | Call GenInv(NCoord,BBT,BBT_inv,NCoord) |
305 | 1 | equemene | |
306 | 1 | equemene | ! ALLOCATE(V(NCoord,NCoord)) |
307 | 1 | equemene | ! tol = 1e-10 |
308 | 1 | equemene | ! ! BBT is destroyed by GINVSE |
309 | 1 | equemene | ! Call GINVSE(BBT,NCoord,NCoord,BBT_inv,NCoord,NCoord,NCoord,tol,V,NCoord,FAIL,IOOUT) |
310 | 1 | equemene | ! DEALLOCATE(V) |
311 | 1 | equemene | ! IF(Fail) Then |
312 | 1 | equemene | ! WRITE (*,*) "Failed in generalized inverse, L220, ConvertBaker...f90" |
313 | 1 | equemene | ! STOP |
314 | 1 | equemene | ! END IF |
315 | 1 | equemene | |
316 | 1 | equemene | !Print *, 'BBT_inv=' |
317 | 1 | equemene | DO J=1,NCoord |
318 | 1 | equemene | ! WRITE(*,'(50(1X,F10.2))') BBT_inv(:,J) |
319 | 1 | equemene | ! Print *, BBT_inv(:,J) |
320 | 1 | equemene | END DO |
321 | 1 | equemene | !Print *, 'End of BBT_inv' |
322 | 1 | equemene | |
323 | 1 | equemene | ! Calculation of (B^T)^-1 = (BB^T)^-1B: |
324 | 1 | equemene | BTransInv_local = 0.d0 |
325 | 1 | equemene | DO I=1,3*Nat |
326 | 1 | equemene | DO J=1,NCoord |
327 | 1 | equemene | BTransInv_local(:,I)=BTransInv_local(:,I)+BBT_inv(:,J)*BMat_BakerT(I,J) |
328 | 1 | equemene | END DO |
329 | 1 | equemene | END DO |
330 | 1 | equemene | |
331 | 1 | equemene | Print *, 'BMat_BakerT=' |
332 | 1 | equemene | DO J=1,3*Nat |
333 | 1 | equemene | WRITE(*,'(50(1X,F12.8))') BMat_BakerT(:,J) |
334 | 1 | equemene | END DO |
335 | 1 | equemene | |
336 | 1 | equemene | WRITE(*,*) "BTransInv_local Trans corresponding to GeomCart" |
337 | 1 | equemene | DO J=1,3*Nat |
338 | 1 | equemene | WRITE(*,'(50(1X,F12.8))') BTransInv_local(:,J) |
339 | 1 | equemene | END DO |
340 | 1 | equemene | |
341 | 1 | equemene | if (FAloBPrimT) DEALLOCATE(Bprimt) |
342 | 1 | equemene | if (FAloBBT) DEALLOCATE(BBt) |
343 | 1 | equemene | if (FAloBBTInv) DEALLOCATE(BBt_inv) |
344 | 1 | equemene | |
345 | 1 | equemene | Grad=0.d0 |
346 | 1 | equemene | DO I=1, 3*Nat |
347 | 1 | equemene | ! Here Grad is gradient in Baker coordinates. GradCart(3*Nat) |
348 | 1 | equemene | Grad(:) = Grad(:) + BTransInv_local(:,I)*GradCart(I) !BTransInv_local = BTransInv in Opt_Geom.f90 |
349 | 1 | equemene | END DO |
350 | 1 | equemene | !WRITE(*,*) 'In Egrad.f90, L263, Gradient:' |
351 | 1 | equemene | !WRITE(*,'(12(1X,F8.4))') Grad(1:NCoord) |
352 | 1 | equemene | CASE ('MIXED') |
353 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad.f90, instead of Egrad_baker.f90. Stop" |
354 | 1 | equemene | CASE ("CART","HYBRID") |
355 | 1 | equemene | WRITE (*,*) "Coord=",TRIM(COORD),": use Egrad.f90, instead of Egrad_baker.f90. Stop" |
356 | 1 | equemene | CASE DEFAULT |
357 | 1 | equemene | WRITE(*,*) "Coord=",AdjustL(Trim(COORD))," not yet implemented in Egradients. STOP" |
358 | 1 | equemene | STOP |
359 | 1 | equemene | END SELECT |
360 | 1 | equemene | |
361 | 1 | equemene | DEALLOCATE(GradCart) |
362 | 1 | equemene | DEALLOCATE(x,y,z) |
363 | 1 | equemene | IF (Debug) WRITE(*,*) 'DBG EGRAD_baker grad',grad(1:NCoord) |
364 | 1 | equemene | |
365 | 1 | equemene | if (debug) WRITE(*,*) '================ Egrad_baker Over ====================' |
366 | 1 | equemene | |
367 | 1 | equemene | !999 CONTINUE |
368 | 1 | equemene | !if (.NOT.Ftmp) WRITE(*,*) 'We should not be here !!!!' |
369 | 1 | equemene | !STOP |
370 | 1 | equemene | ! ====================================================================== |
371 | 1 | equemene | end subroutine Egrad_baker |