source: opeflx.f

Last change on this file was 38d77eb, checked in by baerbaer <baerbaer@…>, 14 years ago

Redirected standard out to logString.

SMMP produced a lot of log messages. This became an issue when run in massively
parallel environments. I replaced all writes to standard out to a write to logString.
The next step is to pass this string to a function that writes the messages to a log
file according to their importance and the chosen log level.

git-svn-id: svn+ssh://svn.berlios.de/svnroot/repos/smmp/trunk@34 26dc1dd8-5c4e-0410-9ffe-d298b4865968

  • Property mode set to 100644
File size: 13.7 KB
Line 
1! **************************************************************
2!
3! This file contains the subroutines: opeflx,gdtflx
4!
5! Copyright 2003-2005 Frank Eisenmenger, U.H.E. Hansmann,
6! Shura Hayryan, Chin-Ku
7! Copyright 2007 Frank Eisenmenger, U.H.E. Hansmann,
8! Jan H. Meinke, Sandipan Mohanty
9!
10! **************************************************************
11
12 subroutine opeflx(nml)
13
14! ......................................................................
15! PURPOSE: Calculate internal energy for FLEX dataset and its partial
16! derivatives vs. variables using recursive algorithm from:
17! Noguti T, Go N, J Phys Soc (Japan) v52 3685-3690 1984; Abe H,
18! Braun W, Noguti T, Go N, Comp Chem v8 239-247 1984; Mazur A K,
19! Abagyan R A, J Biomol Struct Dyn v6 815-832, which I modified
20! for atomic forces instead of simple derivatives (see Lavery R,
21! Sklenar H, Zakrzewska K, Pullman B, J Biomol Struct Dyn v3
22! 989-1014 1986)
23!
24! CALLS: gdtflx
25! ......................................................................
26
27 include 'INCL.H'
28
29 integer nml
30
31 double precision enyflx
32
33 integer ims, i1a, i14, i1, i, i1s, i2s, i2, i2a, j, io, ifivr
34 integer iad, ia, ib, ic, ilavr, ijhb, iowh, it, ity, iv, ivw, jty
35 integer jowh, lad, ntlvr
36 double precision cqi,c, ep, e0, deyds, cqiqj, eynw, eyol, deel
37 double precision ex, cth, dhb, deyrp, eel, eyds, ey, eyrp, ez
38 double precision fvr, hhb, gda, gdn, vrn, p, py, px, pz, yfiv
39
40 double precision rij6, rij12, rij, rij2, sr, sr2, vr, xti
41 double precision xb, xj, xfriv, xfiv, xfi, xfji, xij, xi, xsr
42 double precision yb, yfi, yfji, yfriv, yi, yij, yj, zfi, yti
43 double precision zb, zfiv, zfji, zfriv, zij, zi, zj, zti, del
44
45
46 double precision xfat, yfat, zfat, xtat, ytat, ztat, xfvr, yfvr
47 double precision zfvr, xfrvr, yfrvr, zfrvr
48 dimension xfat(mxat),yfat(mxat),zfat(mxat),
49 & xtat(mxat),ytat(mxat),ztat(mxat),
50 & xfvr(mxvr),yfvr(mxvr),zfvr(mxvr),
51 & xfrvr(mxvr),yfrvr(mxvr),zfrvr(mxvr)
52
53
54 eyel=0.d0
55 eyvw=0.d0
56 eyhb=0.d0
57 eyvr=0.d0
58 eysm=0.d0
59
60 ntlvr=nvrml(nml)
61 if (ntlvr.eq.0) then
62 write (logString, '(a,i4)')
63 & ' opeflx> No variables defined in molecule #',nml
64 return
65 endif
66
67 ifivr=ivrml1(nml)
68 ilavr=ifivr+ntlvr-1
69
70 do i=ifivr,ilavr
71 gdeyvr(i)=0.d0
72 xfvr(i)=0.0
73 yfvr(i)=0.0
74 zfvr(i)=0.0
75 xfrvr(i)=0.0
76 yfrvr(i)=0.0
77 zfrvr(i)=0.0
78 enddo
79
80 do i=iatrs1(irsml1(nml)),iatrs2(irsml2(nml))
81 xfat(i)=0.d0
82 yfat(i)=0.d0
83 zfat(i)=0.d0
84
85 xtat(i)=0.d0
86 ytat(i)=0.d0
87 ztat(i)=0.d0
88 enddo
89
90 i1s=imsml1(nml)+nmsml(nml)
91 i1a=iadml1(nml)+nadml(nml)
92
93 do io=ilavr,ifivr,-1 ! ______ loop over variables in desc. order
94
95 iv=iorvr(io) ! index of var.
96 ia=iatvr(iv) ! prim.mv.at
97
98 ib=iowat(ia) ! base
99 xb=xat(ib)
100 yb=yat(ib)
101 zb=zat(ib)
102
103 it=ityvr(iv) ! type
104 ic=iclvr(iv) ! class
105
106 fvr=0.d0
107
108 if (it.eq.3) then ! torsion
109 ex=xtoat(ib)
110 ey=ytoat(ib)
111 ez=ztoat(ib)
112
113 vr=toat(ia)
114
115 e0=e0to(ic)
116 if (e0.ne.0.d0) then
117 vrn=vr*rnto(ic)
118 eyvr=eyvr+e0*(1.d0+sgto(ic)*cos(vrn))
119 fvr=esnto(ic)*sin(vrn)
120 endif
121
122 elseif (it.eq.2) then ! b.angle
123 ex=xbaat(ia)
124 ey=ybaat(ia)
125 ez=zbaat(ia)
126
127 vr=baat(ia)
128
129 elseif (it.eq.1) then ! b.length
130 ex=xtoat(ia)
131 ey=ytoat(ia)
132 ez=ztoat(ia)
133
134 vr=blat(ia)
135
136 endif
137
138! ============================================ Energies & Atomic forces
139
140 xfiv=0.d0
141 yfiv=0.d0
142 zfiv=0.d0
143 xfriv=0.d0
144 yfriv=0.d0
145 zfriv=0.d0
146
147 i2s=i1s-1 ! last m.s per 'iv'
148 i1s=imsvr1(iv) ! 1st m.s
149
150 do ims=i1s,i2s ! __ loop over m.s
151
152 i1=latms1(ims)
153 i2=latms2(ims)
154
155 do i=i1,i2 ! __ loop over atoms i ===================
156
157 ity=ityat(i)
158 cqi=conv*cgat(i)
159
160 xi=xat(i)
161 yi=yat(i)
162 zi=zat(i)
163
164 xfi=xfat(i)
165 yfi=yfat(i)
166 zfi=zfat(i)
167
168 xti=xtat(i)
169 yti=ytat(i)
170 zti=ztat(i)
171
172 do ivw=ivwat1(i),ivwat2(i) ! loop over vdW-domains of 'i'
173 do j=lvwat1(ivw),lvwat2(ivw) ! atoms j
174
175 jty=ityat(j)
176
177 xj=xat(j)
178 yj=yat(j)
179 zj=zat(j)
180
181 xij=xj-xi
182 yij=yj-yi
183 zij=zj-zi
184
185 rij2=xij*xij+yij*yij+zij*zij
186 rij6=rij2**3
187 rij12=rij6*rij6
188 rij=sqrt(rij2)
189
190 cqiqj=cqi*cgat(j)
191
192 if (epsd) then
193
194 sr=slp_f*rij
195 sr2=sr*sr
196 xsr=(plt-1.d0)*exp(-sr)/2.d0
197 ep=plt-(sr2+2.d0*sr+2.d0)*xsr
198
199 eel=cqiqj/(rij*ep)
200 eyel=eyel+eel
201 deel=eel+cqiqj*(slp_f*sr2*xsr)/(ep*ep)
202
203 else
204
205 eel=cqiqj/rij
206 eyel=eyel+eel
207 deel=eel
208
209 endif
210
211
212 eyrp=aij(ity,jty)/rij12
213 eyds=cij(ity,jty)/rij6
214
215 c=(-12.d0*eyrp+6.d0*eyds- deel)/rij
216
217 xij=xij/rij
218 yij=yij/rij
219 zij=zij/rij
220
221 xfji=c*xij
222 yfji=c*yij
223 zfji=c*zij
224
225 ijhb=ihbty(ity,jty)
226 if (ijhb.ne.0.and.rij.le.cohb) then ! HB Possible
227 if (ijhb.gt.0) then ! i=H,j=acceptor
228 iowh=iowat(i)
229 px=xi-xat(iowh)
230 py=yi-yat(iowh)
231 pz=zi-zat(iowh)
232 else ! i=acceptor,j=H
233 jowh=iowat(j)
234 px=xat(jowh)-xj
235 py=yat(jowh)-yj
236 pz=zat(jowh)-zj
237 endif
238
239 p=sqrt(px*px+py*py+pz*pz)
240 px=px/p
241 py=py/p
242 pz=pz/p
243
244 cth=xij*px+yij*py+zij*pz
245
246 if (cth.gt.0.d0) then
247
248 deyrp=(ahb(ity,jty)-aij(ity,jty))/rij12
249 deyds=(chb(ity,jty)-cij(ity,jty))/rij6
250
251 dhb=deyrp-deyds
252 eyhb=eyhb+eyrp-eyds+cth*dhb
253
254 if (ijhb.gt.0) then ! i=H
255 xti=xti +dhb * (zij*py-yij*pz)
256 yti=yti +dhb * (xij*pz-zij*px)
257 zti=zti +dhb * (yij*px-xij*py)
258 else ! j=H
259 xtat(j)=xtat(j) +dhb * (zij*py-yij*pz)
260 ytat(j)=ytat(j) +dhb * (xij*pz-zij*px)
261 ztat(j)=ztat(j) +dhb * (yij*px-xij*py)
262 endif
263
264 dhb=dhb/rij
265 hhb=cth*(7.d0*deyds-13.d0*deyrp)/rij
266
267 xfji=xfji+ dhb*px+ hhb*xij
268 yfji=yfji+ dhb*py+ hhb*yij
269 zfji=zfji+ dhb*pz+ hhb*zij
270! __________________________________________________ No Hydrogen Bond
271 else
272 eyvw=eyvw+eyrp-eyds
273 endif
274 else
275 eyvw=eyvw+eyrp-eyds
276 endif
277
278 xfi=xfi+xfji
279 yfi=yfi+yfji
280 zfi=zfi+zfji
281
282 xfat(j)=xfat(j)-xfji
283 yfat(j)=yfat(j)-yfji
284 zfat(j)=zfat(j)-zfji
285
286 enddo ! ... atoms j
287 enddo ! ... vdW-domains of i
288
289 do i14=i14at1(i),i14at2(i) ! loop over 1-4 partn. of 'i'
290 j=l14at(i14)
291
292 jty=ityat(j)
293
294 xj=xat(j)
295 yj=yat(j)
296 zj=zat(j)
297
298 xij=xj-xi
299 yij=yj-yi
300 zij=zj-zi
301
302 rij2=xij*xij+yij*yij+zij*zij
303 rij6=rij2**3
304 rij12=rij6*rij6
305 rij=sqrt(rij2)
306
307 cqiqj=cqi*cgat(j)
308
309 if (epsd) then
310
311 sr=slp_f*rij
312 sr2=sr*sr
313 xsr=(plt-1.d0)*exp(-sr)/2.d0
314 ep=plt-(sr2+2.d0*sr+2.d0)*xsr
315
316 eel=cqiqj/(rij*ep)
317 eyel=eyel+eel
318 deel=eel+cqiqj*(slp_f*sr2*xsr)/(ep*ep)
319
320 else
321
322 eel=cqiqj/rij
323 eyel=eyel+eel
324 deel=eel
325
326 endif
327
328 eyrp=a14(ity,jty)/rij12
329 eyds=cij(ity,jty)/rij6
330
331 c=(-12.d0*eyrp+6.d0*eyds- deel )/rij
332
333 xij=xij/rij
334 yij=yij/rij
335 zij=zij/rij
336
337 xfji=c*xij
338 yfji=c*yij
339 zfji=c*zij
340
341 ijhb=ihbty(ity,jty)
342 if (ijhb.ne.0.and.rij.le.cohb) then ! HB Possible
343 if (ijhb.gt.0) then ! i=H,j=acceptor
344 iowh=iowat(i)
345 px=xi-xat(iowh)
346 py=yi-yat(iowh)
347 pz=zi-zat(iowh)
348 else ! i=acceptor,j=H
349 jowh=iowat(j)
350 px=xat(jowh)-xj
351 py=yat(jowh)-yj
352 pz=zat(jowh)-zj
353 endif
354
355 p=sqrt(px*px+py*py+pz*pz)
356 px=px/p
357 py=py/p
358 pz=pz/p
359
360 cth=xij*px+yij*py+zij*pz
361
362 if (cth.gt.0.d0) then
363
364 deyrp=(ahb(ity,jty)-a14(ity,jty))/rij12
365 deyds=(chb(ity,jty)-cij(ity,jty))/rij6
366
367 dhb=deyrp-deyds
368 eyhb=eyhb+eyrp-eyds+cth*dhb
369
370 if (ijhb.gt.0) then ! i=H
371 xti=xti -dhb * (yij*pz-zij*py)
372 yti=yti -dhb * (zij*px-xij*pz)
373 zti=zti -dhb * (xij*py-yij*px)
374 else ! j=H
375 xtat(j)=xtat(j) +dhb * (yij*pz-zij*py)
376 ytat(j)=ytat(j) +dhb * (zij*px-xij*pz)
377 ztat(j)=ztat(j) +dhb * (xij*py-yij*px)
378 endif
379
380 dhb=dhb/rij
381 hhb=cth*(7.d0*deyds-13.d0*deyrp)/rij
382
383 xfji=xfji+ dhb*px+ hhb*xij
384 yfji=yfji+ dhb*py+ hhb*yij
385 zfji=zfji+ dhb*pz+ hhb*zij
386! __________________________________________________ No Hydrogen Bond
387 else
388 eyvw=eyvw+eyrp-eyds
389 endif
390 else
391 eyvw=eyvw+eyrp-eyds
392 endif
393
394 xfi=xfi+xfji
395 yfi=yfi+yfji
396 zfi=zfi+zfji
397
398 xfat(j)=xfat(j)-xfji
399 yfat(j)=yfat(j)-yfji
400 zfat(j)=zfat(j)-zfji
401
402 enddo ! ... 1-4-partners of i
403
404 xfat(i)=xfi
405 yfat(i)=yfi
406 zfat(i)=zfi
407 xtat(i)=xti
408 ytat(i)=yti
409 ztat(i)=zti
410
411 xfiv=xfiv + xfi
412 yfiv=yfiv + yfi
413 zfiv=zfiv + zfi
414
415 xfriv=xfriv + yfi*zi-zfi*yi + xti
416 yfriv=yfriv + zfi*xi-xfi*zi + yti
417 zfriv=zfriv + xfi*yi-yfi*xi + zti
418
419 enddo ! ... atoms i
420 enddo ! ... m.s.
421
422 i2a=i1a-1 ! last 'added' var.
423 i1a=iadvr1(iv) ! 1st 'added' var.
424
425 do iad=i1a,i2a
426 lad=ladvr(iad)
427 xfiv=xfiv+xfvr(lad)
428 yfiv=yfiv+yfvr(lad)
429 zfiv=zfiv+zfvr(lad)
430 xfriv=xfriv+xfrvr(lad)
431 yfriv=yfriv+yfrvr(lad)
432 zfriv=zfriv+zfrvr(lad)
433 enddo
434
435 xfvr(iv)=xfiv
436 yfvr(iv)=yfiv
437 zfvr(iv)=zfiv
438 xfrvr(iv)=xfriv
439 yfrvr(iv)=yfriv
440 zfrvr(iv)=zfriv
441
442 if (it.eq.3.or.it.eq.2) then ! torsion,b.angle
443
444 gdeyvr(iv)= (ey*zb-ez*yb)*xfiv+(ez*xb-ex*zb)*yfiv+
445 & (ex*yb-ey*xb)*zfiv
446 & +ex*xfriv+ey*yfriv+ez*zfriv -fvr
447
448 elseif (it.eq.1) then ! b.length
449
450 gdeyvr(iv)= -(ex*xfiv+ey*yfiv+ez*zfiv) -fvr
451
452 endif
453
454 if (tesgrd) call gdtflx(nml,iv) ! grad.-test
455
456 enddo ! ... variables
457
458 eysm= eyel+eyvw+eyhb+eyvr
459
460 return
461 end
462! *****************************
463 subroutine gdtflx(nml,iv)
464
465! .....................................................................
466! PURPOSE: calculate partial derivative of internal energy for molecule
467! 'nml' vs. variable 'iv' NUMERICALLY and compare with
468! its value obtained analytically
469!
470! CALLS: setvar, enyflx
471! .....................................................................
472
473 include 'INCL.H'
474
475 integer nml, iv, it, i
476 double precision enyflx
477
478 double precision del
479 parameter (del=1.d-7)
480
481 double precision eynw, eyol, gda, gdn, ovr, vlvrx
482 dimension vlvrx(mxvr)
483
484! ____________________________ get & save values of variables
485 do i=1,ivrml1(ntlml)+nvrml(ntlml)-1
486 it=ityvr(i) ! type
487 if (it.eq.3) then ! torsion
488 vlvrx(i)=toat(iatvr(i))
489 elseif (it.eq.2) then ! b.angle
490 vlvrx(i)=baat(iatvr(i))
491 elseif (it.eq.1) then ! b.length
492 vlvrx(i)=blat(iatvr(i))
493 endif
494 enddo
495
496 ovr=vlvrx(iv)
497 eyol=enyflx(nml)
498
499 vlvrx(iv)=ovr+del ! change variable 'iv' by 'del'
500 call setvar(nml,vlvrx)
501 eynw=enyflx(nml) ! new energy
502
503 gdn=(eynw-eyol)/del ! numerical derivative
504 gda=gdeyvr(iv) ! analytical der.
505
506 write (logString, '(1x,2a,2(e12.6,a))') nmvr(iv),': ',gda,' (',
507 & abs(gda-gdn),')'
508
509! _________________________ restore
510 vlvrx(iv)=ovr
511 call setvar(nml,vlvrx)
512
513 return
514 end
515
Note: See TracBrowser for help on using the repository browser.