source: hbond.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: 10.7 KB
Line 
1!**************************************************************
2!
3! This file contains the subroutines: hbond,chhb,ishybd,
4! ishybdo,nursat,interhbond
5!
6! Copyright 2003-2005 Frank Eisenmenger, U.H.E. Hansmann,
7! Shura Hayryan, Chin-Ku
8! Copyright 2007 Frank Eisenmenger, U.H.E. Hansmann,
9! Jan H. Meinke, Sandipan Mohanty
10!
11! **************************************************************
12
13
14 subroutine hbond(nml,mhb,ipr)
15! .................................................................
16! PURPOSE: find hydrogen bonds in molecule 'nml'
17!
18! prints HBonds, if ipr > 0
19!
20! OUTPUT: mhb - number of hyd.bds. of type i->i+4
21!
22! to INCL.H:
23!
24! ntyhb - number of different types of hyd. bds. found
25! nutyhb - number of hyd.bds. found for each type
26! ixtyhb - index for each type of hyd. bd. composed as
27! (atom idx. of H) * 1000 + atm.idx. of acceptor
28!
29! CALLS: chhb,ishybd (ishybdo),nursat
30!
31!................................................................
32
33 include 'INCL.H'
34
35 integer nml
36
37 integer mhb
38
39 integer nursat
40 double precision atbase
41
42 integer ipr, i2s, i, i14, i1, i1s, i2, ia, ixhb, id, ifivr, ii, ih
43 integer ims, io, iv, ivw, ix, k, jj, j, n, na, nd, ntlvr, iat
44
45!f2py intent(out) mhb
46 parameter (atbase=mxat)
47 logical ishb
48
49
50 do i=1,mxtyhb
51 nutyhb(i) = 0
52 ixtyhb(i) = 0
53 enddo
54 ntyhb=0
55 if (nml.eq.0) then
56 ntlvr = nvr
57 ifivr = ivrml1(1)
58 else
59 ntlvr=nvrml(nml)
60 if (ntlvr.eq.0) then
61 write (logString, '(a,i4)')
62 & ' hbond> No variables defined in molecule #',nml
63 return
64 endif
65
66 ifivr=ivrml1(nml)
67! Index of last moving set
68 i1s=imsml1(nml)+nmsml(nml)
69 endif
70! Loop over all variables
71 do io=ifivr+ntlvr-1,ifivr,-1
72! Get index of variable
73 iv=iorvr(io)
74! Index of next to last moving set
75 i2s=i1s-1
76! Index of moving set belonging to iv
77 i1s=imsvr1(iv)
78! Loop over all moving sets between the one belonging to iv and the
79! next to last one
80 do ims=i1s,i2s
81! First atom in moving set
82 i1=latms1(ims)
83! Last atom in moving set
84 i2=latms2(ims)
85! Loop over all atoms in moving set.
86 do i=i1,i2
87! Loop over van der Waals domains of atom i
88 do ivw=ivwat1(i),ivwat2(i)
89! Loop over atoms in van der Waals domain.
90 do j=lvwat1(ivw),lvwat2(ivw)
91
92 call ishybd(i,j,ishb,ih,ia) ! Thornton criteria
93
94 if (ishb) then
95
96 ixhb=ih*atbase+ia
97
98 do k=1,ntyhb
99 if (ixhb.eq.ixtyhb(k)) then
100 nutyhb(k)=nutyhb(k)+1
101 goto 1
102 endif
103 enddo
104
105 if (ntyhb.lt.mxtyhb) then
106 ntyhb=ntyhb+1
107 nutyhb(ntyhb)=1
108 ixtyhb(ntyhb)=ixhb
109 else
110 write (logString, *)
111 & ' hbond> increase parameter MXTYHB'
112 stop
113 endif
114
115 endif ! have h.b.
116
117 1 enddo ! ... atoms j
118 enddo ! ... vdW-domains of i
119
120 do i14=i14at1(i),i14at2(i) ! over 1-4 partners of 'i'
121 j=l14at(i14)
122
123 call ishybd(i,j,ishb,ih,ia) ! Thornton criteria
124! call ishybdo(i,j,ishb,ih,ia)
125
126 if (ishb) then
127
128 ixhb=ih*atbase+ia
129
130 do k=1,ntyhb
131 if (ixhb.eq.ixtyhb(k)) then
132 nutyhb(k)=nutyhb(k)+1
133 goto 2
134 endif
135 enddo
136
137 if (ntyhb.lt.mxtyhb) then
138 ntyhb=ntyhb+1
139 nutyhb(ntyhb)=1
140 ixtyhb(ntyhb)=ixhb
141 else
142 write (logString, *)
143 & ' hbond> increase parameter MXTYHB'
144 stop
145 endif
146
147 endif ! have h.b.
148
149 2 enddo ! ... 1-4-partners of i
150
151 enddo ! ... atoms i
152 enddo ! ... m.s.
153
154 enddo ! ... variables
155
156 mhb=0
157
158! do inhb=1,ntyhb
159! mhb = mhb+nutyhb(inhb)
160! enddo
161
162 if (ipr.gt.0) write (logString, '(1x,a,/)')
163 & ' hbond> Hydrogen Bonds:'
164
165 if (ntyhb.gt.0) then
166
167 ii = 0
168 do i=1,ntyhb
169 jj = nutyhb(i)
170 do j = 1,jj
171 ii = ii + 1
172 ix = ixtyhb(ii)
173
174 id =ix / atbase ! donor atom
175 nd = nursat(id) ! & residue
176
177 ia = ix - id * atbase ! acceptor atom
178 na = nursat(ia)
179
180 n = nd - na
181
182 if (n.gt.4) mhb = mhb+1 ! only count these
183
184 if (ipr.gt.0) then
185
186 if (n.gt.0) then
187 write(*,'(1x,i3,a2,a4,a3,i3,1x,a4,a7,a4,a3,i3,1x,a4,a9,
188 & i2)')
189 & ii,') ',nmat(ia),' ( ',na,seq(na),' ) <-- ',nmat(id),
190 & ' ( ', nd,seq(nd),' ) = i +',n
191 else
192 write(*,'(1x,i3,a2,a4,a3,i3,1x,a4,a7,a4,a3,i3,1x,a4,a9,
193 & i2)')
194 & ii,') ',nmat(ia),' ( ',na,seq(na),' ) <-- ',nmat(id),
195 & ' ( ', nd,seq(nd),' ) = i -',abs(n)
196 endif
197
198 call chhb(ia,id)
199 endif
200
201 enddo
202 enddo
203 endif
204
205 return
206 end
207! .....................................................................
208! Calculates hydrogen bonds between different chains.
209!
210! @return number of intermolecular hydrogen bonds. Returns 0 if only
211! one molecule is present. The value is returned in the
212! variable mhb.
213!
214! @author Jan H. Meinke <j.meinke@fz-juelich.de>
215!
216! .....................................................................
217 subroutine interhbond(mhb)
218
219 include 'INCL.H'
220
221!f2py intent(out) mhb
222
223 logical ishb
224
225 integer*4 mhb
226 integer iml, ires, jml, jat, jres
227
228 integer ia, iat, ih
229
230 mhb = 0
231 do iml = 1, ntlml
232 do jml = iml + 1, ntlml
233 mmhb(iml, jml) = 0
234 do ires= irsml1(iml), irsml2(iml)
235 do jres= irsml1(jml), irsml2(jml)
236 do iat = iatrs1(ires), iatrs2(ires)
237 do jat = iatrs1(jres), iatrs2(jres)
238 call ishybd(iat,jat,ishb,ih,ia)
239 if (ishb) then
240 mhb = mhb + 1
241 mmhb(iml, jml) = mmhb(iml, jml) + 1
242 endif
243 enddo ! jat
244 enddo ! iat
245 enddo ! jres
246 enddo ! ires
247 mmhb(jml, iml) = mmhb(iml, jml)
248 enddo ! jml
249 enddo ! iml
250
251 end ! subroutine interhbond
252! ************************
253 subroutine chhb(i,j)
254
255 include 'INCL.H'
256 integer i
257
258 integer j
259
260 integer ih
261
262 integer ia, ib, id, ihb
263
264 double precision valang
265
266 double precision cdah, cdad, adab, adha, ahab, dah, dad
267
268 ihb = ihbty(ityat(i),ityat(j))
269
270 if (ihb.gt.0) then
271 ih=i
272 ia=j
273 else
274 ih=j
275 ia=i
276 endif
277
278 dah=sqrt((xat(ih)-xat(ia))**2+(yat(ih)-yat(ia))**2+
279 & (zat(ih)-zat(ia))**2)
280
281 id=iowat(ih)
282
283 dad=sqrt((xat(id)-xat(ia))**2+(yat(id)-yat(ia))**2+
284 & (zat(id)-zat(ia))**2)
285 adha=valang(id,ih,ia)*crd
286
287 ib=iowat(ia)
288
289 ahab=valang(ih,ia,ib)*crd
290 adab=valang(id,ia,ib)*crd
291
292 write (logString, *) ' '
293 write (logString, *) ' Dah: ',dah,' Dad: ',dad
294 write (logString, *) ' Adha: ',adha,' Ahab: ',adha,' Adab: ',adab
295 write (logString, *) ' '
296
297 return
298 end
299! *************************************
300 subroutine ishybd(i,j,ishb,ih,ia)
301
302
303! ..........................................................
304! PURPOSE: checks for hydrogen bond between atoms 'i' & 'j'
305! according to geometric criteria
306!
307! OUTPUT: logical 'ishb' - true, if have Hydrogen bond
308! ih - index of Hydrogen atom
309! ia - index of Acceptor atom
310!
311! [I.K.McDonald,J.M.Thornton,Satisfying hydrogen bond
312! potential in proteins.J.Mol.Biol.238(5),777-793 (1994)]
313!
314! D: hydrogen(=H) donor, A: acceptor, B: atom bound to A
315!
316! Dis_HA <= 2.5 & Dis_DA <= 3.9 & Angle(D-H-A) > 90 &
317! Angle(H-A-B) > 90 & Angle(D-A-B) > 90
318! ..........................................................
319
320 include 'INCL.H'
321 double precision cdah, cang, cahb, cdad
322 integer i
323
324 integer j
325
326 integer ih
327
328 integer ia, ib, id, ihb
329
330 double precision valang
331
332 parameter (cdad=3.9d0,
333 & cdah=2.5d0,
334 & cang=110.d0)
335! # cang=90.d0)
336
337 logical ishb
338
339 cahb = cang * cdr
340 ishb = .false.
341
342
343 if (i.le.0.or.j.le.0) return
344
345 ihb = ihbty(ityat(i),ityat(j))
346
347 if (ihb.eq.0) then
348 return
349 elseif (ihb.gt.0) then
350 ih=i
351 ia=j
352 else
353 ih=j
354 ia=i
355 endif
356
357 if (sqrt((xat(ih)-xat(ia))**2+(yat(ih)-yat(ia))**2+
358 & (zat(ih)-zat(ia))**2).gt.cdah) return
359
360 id=iowat(ih)
361
362 if (id.le.0.or.sqrt((xat(id)-xat(ia))**2+(yat(id)-yat(ia))**2+
363 & (zat(id)-zat(ia))**2).gt.cdad
364 & .or.valang(id,ih,ia).lt.cahb) return
365
366 ib=iowat(ia)
367
368 if (ib.gt.0.and.valang(ih,ia,ib).ge.cahb
369 & .and.valang(id,ia,ib).ge.cahb) ishb=.true.
370
371 return
372 end
373! **************************************
374
375 subroutine ishybdo(i,j,ishb,ih,ia)
376
377! ..........................................................
378! PURPOSE: checks for hydrogen bond between atoms 'i' & 'j'
379! according to geometric criteria
380!
381! OUTPUT: logical 'ishb' - true, if have Hydrogen bond
382! ih - index of Hydrogen atom
383! ia - index of Acceptor atom
384!
385! D: hydrogen(=H) donor, A: acceptor
386!
387! Dis_AH <= 2.5 & Angle(D-H-A) >= 160
388! ...........................................................
389
390 include 'INCL.H'
391 integer i
392
393 integer j
394
395 integer ih
396
397 double precision cahb, cdah, cang, valang
398
399 integer ia, ihb, id
400 parameter (cdah=2.5d0,
401 & cang=140.d0)
402
403 logical ishb
404
405 cahb = cang * cdr
406
407
408 ishb = .false.
409
410 if (i.le.0.or.j.le.0) return
411
412 ihb = ihbty(ityat(i),ityat(j))
413
414 if (ihb.eq.0) then
415 return
416 elseif (ihb.gt.0) then
417 ih=i
418 ia=j
419 else
420 ih=j
421 ia=i
422 endif
423
424 if (sqrt((xat(ih)-xat(ia))**2+(yat(ih)-yat(ia))**2+
425 & (zat(ih)-zat(ia))**2).gt.cdah) return
426
427 id=iowat(ih)
428
429 if (id.gt.0.and.valang(id,ih,ia).ge.cahb) ishb=.true.
430
431 return
432 end
433
Note: See TracBrowser for help on using the repository browser.