source: enylun.f@ c076b43

Last change on this file since c076b43 was c076b43, checked in by sandipan <sandipan@…>, 16 years ago

Updated version of Lund force field functions

The Lund force field functions were debugged. A new library file
lib.lun was added to describe the geometrical parameters of
protein molecules consistent with the force field. With this
library file defining the bond lengths and bond angles, energy
values calculated with SMMP with the Lund force field can be
directly compared with those calculated with PROFASI.

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

  • Property mode set to 100644
File size: 30.3 KB
Line 
1! *******************************************************************
2! SMMP version of Anders Irback's force field, to be called the Lund
3! force field. This file contains the function enylun, which in turn
4! calls all the terms in the energy function. The terms Bias (ebias),
5! Hydrogen bonds (ehbmm and ehbms), Hydrophobicity (ehp) and the
6! Excluded volume (eexvol and eloexv) are also implemented in this
7! file.
8!
9! Copyright 2007 Frank Eisenmenger, U.H.E. Hansmann,
10! Jan H. Meinke, Sandipan Mohanty
11!
12 subroutine set_local_constr(ires,frst,lst,root)
13 include 'INCL.H'
14 include 'incl_lund.h'
15 integer ires,frst,lst,root
16 do iat1=iCa(ires)+frst,iCa(ires)+lst
17 do iat2=iat1+1,iCa(ires)+lst
18 matcon(iat2-iat1,iat1)=0
19 matcon(iat1-iat2,iat2)=0
20 enddo
21 enddo
22 iat1=iCa(ires)+root
23 do iat2=iCa(ires)+root,iCa(ires)+lst
24 matcon(iat2-iat1,iat1)=0
25 matcon(iat1-iat2,iat2)=0
26 enddo
27 end subroutine set_local_constr
28
29 subroutine init_lundff
30 include 'INCL.H'
31 include 'incl_lund.h'
32 character mynm*4
33 logical prlvr
34
35 print *,'initializing Lund forcefield'
36! Some parameters in the Lund force field.
37! The correspondence between internal energy scale and kcal/mol
38 eunit=1.3315d0
39! eunit=1.0
40! Bias
41 kbias=100.0d0*eunit
42! print *,'Bias'
43! Hydrogen bonds
44 epshb1=3.1d0*eunit
45 epshb2=2.0d0*eunit
46 sighb=2.0d0
47 cthb=4.5d0
48 cthb2=cthb*cthb
49 powa=0.5d0
50 powb=0.5d0
51 blhb=-30.0d0*(((sighb/cthb)**10-(sighb/cthb)**12))/cthb2
52 alhb=-(5*((sighb/cthb)**12)-6*((sighb/cthb)**10))-blhb*cthb2
53 sighb2=sighb*sighb
54 cdon=1.0d0
55 cacc=(1.0d0/1.23d0)**powb
56 csacc=(1.0d0/1.25d0)**powb
57! print *,'Hydrogen bonds'
58! Hydrophobicity
59! print *,'Hydrophobicity with nhptyp = ',nhptyp
60
61 hpstrg(1)=0.0d0*eunit
62 hpstrg(2)=0.1d0*eunit
63 hpstrg(3)=0.1d0*eunit
64 hpstrg(4)=0.1d0*eunit
65 hpstrg(5)=0.9d0*eunit
66 hpstrg(6)=2.8d0*eunit
67 hpstrg(7)=0.1d0*eunit
68 hpstrg(8)=2.8d0*eunit
69 hpstrg(9)=3.2d0*eunit
70
71 do i=1,mxrs
72 do j=1,6
73 ihpat(i,j)=0
74 enddo
75 nhpat(i)=0
76 enddo
77 do i=irsml1(1),irsml2(ntlml)
78 mynm=seq(i)
79 call tolost(mynm)
80 if ((mynm.eq.'pro').or.(mynm.eq.'cpro')
81 & .or.(mynm.eq.'cpru').or.(mynm.eq.'prou')
82 & .or.(mynm.eq.'pron').or.(mynm.eq.'pro+')) then
83 prlvr=.true. ! residue i is a proline variant
84 print *, 'proline variant ',mynm,i
85 else
86 prlvr=.false.
87 endif
88
89 if (mynm.eq.'ala') then
90 nhpat(i)=1
91 ihpat(i,1)=iCa(i)+2
92 else if (mynm.eq.'val') then
93 nhpat(i)=3
94 ihpat(i,1)=iCa(i)+2
95 ihpat(i,2)=iCa(i)+4
96 ihpat(i,3)=iCa(i)+8
97 else if (prlvr) then
98 nhpat(i)=3
99 ihpat(i,1)=iCa(i)+2
100 ihpat(i,2)=iCa(i)+5
101 ihpat(i,3)=iCa(i)+8
102 else if (mynm.eq.'leu') then
103 nhpat(i)=4
104 ihpat(i,1)=iCa(i)+2
105 ihpat(i,2)=iCa(i)+5
106 ihpat(i,3)=iCa(i)+7
107 ihpat(i,4)=iCa(i)+11
108 else if (mynm.eq.'ile') then
109 nhpat(i)=4
110 ihpat(i,1)=iCa(i)+2
111 ihpat(i,2)=iCa(i)+4
112 ihpat(i,3)=iCa(i)+8
113 ihpat(i,4)=iCa(i)+11
114 else if (mynm.eq.'met') then
115 nhpat(i)=4
116 ihpat(i,1)=iCa(i)+2
117 ihpat(i,2)=iCa(i)+5
118 ihpat(i,3)=iCa(i)+8
119 ihpat(i,4)=iCa(i)+9
120 else if (mynm.eq.'phe') then
121 nhpat(i)=6
122 ihpat(i,1)=iCa(i)+5
123 ihpat(i,2)=iCa(i)+6
124 ihpat(i,3)=iCa(i)+8
125 ihpat(i,4)=iCa(i)+10
126 ihpat(i,5)=iCa(i)+12
127 ihpat(i,6)=iCa(i)+14
128 else if (mynm.eq.'tyr') then
129 nhpat(i)=6
130 ihpat(i,1)=iCa(i)+5
131 ihpat(i,2)=iCa(i)+6
132 ihpat(i,3)=iCa(i)+8
133 ihpat(i,4)=iCa(i)+10
134 ihpat(i,5)=iCa(i)+13
135 ihpat(i,6)=iCa(i)+15
136 else if (mynm.eq.'trp') then
137 nhpat(i)=6
138 ihpat(i,1)=iCa(i)+10
139 ihpat(i,2)=iCa(i)+11
140 ihpat(i,3)=iCa(i)+13
141 ihpat(i,4)=iCa(i)+15
142 ihpat(i,5)=iCa(i)+17
143 ihpat(i,6)=iCa(i)+19
144 endif
145 enddo
146! print *,'Hydrophobicity'
147
148! Excluded volume and local pair excluded volume terms
149 exvk=0.1*eunit
150 exvcut=4.3
151 exvcut2=exvcut*exvcut
152
153 sigsa(1)=1.0d0 ! hydrogen
154 sigsa(2)=1.0d0 ! hydrogen
155 sigsa(3)=1.0d0 ! hydrogen
156 sigsa(4)=1.0d0 ! hydrogen
157 sigsa(5)=1.0d0 ! hydrogen
158 sigsa(6)=1.0d0 ! hydrogen
159 sigsa(7)=1.75d0 ! carbon
160 sigsa(8)=1.75d0 ! carbon
161 sigsa(9)=1.75d0 ! carbon
162 sigsa(10)=1.42d0 ! oxygen
163 sigsa(11)=1.42d0 ! oxygen
164 sigsa(12)=1.42d0 ! oxygen
165 sigsa(13)=1.55d0 ! nitrogen
166 sigsa(14)=1.55d0 ! nitrogen
167 sigsa(15)=1.55d0 ! nitrogen
168 sigsa(16)=1.77d0 ! sulfur
169 sigsa(17)=1.0d0 ! hydrogen
170 sigsa(18)=1.75d0 ! carbon
171
172 do i=1,mxtyat
173 do j=1,mxtyat
174 sig2lcp(i,j)=(sigsa(i)+sigsa(j))**2
175 asalcp(i,j)=-7*((sig2lcp(i,j)/exvcut2)**6.0d0)
176 bsalcp(i,j)=6*((sig2lcp(i,j)/exvcut2)**6.0d0)/exvcut2
177 enddo
178 enddo
179! print *,'Local pair excluded volume constants'
180
181 exvlam=0.75d0
182 exvcutg=exvcut*exvlam
183 exvcutg2=exvcutg*exvcutg
184
185 do i=1,mxtyat
186 do j=1,mxtyat
187 sig2exv(i,j)=(exvlam*(sigsa(i)+sigsa(j)))**2
188 asaexv(i,j)=-7*((sig2exv(i,j)/exvcutg2)**6.0)
189 bsaexv(i,j)=6*((sig2exv(i,j)/exvcutg2)**6.0)/exvcutg2
190 enddo
191 enddo
192! print *,'General excluded volume constants'
193
194! Initialization of the connections matrix matcon(i,j). The index
195! i runs from -mxconr to +mxconr, and j from 1 to mxat.
196! matcon(i2-i1,i1) = 0, if the distance between atoms i1 and i2 is fixed
197! = 2, if atoms i1 and i2 are separated by 3 covalent
198! bonds and their distance can change
199! = 1, for all other pairs
200! if abs(i2-i1) > mxconr, the atoms are assumed to be separated by
201! many bonds, and with no restriction on their distances. On a protein
202! molecule made of natural amino acids, atoms with indices separated
203! by more than 35 can not be connected by three covalent bonds.
204
205 do i=1,mxat
206 do j=-mxconr,mxconr
207 matcon(j,i)=1
208 enddo
209 matcon(0,i)=0
210 enddo
211! continued...
212 do iml=1,ntlml
213 do iat1=iatrs1(irsml1(iml)),iatrs2(irsml2(iml))
214 do j=1,nbdat(iat1)
215 iat2=ibdat(j,iat1)
216 matcon(iat2-iat1,iat1)=0 ! directly bonded atoms
217 matcon(iat1-iat2,iat2)=0 !
218 do k=1,nbdat(iat2)
219 iat3=ibdat(k,iat2)
220 matcon(iat3-iat1,iat1)=0 !
221 matcon(iat1-iat3,iat3)=0 ! 2 covalent bonds
222 do l=1,nbdat(iat3)
223 iat4=ibdat(l,iat3)
224 if (iat4.ne.iat2) then
225 matcon(iat4-iat1,iat1)=2 !
226 matcon(iat1-iat4,iat4)=2 ! 3 covalent bonds
227 endif
228 do m=1,nbdat(iat2)
229 iat3p=ibdat(m,iat2)
230 if (iat3p.ne.iat3) then
231 matcon(iat4-iat3p,iat3p)=2 ! 3 covalent bonds
232 matcon(iat3p-iat4,iat4)=2 !
233 endif
234 enddo
235 enddo
236 enddo
237 do k=1,nbdat(iat1)
238 iat3=ibdat(k,iat1)
239 matcon(iat3-iat2,iat2)=0 ! also 2 bonds iat2-iat1-iat3
240 matcon(iat2-iat3,iat3)=0 !
241 enddo
242 enddo
243 enddo
244
245! print *,'going to initialize connections for first residue'
246! print *,'iN,iCa,iC =',iN(irsml1(iml)),
247! # iCa(irsml1(iml)),iC(irsml1(iml))
248 do iat1=iN(irsml1(iml))+1,iCa(irsml1(iml))-1
249! print *,'connections for iat1 = ',iat1
250 matcon(iat1-iN(irsml1(iml)),iN(irsml1(iml)))=0
251 matcon(iN(irsml1(iml))-iat1,iat1)=0
252 matcon(iat1-iCa(irsml1(iml)),iCa(irsml1(iml)))=0
253 matcon(iCa(irsml1(iml))-iat1,iat1)=0
254 do iat2=iat1,iCa(irsml1(iml))-1
255 matcon(iat2-iat1,iat1)=0
256 matcon(iat1-iat2,iat2)=0
257 enddo
258 matcon(iat1-iCa(irsml1(iml))-1,iCa(irsml1(iml))+1)=2
259 matcon(iCa(irsml1(iml))+1-iat1,iat1)=2
260 matcon(iat1-iCa(irsml1(iml))-2,iCa(irsml1(iml))+2)=2
261 matcon(iCa(irsml1(iml))+2-iat1,iat1)=2
262 matcon(iat1-iC(irsml1(iml)),iC(irsml1(iml)))=2
263 matcon(iC(irsml1(iml))-iat1,iat1)=2
264 enddo
265
266! Below: for certain residues, some atoms separated by 3 or more bonds
267! do not change distance. So, the connection matrix term for such pairs
268! should be zero.
269
270 do irs=irsml1(iml),irsml2(iml)
271 if (irs.eq.irsml1(iml)) then
272 iatoff=1
273 else
274 iatoff=0
275 endif
276 if (irs.eq.irsml2(iml)) then
277 iatmrg=2
278 else
279 iatmrg=0
280 endif
281 mynm=seq(irs)
282 call tolost(mynm)
283 if ((mynm.eq.'pro')) then
284 prlvr=.true. ! residue i is a proline variant
285 else
286 prlvr=.false.
287 endif
288 if ((mynm.eq.'asn')) then
289 call set_local_constr(irs,5,9,2)
290 else if ((mynm.eq.'gln')) then
291 call set_local_constr(irs,8,12,5)
292 else if ((mynm.eq.'arg')) then
293 call set_local_constr(irs,11,19,8)
294 else if ((mynm.eq.'his')) then
295 call set_local_constr(irs,5,12,2)
296 else if (mynm.eq.'phe') then
297 call set_local_constr(irs,5,15,2)
298 else if (mynm.eq.'tyr') then
299 call set_local_constr(irs,5,16,2)
300 iat1=iCa(irs)+12 ! H_h
301 iat2=iCa(irs)+13 ! C_e2
302 iat3=iCa(irs)+8 ! C_e1
303 matcon(iat2-iat1,iat1)=2
304 matcon(iat1-iat2,iat2)=2
305 matcon(iat3-iat1,iat1)=2
306 matcon(iat1-iat3,iat3)=2
307 else if (mynm.eq.'trp') then
308 call set_local_constr(irs,5,19,2)
309 else if (prlvr) then
310! Proline. Many more distances are fixed because of the fixed
311! phi angle
312 call set_local_constr(irs,-3,11,-3)
313 matcon(iCa(irs-1)-iCa(irs)-8,iCa(irs)+8)=0
314 matcon(iCa(irs)+8-iCa(irs-1),iCa(irs-1))=0
315 endif
316 enddo
317
318! Distances fixed because of the constant omega angle
319 do irs=irsml1(iml),irsml2(iml)-1
320 matcon(iCa(irs+1)-iC(irs)-1,iC(irs)+1)=0 ! O_i -- Ca_{i+1}
321 matcon(-iCa(irs+1)+iC(irs)+1,iCa(irs+1))=0 ! O_i -- Ca_{i+1}
322
323 matcon(iN(irs+1)-iC(irs),iC(irs)+1)=0 ! O_i -- H_{i+1}
324 matcon(-iN(irs+1)+iC(irs),iN(irs+1)+1)=0 ! O_i -- H_{i+1}
325
326 matcon(iCa(irs+1)-iCa(irs),iCa(irs))=0 ! Ca_i -- Ca_{i+1}
327 matcon(-iCa(irs+1)+iCa(irs),iCa(irs+1))=0 ! Ca_i -- Ca_{i+1}
328
329 matcon(iN(irs+1)+1-iCa(irs),iCa(irs))=0 ! Ca_i -- H_{i+1}
330 matcon(-iN(irs+1)-1+iCa(irs),iN(irs+1)+1)=0 ! Ca_i -- H_{i+1}
331 enddo
332
333 enddo
334! finished initializing matrix conmat
335! print *,'Connections matrix'
336! Local pair excluded volume
337 do i=1,mxml
338 ilpst(i)=1
339 ilpnd(i)=0
340 enddo
341 do i=1,50*mxrs
342 lcp1(i)=0
343 lcp2(i)=0
344 enddo
345 ilp=0
346 do iml=1,ntlml
347 do iat1=iatrs1(irsml1(iml)),iatrs2(irsml2(iml))
348 do iat2=iat1+1,iatrs2(irsml2(iml))
349 if ((iat2-iat1.le.mxconr).and.
350 & matcon(iat2-iat1,iat1).eq.2) then
351 ilp=ilp+1
352 lcp1(ilp)=iat1
353 lcp2(ilp)=iat2
354 endif
355 enddo
356 enddo
357
358 ilpnd(iml)=ilp
359 if (iml.lt.ntlml) then
360 ilpst(iml+1)=ilp+1
361 endif
362 print *,'molecule ',iml,' lc pair range ',ilpst(iml),ilpnd(iml)
363! print *,'local pair list'
364 do lci=ilpst(iml),ilpnd(iml)
365 iat1=lcp1(lci)
366 iat2=lcp2(lci)
367! print *,lci,iat1,iat2,ityat(iat1),ityat(iat2)
368 enddo
369 enddo
370
371 print *,'finished initializing Lund force field'
372 end
373
374 integer function ihptype(iaa)
375 include 'INCL.H'
376 character mynm*4
377 mynm=seq(iaa)
378 ityp=-1
379 call tolost(mynm)
380 if (mynm.eq.'ala') then
381 ityp=1
382 else if ((mynm.eq.'val').or.(mynm.eq.'leu').or.(mynm.eq.'ile')
383 & .or.(mynm.eq.'met').or.(mynm.eq.'pro')) then
384 ityp=2
385 else if ((mynm.eq.'phe').or.(mynm.eq.'tyr').or.(mynm.eq.'trp'))
386 & then
387 ityp=3
388 endif
389 ihptype=ityp
390 return
391 end
392
393 real*8 function ebiasrs(irsd)
394 include 'INCL.H'
395 include 'incl_lund.h'
396 dimension q1(2),q2(2)
397 data q1/-0.2,0.2/
398 data q2/0.42,-0.42/
399
400 et=0.0
401 do i=0,1
402 iat1=iN(irsd)+i
403 do j=0,1
404 iat2=iC(irsd)+j
405 xij=xat(iat1)-xat(iat2)
406 yij=yat(iat1)-yat(iat2)
407 zij=zat(iat1)-zat(iat2)
408 et=et+q1(i+1)*q2(j+1)/sqrt(xij*xij+yij*yij+zij*zij)
409 enddo
410 enddo
411 ebiasrs=kbias*et
412 return
413 end
414! Evaluates backbone backbone hydrogen bond strength for residues
415! i and j, taking the donor from residue i and acceptor from residue j
416 real*8 function ehbmmrs(i,j)
417 include 'INCL.H'
418 include 'incl_lund.h'
419 double precision r2,r4,r6,dx,dy,dz,ca,cb
420 integer d1,d2,a1,a2
421 d1=iN(i)
422 d2=d1+1
423 a1=iC(j)+1 ! dipoles are numbered from -ve to +ve charge
424 a2=iC(j) ! so, acceptor 1 is the Oxygen, and a2, the carbon
425 dx=xat(a1)-xat(d2)
426 dy=yat(a1)-yat(d2)
427 dz=zat(a1)-zat(d2)
428 r2=dx*dx+dy*dy+dz*dz
429 if (r2.gt.cthb2) then
430! print *,'hbmm = 0 ',cthb2,r2,a1,a2,d1,d2
431! print *,'a1,a2,d1,d2,r2 = ',a1,a2,d1,d2,r2,sighb2,cthb
432 ehbmmrs=0
433 return
434 endif
435 ca=(xat(d2)-xat(d1))*dx+(yat(d2)-yat(d1))*dy+(zat(d2)-zat(d1))*dz
436 cb=(xat(a2)-xat(a1))*dx+(yat(a2)-yat(a1))*dy+(zat(a2)-zat(a1))*dz
437 if (powa.gt.0.and.ca.le.0) then
438! print *,'hbmm, returning 0 because of angle a'
439 ehbmmrs=0
440 return
441 endif
442 if (powb.gt.0.and.cb.le.0) then
443! print *,'hbmm, returning 0 because of angle b'
444 ehbmmrs=0
445 return
446 endif
447 r6=sighb2/r2
448 r4=r6*r6
449 r6=r6*r4
450 rdon2=(xat(d2)-xat(d1))**2+(yat(d2)-yat(d1))**2+
451 & (zat(d2)-zat(d1))**2
452 racc2=(xat(a2)-xat(a1))**2+(yat(a2)-yat(a1))**2+
453 & (zat(a2)-zat(a1))**2
454 evlu=((ca*ca/(r2*rdon2))**(0.5*powa))*
455 & ((cb*cb/(r2*racc2))**(0.5*powb))
456 evlu=evlu*(r6*(5*r6-6*r4)+alhb+blhb*r2)
457! print *,'found hbmm contribution ',i,j,epshb1,evlu
458 ehbmmrs=epshb1*evlu
459 return
460 end
461 real*8 function enylun(nml)
462! nml = 1 .. ntlml. No provision exists to handle out of range values
463! for nml inside this function.
464 include 'INCL.H'
465 include 'incl_lund.h'
466 character mynm*4
467 logical prlvr
468 eyhb=0.0 ! backbone-backbone and sidechain-backbone HB
469 eyel=0.0 ! Bias, or local electrostatic term
470 eysl=0.0 ! Hydrophobicity, replaces the solvent term of SMMP
471 eyvr=0.0 ! Local pair excluded volume, in a sense a variable potential
472 eyvw=0.0 ! atom-atom repulsion, excluded volume
473! atom-atom repulsion is calculated on a system wide basis, instead of
474! molecule by molecule for efficiency. Look into function exvlun.
475
476 istres=irsml1(nml)
477 indres=irsml2(nml)
478
479! First, all terms that can be calculated on a residue by residue basis
480 do i=istres,indres
481 mynm=seq(i)
482 call tolost(mynm)
483 if ((mynm.eq.'pro').or.(mynm.eq.'cpro')
484 & .or.(mynm.eq.'cpru').or.(mynm.eq.'prou')
485 & .or.(mynm.eq.'pron').or.(mynm.eq.'pro+')) then
486 prlvr=.true. ! residue i is a proline variant
487 else
488 prlvr=.false.
489 endif
490
491! Bias, or local electrostatic term. Excluded from the list are
492! residues at the ends of the chain, glycine and all proline variants
493 if ((i.ne.istres).and.(i.ne.indres).and.
494 & .not.prlvr.and.mynm.ne.'gly') then
495 eyel=eyel+ebiasrs(i)
496 endif
497! Backbone--backbone hydrogen bonds
498 shbm1=1.0
499 shbm2=1.0
500 if ((i.eq.istres).or.(i.eq.indres)) shbm1=0.5
501! Residue i contributes the donor, and j, the acceptor, so both i and
502! j run over the whole set of amino acids.
503! No terms for residue i, if it is a proline variant.
504 if (.not.prlvr) then
505 do j=istres,indres
506 if ((j.lt.(i-2)).or.(j.gt.(i+1))) then
507 shbm2=1.0
508 if ((j.eq.istres).or.(j.eq.indres)) shbm2=0.5
509 etmp=ehbmmrs(i,j)
510 eyhb=eyhb+shbm1*shbm2*etmp
511 endif
512 enddo
513 endif
514! Hydrophobicity, only if residue i is hydrophobic to start with
515 ihpi=ihptype(i)
516 if (ihpi.ge.0) then
517! Unlike hydrogen bonds, the hydrophobicity potential is symmetric
518! in i and j. So, the loop for j runs from i+1 to the end.
519
520 do j=i+1,indres
521 ihpj=ihptype(j)
522 if (ihpj.ge.0) then
523 etmp=ehp(i,j,ihpi,ihpj)
524 if (j.eq.(i+1)) etmp=0
525 if (j.eq.(i+2)) etmp=0.5*etmp
526! print *, 'hydrophobicity contribution ',etmp,i,j
527 eysl=eysl+etmp
528 endif
529 enddo
530 endif
531 enddo
532
533! Terms that are not calculated residue by residue ...
534
535! Local pair or third-neighbour excluded volume
536! Numerically this is normally the term with largest positive
537! contribution to the energy in an equilibrated stystem.
538
539 i1=ilpst(nml)
540 i2=ilpnd(nml)
541 etmp=0.0
542 do i=i1,i2
543 etmp1=0.0
544 iat1=lcp1(i)
545 iat2=lcp2(i)
546 iatt1=ityat(iat1)
547 iatt2=ityat(iat2)
548 xij=xat(iat1)-xat(iat2)
549 yij=yat(iat1)-yat(iat2)
550 zij=zat(iat1)-zat(iat2)
551 r2=(xij*xij + yij*yij + zij*zij)
552 if (r2.le.exvcut2) then
553 r6=sig2lcp(iatt1,iatt2)/r2
554 r6=r6*r6*r6
555 etmp1=(r6*r6+asalcp(iatt1,iatt2)+bsalcp(iatt1,iatt2)*r2)
556 if (etmp1.ge.0) then
557! print *,'local pair contribution ',iat1,iat2,iatt1,
558! & iatt2,sqrt(sig2lcp(iatt1,iatt2)),etmp1
559 endif
560 etmp=etmp+etmp1
561 endif
562! print *,'pair : ',iat1,iat2,' contribution ',etmp1
563! print *,exvcut2,r2
564 enddo
565 eyvr=exvk*etmp
566 eysm=eyel+eyhb+eyvr+eysl
567 enylun=eysm
568 end
569 real*8 function eninlun()
570 include 'INCL.H'
571 eysmi=0.0
572 eyeli=0.0
573 eyvwi=0.0
574 eyhbi=0.0
575 eninlun=0.0
576 return
577 end
578 real*8 function ehp(i1,i2,ihp1,ihp2)
579 include 'INCL.H'
580 include 'incl_lund.h'
581 dimension r2min(12)
582
583 b2=20.25
584 a2=12.25
585! ihp1=ihptype(i1)
586! ihp2=ihptype(i2)
587 if ((ihp1.le.0).or.(ihp2.le.0)) then
588 ehp=0.0
589 return
590 endif
591 ni=nhpat(i1)
592 nj=nhpat(i2)
593 do i=1,ni+nj
594 r2min(i)=100000000 ! quite far
595 enddo
596 do i=1,ni
597 k=ihpat(i1,i)
598 do j=1,nj
599 l=ihpat(i2,j)
600 xij=xat(k)-xat(l)
601 yij=yat(k)-yat(l)
602 zij=zat(k)-zat(l)
603 dtmp=(xij*xij + yij*yij + zij*zij)
604 if (dtmp.le.r2min(i)) r2min(i)=dtmp
605 if (dtmp.le.r2min(ni+j)) r2min(ni+j)=dtmp
606 enddo
607 enddo
608 ssum=0
609 do i=1,ni+nj
610 if (r2min(i).le.b2) then
611 if (r2min(i).lt.a2) then
612 ssum=ssum+1
613 else
614 ssum=ssum+(b2-r2min(i))/(b2-a2)
615 endif
616 endif
617! hpstrgth=hpstrg((ihp1-1)*nhptyp+ihp2)
618! print *, 'hp diagnosis ',ni,nj,ssum/(ni+nj),hpstrgth
619 enddo
620 ehp=-hpstrg((ihp1-1)*nhptyp+ihp2)*ssum/(ni+nj)
621 return
622 end
623
624 real*8 function exvlun(nml)
625 include 'INCL.H'
626 include 'incl_lund.h'
627! For multi-chain systems it makes little sense to split the calculation
628! of this term into an 'interaction part' and a contribution from
629! individual molecules. So, normally this should always be called with
630! argument nml=0. Only for diagnostic reasons, you might want to find
631! the contribution from one molecule in a multi-chain system assuming
632! there was no other molecule.
633 dimension isort(mxat),ngbr(mxat),locccl(mxat),incell(mxcell)
634 dimension icell(mxat)
635 integer :: jx1, jy1, jz1
636
637 if (nml.eq.0) then
638 istat=iatrs1(irsml1(1))
639 indat=iatrs2(irsml2(ntlml))
640 else
641 istat=iatrs1(irsml1(nml))
642 indat=iatrs2(irsml2(nml))
643 endif
644
645 eyvw=0.0
646! The beginning part of this implementation is very similar to the
647! assignment of cells to the atoms during calculation of solvent
648! accessible surface area. So, much of that part is similar. But
649! unlike the accessible surface calculations, this term is symmetric
650! in any two participating atoms. So, the part after the assignment
651! of cells differs even in the structure of the code.
652
653 do i=1,mxcell
654 incell(i)=0
655 enddo
656! print *,'evaluating general excluded volume :',istat,',',indat
657! Find minimal containing box
658 xmin=xat(istat)
659 ymin=yat(istat)
660 zmin=zat(istat)
661 xmax=xmin
662 ymax=ymin
663 zmax=zmin
664
665 do i=istat,indat
666 if (xat(i).le.xmin) then
667 xmin=xat(i)
668 else if (xat(i).ge.xmax) then
669 xmax=xat(i)
670 endif
671 if (yat(i).le.ymin) then
672 ymin=yat(i)
673 else if (yat(i).ge.ymax) then
674 ymax=yat(i)
675 endif
676 if (zat(i).le.zmin) then
677 zmin=zat(i)
678 else if (zat(i).ge.zmax) then
679 zmax=zat(i)
680 endif
681 enddo
682
683 sizex=xmax-xmin
684 sizey=ymax-ymin
685 sizez=zmax-zmin
686! Number of cells along each directions that fit into the box.
687 ndx=int(sizex/exvcutg)+1
688 ndy=int(sizey/exvcutg)+1
689 ndz=int(sizez/exvcutg)+1
690
691 nxy=ndx*ndy
692 ncell=nxy*ndz
693! print *,'Number of cells along x,y,z = ',ndx,',',ndy,',',ndz
694 if (ncell.ge.mxcell) then
695 print *,'exvlun> required number of cells',ncell,
696 & ' exceeded the limit ',mxcell
697 print *,'recompile with a higher mxcell.'
698 stop
699 endif
700! Expand box to contain an integral number of cells along each direction
701 shiftx=(dble(ndx)*exvcutg-sizex)/2.0
702 shifty=(dble(ndy)*exvcutg-sizey)/2.0
703 shiftz=(dble(ndz)*exvcutg-sizez)/2.0
704 xmin=xmin-shiftx
705 ymin=ymin-shifty
706 zmin=zmin-shiftz
707 xmax=xmax+shiftx
708 ymax=ymax+shifty
709 zmax=zmax+shiftz
710
711! Set occupied cells to zero. Note that the maximum number of occupied
712! cells is the number of atoms in the system.
713 nocccl=0
714 do i=1,mxat
715 locccl(i)=0
716 enddo
717
718! Put atoms in cells
719 do j=istat,indat
720 mx=min(int(max((xat(j)-xmin)/exvcutg,0.0d0)),ndx-1)
721 my=min(int(max((yat(j)-ymin)/exvcutg,0.0d0)),ndy-1)
722 mz=min(int(max((zat(j)-zmin)/exvcutg,0.0d0)),ndz-1)
723 icellj=mx+my*ndx+mz*nxy+1
724 icell(j)=icellj
725 if (icellj.gt.mxcell) then
726 print *,'exvlun> bad cell index ',icellj,' for atom ',j
727 stop
728 else
729 if (incell(icellj).eq.0) then
730! previously unoccupied cell
731 nocccl=nocccl+1
732 locccl(nocccl)=icellj
733 endif
734 incell(icellj)=incell(icellj)+1
735 endif
736 enddo
737! print *,'finished assigning cells. nocccl = ',nocccl
738! Cummulative occupancy of i'th cell
739 do i=1,ncell
740 incell(i+1)=incell(i+1)+incell(i)
741 enddo
742! print *,'finished making cumulative cell sums'
743! Sorting atoms by their cell index
744 do i=istat,indat
745 j=icell(i)
746 jj=incell(j)
747 isort(jj)=i
748 incell(j)=jj-1
749 enddo
750! print *,'sorted atoms by cell index'
751 etmp=0.0
752 do icl=1,nocccl
753! loop through occupied cells
754 lcell=locccl(icl)
755 ix=mod(lcell-1,ndx)
756 iz = (lcell - 1) / nxy
757 iy = ((lcell - 1) - iz * nxy) / ndx
758
759c print *,'icl=',icl,'absolute index of cell = ',lcell
760c print *,'iz,iy,ix = ',iz,iy,ix
761c find all atoms in current cell and all its forward-going neighbours
762 nex=ix+1
763 ney=iy+1
764 nez=iz+1
765 nsame=0
766 nngbr=0
767 do jx=ix,nex
768 if (jx.ge.ndx) then
769 jx1=0
770 else
771 jx1=jx
772 endif
773 do jy=iy,ney
774 if (jy.ge.ndy) then
775 jy1=0
776 else
777 jy1=jy
778 endif
779 do jz=iz,nez
780 if (jz.ge.ndz) then
781 jz1=0
782 else
783 jz1=jz
784 endif
785 jcl=jx1 + ndx*jy1 + nxy*jz1 + 1
786! write(*,*)'jcl,jx1,jy1,jz1:', jcl,jx1,jy1,jz1
787 do ii=incell(jcl)+1,incell(jcl+1)
788! count the total number of neighbours
789 nngbr=nngbr+1
790 if (jx.eq.ix.and.jy.eq.iy.and.jz.eq.iz) then
791! count how many neighbours are from the same cell
792 nsame=nsame+1
793 endif
794 ngbr(nngbr)=isort(ii)
795 enddo
796 enddo
797 enddo
798 enddo
799! A few more cells need to be searched, so that we cover 13 of the 26
800! neighbouring cells.
801! 1
802 jx=ix+1
803 if (jx.ge.ndx) jx=0
804 jy=iy
805 jz=iz-1
806 if (jz.lt.0) jz=ndz-1
807 jcl=jx+ndx*jy+nxy*jz+1
808 do ii=incell(jcl)+1,incell(jcl+1)
809 nngbr=nngbr+1
810 ngbr(nngbr)=isort(ii)
811 enddo
812! 2
813 jx=ix
814 jy=iy-1
815 if (jy.lt.0) jy =ndy-1
816 jz=iz+1
817 if (jz.ge.ndz) jz=0
818 jcl=jx+ndx*jy+nxy*jz+1
819 do ii=incell(jcl)+1,incell(jcl+1)
820 nngbr=nngbr+1
821 ngbr(nngbr)=isort(ii)
822 enddo
823! 3
824 jx=ix-1
825 if (jx.lt.0)jx =ndx-1
826 jy=iy+1
827 if (jy.ge.ndy) jy=0
828 jz=iz
829 jcl=jx+ndx*jy+nxy*jz+1
830 do ii=incell(jcl)+1,incell(jcl+1)
831 nngbr=nngbr+1
832 ngbr(nngbr)=isort(ii)
833 enddo
834! 4
835 jx=ix+1
836 if (jx.ge.ndx) jx=0
837 jy=iy+1
838 if (jy.ge.ndy) jy=0
839 jz=iz-1
840 if (jz.lt.0) jz=ndz-1
841 jcl=jx+ndx*jy+nxy*jz+1
842 do ii=incell(jcl)+1,incell(jcl+1)
843 nngbr=nngbr+1
844 ngbr(nngbr)=isort(ii)
845 enddo
846! 5
847 jx=ix+1
848 if (jx.ge.ndx) jx = 0
849 jy=iy-1
850 if (jy.lt.0) jy=ndy-1
851 jz=iz+1
852 if (jz.ge.ndz) jz=0
853 jcl=jx+ndx*jy+nxy*jz+1
854 do ii=incell(jcl)+1,incell(jcl+1)
855 nngbr=nngbr+1
856 ngbr(nngbr)=isort(ii)
857 enddo
858! 6
859 jx=ix+1
860 if (jx.ge.ndx) jx=0
861 jy=iy-1
862 if (jy.lt.0) jy=ndy-1
863 jz=iz-1
864 if (jz.lt.0) jz=ndy-1
865 jcl=jx+ndx*jy+nxy*jz+1
866 do ii=incell(jcl)+1,incell(jcl+1)
867 nngbr=nngbr+1
868 ngbr(nngbr)=isort(ii)
869 enddo
870
871! print *,'atoms in same cell ',nsame
872! print *,'atoms in neighbouring cells ',nngbr
873 do i1=1,nsame
874! Over all atoms from the original cell
875 iat1=ngbr(i1)
876 do i2=i1,nngbr
877! Over all atoms in the original+neighbouring cells
878 iat2=ngbr(i2)
879 xij=xat(iat1)-xat(iat2)
880 yij=yat(iat1)-yat(iat2)
881 zij=zat(iat1)-zat(iat2)
882 r2=(xij*xij+yij*yij+zij*zij)
883
884 if (r2.le.exvcutg2) then
885 if (abs(iat2-iat1).gt.mxconr.or.
886 & matcon(iat2-iat1,iat1).eq.1) then
887 iatt1=ityat(iat1)
888 iatt2=ityat(iat2)
889 r6=sig2exv(iatt1,iatt2)/r2
890 r6=r6*r6*r6
891 etmp1=r6*r6+asaexv(iatt1,iatt2)
892 & +bsaexv(iatt1,iatt2)*r2
893 etmp=etmp+etmp1
894! if (etmp1.ge.2000) then
895! print *,'contribution ',iat1,iat2,etmp1
896! call outpdb(1,'EXAMPLES/clash.pdb')
897! stop
898! endif
899 endif
900 endif
901 enddo
902 enddo
903 enddo
904! irs=1
905! do iat=iatrs1(irs),iatrs2(irs)
906! do j=-mxconr,mxconr
907! print *,iat,j,':',matcon(j,iat)
908! enddo
909! enddo
910! irs=irsml2(1)
911! do iat=iatrs1(irs),iatrs2(irs)
912! do j=-mxconr,mxconr
913! print *,iat,j,':',matcon(j,iat)
914! enddo
915! enddo
916
917 eyvw=exvk*etmp
918 exvlun=eyvw
919 return
920 end
921
922 real*8 function exvbrfc()
923! Brute force excluded volume evaluation
924 include 'INCL.H'
925 include 'incl_lund.h'
926
927 etmp=0.0
928 etmp1=0.0
929 print *,'max connection radius is ',mxconr
930 do iat1=iatrs1(irsml1(1)),iatrs2(irsml2(ntlml))
931 do iat2=iat1+1,iatrs2(irsml2(ntlml))
932 xij=xat(iat1)-xat(iat2)
933 yij=yat(iat1)-yat(iat2)
934 zij=zat(iat1)-zat(iat2)
935 r2=(xij*xij+yij*yij+zij*zij)
936
937 if (r2.le.exvcutg2) then
938 if (abs(iat2-iat1).gt.mxconr.or.
939 & matcon(iat2-iat1,iat1).eq.1) then
940 iatt1=ityat(iat1)
941 iatt2=ityat(iat2)
942 r6=sig2exv(iatt1,iatt2)/r2
943 r6=r6*r6*r6
944 etmp1=r6*r6+asaexv(iatt1,iatt2)
945 & +bsaexv(iatt1,iatt2)*r2
946 etmp=etmp+etmp1
947 else
948! print *,'atoms ', iat1,' and ',iat2,' were close',
949! # 'but matcon is ',matcon(iat2-iat1,iat1)
950 endif
951 endif
952 enddo
953 enddo
954 exvbrfc=etmp*exvk
955 return
956 end
Note: See TracBrowser for help on using the repository browser.