dynamic dissulfides with finegrain parallelization cont.
[unres.git] / source / unres / src_MD / ssMD.F
1 c----------------------------------------------------------------------------
2       subroutine check_energies
3 c      implicit none
4
5 c     Includes
6       include 'DIMENSIONS'
7       include 'COMMON.CHAIN'
8       include 'COMMON.VAR'
9       include 'COMMON.IOUNITS'
10       include 'COMMON.SBRIDGE'
11       include 'COMMON.LOCAL'
12       include 'COMMON.GEO'
13
14 c     External functions
15       double precision ran_number
16       external ran_number
17
18 c     Local variables
19       integer i,j,k,l,lmax,p,pmax
20       double precision rmin,rmax
21       double precision eij
22
23       double precision d
24       double precision wi,rij,tj,pj
25
26
27 c      return
28
29       i=5
30       j=14
31
32       d=dsc(1)
33       rmin=2.0D0
34       rmax=12.0D0
35
36       lmax=10000
37       pmax=1
38
39       do k=1,3
40         c(k,i)=0.0D0
41         c(k,j)=0.0D0
42         c(k,nres+i)=0.0D0
43         c(k,nres+j)=0.0D0
44       enddo
45
46       do l=1,lmax
47
48 ct        wi=ran_number(0.0D0,pi)
49 c        wi=ran_number(0.0D0,pi/6.0D0)
50 c        wi=0.0D0
51 ct        tj=ran_number(0.0D0,pi)
52 ct        pj=ran_number(0.0D0,pi)
53 c        pj=ran_number(0.0D0,pi/6.0D0)
54 c        pj=0.0D0
55
56         do p=1,pmax
57 ct           rij=ran_number(rmin,rmax)
58
59            c(1,j)=d*sin(pj)*cos(tj)
60            c(2,j)=d*sin(pj)*sin(tj)
61            c(3,j)=d*cos(pj)
62
63            c(3,nres+i)=-rij
64
65            c(1,i)=d*sin(wi)
66            c(3,i)=-rij-d*cos(wi)
67
68            do k=1,3
69               dc(k,nres+i)=c(k,nres+i)-c(k,i)
70               dc_norm(k,nres+i)=dc(k,nres+i)/d
71               dc(k,nres+j)=c(k,nres+j)-c(k,j)
72               dc_norm(k,nres+j)=dc(k,nres+j)/d
73            enddo
74
75            call dyn_ssbond_ene(i,j,eij)
76         enddo
77       enddo
78
79       call exit(1)
80
81       return
82       end
83
84 C-----------------------------------------------------------------------------
85
86       subroutine dyn_ssbond_ene(resi,resj,eij)
87 c      implicit none
88
89 c     Includes
90       include 'DIMENSIONS'
91       include 'COMMON.SBRIDGE'
92       include 'COMMON.CHAIN'
93       include 'COMMON.DERIV'
94       include 'COMMON.LOCAL'
95       include 'COMMON.INTERACT'
96       include 'COMMON.VAR'
97       include 'COMMON.IOUNITS'
98       include 'COMMON.CALC'
99 #ifndef CLUST
100 #ifndef WHAM
101       include 'COMMON.MD'
102 #endif
103 #endif
104
105 c     External functions
106       double precision h_base
107       external h_base
108
109 c     Input arguments
110       integer resi,resj
111
112 c     Output arguments
113       double precision eij
114
115 c     Local variables
116       logical havebond
117 c      integer itypi,itypj,k,l
118       double precision rrij,ssd,deltat1,deltat2,deltat12,cosphi
119       double precision sig0ij,ljd,sig,fac,e1,e2
120       double precision dcosom1(3),dcosom2(3),ed
121       double precision pom1,pom2
122       double precision ljA,ljB,ljXs
123       double precision d_ljB(1:3)
124       double precision ssA,ssB,ssC,ssXs
125       double precision ssxm,ljxm,ssm,ljm
126       double precision d_ssxm(1:3),d_ljxm(1:3),d_ssm(1:3),d_ljm(1:3)
127       double precision f1,f2,h1,h2,hd1,hd2
128       double precision omega,delta_inv,deltasq_inv,fac1,fac2
129 c-------FIRST METHOD
130       double precision xm,d_xm(1:3)
131 c-------END FIRST METHOD
132 c-------SECOND METHOD
133 c$$$      double precision ss,d_ss(0:3),ljf,d_ljf(0:3)
134 c-------END SECOND METHOD
135
136 c-------TESTING CODE
137       logical checkstop,transgrad
138       common /sschecks/ checkstop,transgrad
139
140       integer icheck,nicheck,jcheck,njcheck
141       double precision echeck(-1:1),deps,ssx0,ljx0
142 c-------END TESTING CODE
143
144
145       i=resi
146       j=resj
147
148       itypi=itype(i)
149       dxi=dc_norm(1,nres+i)
150       dyi=dc_norm(2,nres+i)
151       dzi=dc_norm(3,nres+i)
152       dsci_inv=vbld_inv(i+nres)
153
154       itypj=itype(j)
155       xj=c(1,nres+j)-c(1,nres+i)
156       yj=c(2,nres+j)-c(2,nres+i)
157       zj=c(3,nres+j)-c(3,nres+i)
158       dxj=dc_norm(1,nres+j)
159       dyj=dc_norm(2,nres+j)
160       dzj=dc_norm(3,nres+j)
161       dscj_inv=vbld_inv(j+nres)
162
163       chi1=chi(itypi,itypj)
164       chi2=chi(itypj,itypi)
165       chi12=chi1*chi2
166       chip1=chip(itypi)
167       chip2=chip(itypj)
168       chip12=chip1*chip2
169       alf1=alp(itypi)
170       alf2=alp(itypj)
171       alf12=0.5D0*(alf1+alf2)
172
173       rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
174       rij=dsqrt(rrij)  ! sc_angular needs rij to really be the inverse
175 c     The following are set in sc_angular
176 c      erij(1)=xj*rij
177 c      erij(2)=yj*rij
178 c      erij(3)=zj*rij
179 c      om1=dxi*erij(1)+dyi*erij(2)+dzi*erij(3)
180 c      om2=dxj*erij(1)+dyj*erij(2)+dzj*erij(3)
181 c      om12=dxi*dxj+dyi*dyj+dzi*dzj
182       call sc_angular
183       rij=1.0D0/rij  ! Reset this so it makes sense
184
185       sig0ij=sigma(itypi,itypj)
186       sig=sig0ij*dsqrt(1.0D0/sigsq)
187
188       ljXs=sig-sig0ij
189       ljA=eps1*eps2rt**2*eps3rt**2
190       ljB=ljA*bb(itypi,itypj)
191       ljA=ljA*aa(itypi,itypj)
192       ljxm=ljXs+(-2.0D0*aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
193
194       ssXs=d0cm
195       deltat1=1.0d0-om1
196       deltat2=1.0d0+om2
197       deltat12=om2-om1+2.0d0
198       cosphi=om12-om1*om2
199       ssA=akcm
200       ssB=akct*deltat12
201       ssC=ss_depth
202      &     +akth*(deltat1*deltat1+deltat2*deltat2)
203      &     +v1ss*cosphi+v2ss*cosphi*cosphi+v3ss*cosphi*cosphi*cosphi
204       ssxm=ssXs-0.5D0*ssB/ssA
205
206 c-------TESTING CODE
207 c$$$c     Some extra output
208 c$$$      ssm=ssC-0.25D0*ssB*ssB/ssA
209 c$$$      ljm=-0.25D0*ljB*bb(itypi,itypj)/aa(itypi,itypj)
210 c$$$      ssx0=ssB*ssB-4.0d0*ssA*ssC
211 c$$$      if (ssx0.gt.0.0d0) then
212 c$$$        ssx0=ssXs+0.5d0*(-ssB+sqrt(ssx0))/ssA
213 c$$$      else
214 c$$$        ssx0=ssxm
215 c$$$      endif
216 c$$$      ljx0=ljXs+(-aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
217 c$$$      write(iout,'(a,4f8.2,2f15.2,3f6.2)')"SSENERGIES ",
218 c$$$     &     ssxm,ljxm,ssx0,ljx0,ssm,ljm,om1,om2,om12
219 c$$$      return
220 c-------END TESTING CODE
221
222 c-------TESTING CODE
223 c     Stop and plot energy and derivative as a function of distance
224       if (checkstop) then
225         ssm=ssC-0.25D0*ssB*ssB/ssA
226         ljm=-0.25D0*ljB*bb(itypi,itypj)/aa(itypi,itypj)
227         if (ssm.lt.ljm .and.
228      &       dabs(rij-0.5d0*(ssxm+ljxm)).lt.0.35d0*(ljxm-ssxm)) then
229           nicheck=1000
230           njcheck=1
231           deps=0.5d-7
232         else
233           checkstop=.false.
234         endif
235       endif
236       if (.not.checkstop) then
237         nicheck=0
238         njcheck=-1
239       endif
240
241       do icheck=0,nicheck
242       do jcheck=-1,njcheck
243       if (checkstop) rij=(ssxm-1.0d0)+
244      &       ((ljxm-ssxm+2.0d0)*icheck)/nicheck+jcheck*deps
245 c-------END TESTING CODE
246
247       if (rij.gt.ljxm) then
248         havebond=.false.
249         ljd=rij-ljXs
250         fac=(1.0D0/ljd)**expon
251         e1=fac*fac*aa(itypi,itypj)
252         e2=fac*bb(itypi,itypj)
253         eij=eps1*eps2rt*eps3rt*(e1+e2)
254         eps2der=eij*eps3rt
255         eps3der=eij*eps2rt
256         eij=eij*eps2rt*eps3rt
257
258         sigder=-sig/sigsq
259         e1=e1*eps1*eps2rt**2*eps3rt**2
260         ed=-expon*(e1+eij)/ljd
261         sigder=ed*sigder
262         eom1=eps2der*eps2rt_om1-2.0D0*alf1*eps3der+sigder*sigsq_om1
263         eom2=eps2der*eps2rt_om2+2.0D0*alf2*eps3der+sigder*sigsq_om2
264         eom12=eij*eps1_om12+eps2der*eps2rt_om12
265      &       -2.0D0*alf12*eps3der+sigder*sigsq_om12
266       else if (rij.lt.ssxm) then
267         havebond=.true.
268         ssd=rij-ssXs
269         eij=ssA*ssd*ssd+ssB*ssd+ssC
270
271         ed=2*akcm*ssd+akct*deltat12
272         pom1=akct*ssd
273         pom2=v1ss+2*v2ss*cosphi+3*v3ss*cosphi*cosphi
274         eom1=-2*akth*deltat1-pom1-om2*pom2
275         eom2= 2*akth*deltat2+pom1-om1*pom2
276         eom12=pom2
277       else
278         omega=v1ss+2.0d0*v2ss*cosphi+3.0d0*v3ss*cosphi*cosphi
279
280         d_ssxm(1)=0.5D0*akct/ssA
281         d_ssxm(2)=-d_ssxm(1)
282         d_ssxm(3)=0.0D0
283
284         d_ljxm(1)=sig0ij/sqrt(sigsq**3)
285         d_ljxm(2)=d_ljxm(1)*sigsq_om2
286         d_ljxm(3)=d_ljxm(1)*sigsq_om12
287         d_ljxm(1)=d_ljxm(1)*sigsq_om1
288
289 c-------FIRST METHOD, DISCONTINUOUS SECOND DERIVATIVE
290         xm=0.5d0*(ssxm+ljxm)
291         do k=1,3
292           d_xm(k)=0.5d0*(d_ssxm(k)+d_ljxm(k))
293         enddo
294         if (rij.lt.xm) then
295           havebond=.true.
296           ssm=ssC-0.25D0*ssB*ssB/ssA
297           d_ssm(1)=0.5D0*akct*ssB/ssA
298           d_ssm(2)=2.0D0*akth*deltat2-om1*omega-d_ssm(1)
299           d_ssm(1)=-2.0D0*akth*deltat1-om2*omega+d_ssm(1)
300           d_ssm(3)=omega
301           f1=(rij-xm)/(ssxm-xm)
302           f2=(rij-ssxm)/(xm-ssxm)
303           h1=h_base(f1,hd1)
304           h2=h_base(f2,hd2)
305           eij=ssm*h1+Ht*h2
306           delta_inv=1.0d0/(xm-ssxm)
307           deltasq_inv=delta_inv*delta_inv
308           fac=ssm*hd1-Ht*hd2
309           fac1=deltasq_inv*fac*(xm-rij)
310           fac2=deltasq_inv*fac*(rij-ssxm)
311           ed=delta_inv*(Ht*hd2-ssm*hd1)
312           eom1=fac1*d_ssxm(1)+fac2*d_xm(1)+h1*d_ssm(1)
313           eom2=fac1*d_ssxm(2)+fac2*d_xm(2)+h1*d_ssm(2)
314           eom12=fac1*d_ssxm(3)+fac2*d_xm(3)+h1*d_ssm(3)
315         else
316           havebond=.false.
317           ljm=-0.25D0*ljB*bb(itypi,itypj)/aa(itypi,itypj)
318           d_ljm(1)=-0.5D0*bb(itypi,itypj)/aa(itypi,itypj)*ljB
319           d_ljm(2)=d_ljm(1)*(0.5D0*eps2rt_om2/eps2rt+alf2/eps3rt)
320           d_ljm(3)=d_ljm(1)*(0.5D0*eps1_om12+0.5D0*eps2rt_om12/eps2rt-
321      +         alf12/eps3rt)
322           d_ljm(1)=d_ljm(1)*(0.5D0*eps2rt_om1/eps2rt-alf1/eps3rt)
323           f1=(rij-ljxm)/(xm-ljxm)
324           f2=(rij-xm)/(ljxm-xm)
325           h1=h_base(f1,hd1)
326           h2=h_base(f2,hd2)
327           eij=Ht*h1+ljm*h2
328           delta_inv=1.0d0/(ljxm-xm)
329           deltasq_inv=delta_inv*delta_inv
330           fac=Ht*hd1-ljm*hd2
331           fac1=deltasq_inv*fac*(ljxm-rij)
332           fac2=deltasq_inv*fac*(rij-xm)
333           ed=delta_inv*(ljm*hd2-Ht*hd1)
334           eom1=fac1*d_xm(1)+fac2*d_ljxm(1)+h2*d_ljm(1)
335           eom2=fac1*d_xm(2)+fac2*d_ljxm(2)+h2*d_ljm(2)
336           eom12=fac1*d_xm(3)+fac2*d_ljxm(3)+h2*d_ljm(3)
337         endif
338 c-------END FIRST METHOD, DISCONTINUOUS SECOND DERIVATIVE
339
340 c-------SECOND METHOD, CONTINUOUS SECOND DERIVATIVE
341 c$$$        ssd=rij-ssXs
342 c$$$        ljd=rij-ljXs
343 c$$$        fac1=rij-ljxm
344 c$$$        fac2=rij-ssxm
345 c$$$
346 c$$$        d_ljB(1)=ljB*(eps2rt_om1/eps2rt-2.0d0*alf1/eps3rt)
347 c$$$        d_ljB(2)=ljB*(eps2rt_om2/eps2rt+2.0d0*alf2/eps3rt)
348 c$$$        d_ljB(3)=ljB*(eps1_om12+eps2rt_om12/eps2rt-2.0d0*alf12/eps3rt)
349 c$$$
350 c$$$        ssm=ssC-0.25D0*ssB*ssB/ssA
351 c$$$        d_ssm(1)=0.5D0*akct*ssB/ssA
352 c$$$        d_ssm(2)=2.0D0*akth*deltat2-om1*omega-d_ssm(1)
353 c$$$        d_ssm(1)=-2.0D0*akth*deltat1-om2*omega+d_ssm(1)
354 c$$$        d_ssm(3)=omega
355 c$$$
356 c$$$        ljm=-0.25D0*bb(itypi,itypj)/aa(itypi,itypj)
357 c$$$        do k=1,3
358 c$$$          d_ljm(k)=ljm*d_ljB(k)
359 c$$$        enddo
360 c$$$        ljm=ljm*ljB
361 c$$$
362 c$$$        ss=ssA*ssd*ssd+ssB*ssd+ssC
363 c$$$        d_ss(0)=2.0d0*ssA*ssd+ssB
364 c$$$        d_ss(2)=akct*ssd
365 c$$$        d_ss(1)=-d_ss(2)-2.0d0*akth*deltat1-om2*omega
366 c$$$        d_ss(2)=d_ss(2)+2.0d0*akth*deltat2-om1*omega
367 c$$$        d_ss(3)=omega
368 c$$$
369 c$$$        ljf=bb(itypi,itypj)/aa(itypi,itypj)
370 c$$$        ljf=9.0d0*ljf*(-0.5d0*ljf)**(1.0d0/3.0d0)
371 c$$$        d_ljf(0)=ljf*2.0d0*ljB*fac1
372 c$$$        do k=1,3
373 c$$$          d_ljf(k)=d_ljm(k)+ljf*(d_ljB(k)*fac1*fac1-
374 c$$$     &         2.0d0*ljB*fac1*d_ljxm(k))
375 c$$$        enddo
376 c$$$        ljf=ljm+ljf*ljB*fac1*fac1
377 c$$$
378 c$$$        f1=(rij-ljxm)/(ssxm-ljxm)
379 c$$$        f2=(rij-ssxm)/(ljxm-ssxm)
380 c$$$        h1=h_base(f1,hd1)
381 c$$$        h2=h_base(f2,hd2)
382 c$$$        eij=ss*h1+ljf*h2
383 c$$$        delta_inv=1.0d0/(ljxm-ssxm)
384 c$$$        deltasq_inv=delta_inv*delta_inv
385 c$$$        fac=ljf*hd2-ss*hd1
386 c$$$        ed=d_ss(0)*h1+d_ljf(0)*h2+delta_inv*fac
387 c$$$        eom1=d_ss(1)*h1+d_ljf(1)*h2+deltasq_inv*fac*
388 c$$$     &       (fac1*d_ssxm(1)-fac2*(d_ljxm(1)))
389 c$$$        eom2=d_ss(2)*h1+d_ljf(2)*h2+deltasq_inv*fac*
390 c$$$     &       (fac1*d_ssxm(2)-fac2*(d_ljxm(2)))
391 c$$$        eom12=d_ss(3)*h1+d_ljf(3)*h2+deltasq_inv*fac*
392 c$$$     &       (fac1*d_ssxm(3)-fac2*(d_ljxm(3)))
393 c$$$
394 c$$$        havebond=.false.
395 c$$$        if (ed.gt.0.0d0) havebond=.true.
396 c-------END SECOND METHOD, CONTINUOUS SECOND DERIVATIVE
397
398       endif
399
400       if (havebond) then
401 #ifndef CLUST
402 #ifndef WHAM
403 c        if (dyn_ssbond_ij(i,j).eq.1.0d300) then
404 c          write(iout,'(a15,f12.2,f8.1,2i5)')
405 c     &         "SSBOND_E_FORM",totT,t_bath,i,j
406 c        endif
407 #endif
408 #endif
409         dyn_ssbond_ij(i,j)=eij
410       else if (.not.havebond .and. dyn_ssbond_ij(i,j).lt.1.0d300) then
411         dyn_ssbond_ij(i,j)=1.0d300
412 #ifndef CLUST
413 #ifndef WHAM
414 c        write(iout,'(a15,f12.2,f8.1,2i5)')
415 c     &       "SSBOND_E_BREAK",totT,t_bath,i,j
416 #endif
417 #endif
418       endif
419
420 c-------TESTING CODE
421       if (checkstop) then
422         if (jcheck.eq.0) write(iout,'(a,3f15.8,$)')
423      &       "CHECKSTOP",rij,eij,ed
424         echeck(jcheck)=eij
425       endif
426       enddo
427       if (checkstop) then
428         write(iout,'(f15.8)')(echeck(1)-echeck(-1))*0.5d0/deps
429       endif
430       enddo
431       if (checkstop) then
432         transgrad=.true.
433         checkstop=.false.
434       endif
435 c-------END TESTING CODE
436
437       do k=1,3
438         dcosom1(k)=(dc_norm(k,nres+i)-om1*erij(k))/rij
439         dcosom2(k)=(dc_norm(k,nres+j)-om2*erij(k))/rij
440       enddo
441       do k=1,3
442         gg(k)=ed*erij(k)+eom1*dcosom1(k)+eom2*dcosom2(k)
443       enddo
444       do k=1,3
445         gvdwx(k,i)=gvdwx(k,i)-gg(k)
446      &       +(eom12*(dc_norm(k,nres+j)-om12*dc_norm(k,nres+i))
447      &       +eom1*(erij(k)-om1*dc_norm(k,nres+i)))*dsci_inv
448         gvdwx(k,j)=gvdwx(k,j)+gg(k)
449      &       +(eom12*(dc_norm(k,nres+i)-om12*dc_norm(k,nres+j))
450      &       +eom2*(erij(k)-om2*dc_norm(k,nres+j)))*dscj_inv
451       enddo
452 cgrad      do k=i,j-1
453 cgrad        do l=1,3
454 cgrad          gvdwc(l,k)=gvdwc(l,k)+gg(l)
455 cgrad        enddo
456 cgrad      enddo
457
458       do l=1,3
459         gvdwc(l,i)=gvdwc(l,i)-gg(l)
460         gvdwc(l,j)=gvdwc(l,j)+gg(l)
461       enddo
462
463       return
464       end
465
466 C-----------------------------------------------------------------------------
467
468       double precision function h_base(x,deriv)
469 c     A smooth function going 0->1 in range [0,1]
470 c     It should NOT be called outside range [0,1], it will not work there.
471       implicit none
472
473 c     Input arguments
474       double precision x
475
476 c     Output arguments
477       double precision deriv
478
479 c     Local variables
480       double precision xsq
481
482
483 c     Two parabolas put together.  First derivative zero at extrema
484 c$$$      if (x.lt.0.5D0) then
485 c$$$        h_base=2.0D0*x*x
486 c$$$        deriv=4.0D0*x
487 c$$$      else
488 c$$$        deriv=1.0D0-x
489 c$$$        h_base=1.0D0-2.0D0*deriv*deriv
490 c$$$        deriv=4.0D0*deriv
491 c$$$      endif
492
493 c     Third degree polynomial.  First derivative zero at extrema
494       h_base=x*x*(3.0d0-2.0d0*x)
495       deriv=6.0d0*x*(1.0d0-x)
496
497 c     Fifth degree polynomial.  First and second derivatives zero at extrema
498 c$$$      xsq=x*x
499 c$$$      h_base=x*xsq*(6.0d0*xsq-15.0d0*x+10.0d0)
500 c$$$      deriv=x-1.0d0
501 c$$$      deriv=deriv*deriv
502 c$$$      deriv=30.0d0*xsq*deriv
503
504       return
505       end
506
507 c----------------------------------------------------------------------------
508
509       subroutine dyn_set_nss
510 c     Adjust nss and other relevant variables based on dyn_ssbond_ij
511 c      implicit none
512
513 c     Includes
514       include 'DIMENSIONS'
515 #ifdef MPI
516       include "mpif.h"
517 #endif
518       include 'COMMON.SBRIDGE'
519       include 'COMMON.CHAIN'
520       include 'COMMON.IOUNITS'
521       include 'COMMON.SETUP'
522 #ifndef CLUST
523 #ifndef WHAM
524       include 'COMMON.MD'
525 #endif
526 #endif
527
528 c     Local variables
529       double precision emin
530       integer i,j,imin
531       integer diff,allflag(maxdim),allnss,
532      &     allihpb(maxdim),alljhpb(maxdim),
533      &     newnss,newihpb(maxdim),newjhpb(maxdim)
534       logical found
535       integer i_newnss(max_fg_procs),displ(max_fg_procs)
536       integer g_newihpb(maxdim),g_newjhpb(maxdim),g_newnss
537
538       allnss=0
539       do i=1,nres-1
540         do j=i+1,nres
541           if (dyn_ssbond_ij(i,j).lt.1.0d300) then
542             allnss=allnss+1
543             allflag(allnss)=0
544             allihpb(allnss)=i
545             alljhpb(allnss)=j
546           endif
547         enddo
548       enddo
549
550 cmc      write(iout,*)"ALLNSS ",allnss,(allihpb(i),alljhpb(i),i=1,allnss)
551
552  1    emin=1.0d300
553       do i=1,allnss
554         if (allflag(i).eq.0 .and.
555      &       dyn_ssbond_ij(allihpb(i),alljhpb(i)).lt.emin) then
556           emin=dyn_ssbond_ij(allihpb(i),alljhpb(i))
557           imin=i
558         endif
559       enddo
560       if (emin.lt.1.0d300) then
561         allflag(imin)=1
562         do i=1,allnss
563           if (allflag(i).eq.0 .and.
564      &         (allihpb(i).eq.allihpb(imin) .or.
565      &         alljhpb(i).eq.allihpb(imin) .or.
566      &         allihpb(i).eq.alljhpb(imin) .or.
567      &         alljhpb(i).eq.alljhpb(imin))) then
568             allflag(i)=-1
569           endif
570         enddo
571         goto 1
572       endif
573
574 cmc      write(iout,*)"ALLNSS ",allnss,(allihpb(i),alljhpb(i),i=1,allnss)
575
576       newnss=0
577       do i=1,allnss
578         if (allflag(i).eq.1) then
579           newnss=newnss+1
580           newihpb(newnss)=allihpb(i)
581           newjhpb(newnss)=alljhpb(i)
582         endif
583       enddo
584
585       if (nfgtasks.gt.1)then
586
587         call MPI_Reduce(newnss,g_newnss,1,
588      &    MPI_INTEGER,MPI_SUM,king,FG_COMM,IERR)
589         call MPI_Gather(newnss,1,MPI_INTEGER,
590      &                  i_newnss,1,MPI_INTEGER,king,FG_COMM,IERR)
591         displ(0)=0
592         do i=1,nfgtasks-1,1
593           displ(i)=i_newnss(i-1)+displ(i-1)
594         enddo
595         call MPI_Gatherv(newihpb,newnss,MPI_INTEGER,
596      &                   g_newihpb,i_newnss,displ,MPI_INTEGER,
597      &                   king,FG_COMM,IERR)     
598         call MPI_Gatherv(newjhpb,newnss,MPI_INTEGER,
599      &                   g_newjhpb,i_newnss,displ,MPI_INTEGER,
600      &                   king,FG_COMM,IERR)     
601         if(fg_rank.eq.0) then
602 c         print *,'g_newnss',g_newnss
603 c         print *,'g_newihpb',(g_newihpb(i),i=1,g_newnss)
604 c         print *,'g_newjhpb',(g_newjhpb(i),i=1,g_newnss)
605          newnss=g_newnss  
606          do i=1,newnss
607           newihpb(i)=g_newihpb(i)
608           newjhpb(i)=g_newjhpb(i)
609          enddo
610         endif
611       endif
612
613       diff=newnss-nss
614
615 cmc      write(iout,*)"NEWNSS ",newnss,(newihpb(i),newjhpb(i),i=1,newnss)
616
617       do i=1,nss
618         found=.false.
619         do j=1,newnss
620           if (idssb(i).eq.newihpb(j) .and.
621      &         jdssb(i).eq.newjhpb(j)) found=.true.
622         enddo
623 #ifndef CLUST
624 #ifndef WHAM
625         if (.not.found) write(iout,'(a15,f12.2,f8.1,2i5)')
626      &       "SSBOND_BREAK",totT,t_bath,idssb(i),jdssb(i)
627 #endif
628 #endif
629       enddo
630
631       do i=1,newnss
632         found=.false.
633         do j=1,nss
634           if (newihpb(i).eq.idssb(j) .and.
635      &         newjhpb(i).eq.jdssb(j)) found=.true.
636         enddo
637 #ifndef CLUST
638 #ifndef WHAM
639         if (.not.found) write(iout,'(a15,f12.2,f8.1,2i5)')
640      &       "SSBOND_FORM",totT,t_bath,newihpb(i),newjhpb(i)
641 #endif
642 #endif
643       enddo
644
645       nss=newnss
646       do i=1,nss
647         idssb(i)=newihpb(i)
648         jdssb(i)=newjhpb(i)
649       enddo
650
651       return
652       end
653
654 c----------------------------------------------------------------------------
655
656 #ifdef WHAM
657       subroutine read_ssHist
658       implicit none
659
660 c     Includes
661       include 'DIMENSIONS'
662       include "DIMENSIONS.FREE"
663       include 'COMMON.FREE'
664
665 c     Local variables
666       integer i,j
667       character*80 controlcard
668
669       do i=1,dyn_nssHist
670         call card_concat(controlcard,.true.)
671         read(controlcard,*)
672      &       dyn_ssHist(i,0),(dyn_ssHist(i,j),j=1,2*dyn_ssHist(i,0))
673       enddo
674
675       return
676       end
677 #endif
678
679 c----------------------------------------------------------------------------
680
681
682 C-----------------------------------------------------------------------------
683 C-----------------------------------------------------------------------------
684 C-----------------------------------------------------------------------------
685 C-----------------------------------------------------------------------------
686 C-----------------------------------------------------------------------------
687 C-----------------------------------------------------------------------------
688 C-----------------------------------------------------------------------------
689
690 c$$$c-----------------------------------------------------------------------------
691 c$$$
692 c$$$      subroutine ss_relax(i_in,j_in)
693 c$$$      implicit none
694 c$$$
695 c$$$c     Includes
696 c$$$      include 'DIMENSIONS'
697 c$$$      include 'COMMON.VAR'
698 c$$$      include 'COMMON.CHAIN'
699 c$$$      include 'COMMON.IOUNITS'
700 c$$$      include 'COMMON.INTERACT'
701 c$$$
702 c$$$c     Input arguments
703 c$$$      integer i_in,j_in
704 c$$$
705 c$$$c     Local variables
706 c$$$      integer i,iretcode,nfun_sc
707 c$$$      logical scfail
708 c$$$      double precision var(maxvar),e_sc,etot
709 c$$$
710 c$$$
711 c$$$      mask_r=.true.
712 c$$$      do i=nnt,nct
713 c$$$        mask_side(i)=0
714 c$$$      enddo
715 c$$$      mask_side(i_in)=1
716 c$$$      mask_side(j_in)=1
717 c$$$
718 c$$$c     Minimize the two selected side-chains
719 c$$$      call overlap_sc(scfail)  ! Better not fail!
720 c$$$      call minimize_sc(e_sc,var,iretcode,nfun_sc)
721 c$$$
722 c$$$      mask_r=.false.
723 c$$$
724 c$$$      return
725 c$$$      end
726 c$$$
727 c$$$c-------------------------------------------------------------
728 c$$$
729 c$$$      subroutine minimize_sc(etot_sc,iretcode,nfun)
730 c$$$c     Minimize side-chains only, starting from geom but without modifying
731 c$$$c     bond lengths.
732 c$$$c     If mask_r is already set, only the selected side-chains are minimized,
733 c$$$c     otherwise all side-chains are minimized keeping the backbone frozen.
734 c$$$      implicit none
735 c$$$
736 c$$$c     Includes
737 c$$$      include 'DIMENSIONS'
738 c$$$      include 'COMMON.IOUNITS'
739 c$$$      include 'COMMON.VAR'
740 c$$$      include 'COMMON.CHAIN'
741 c$$$      include 'COMMON.GEO'
742 c$$$      include 'COMMON.MINIM'
743 c$$$      integer icall
744 c$$$      common /srutu/ icall
745 c$$$
746 c$$$c     Output arguments
747 c$$$      double precision etot_sc
748 c$$$      integer iretcode,nfun
749 c$$$
750 c$$$c     External functions/subroutines
751 c$$$      external func_sc,grad_sc,fdum
752 c$$$
753 c$$$c     Local variables
754 c$$$      integer liv,lv
755 c$$$      parameter (liv=60,lv=(77+maxvar*(maxvar+17)/2)) 
756 c$$$      integer iv(liv)
757 c$$$      double precision rdum(1)
758 c$$$      double precision d(maxvar),v(1:lv),x(maxvar),xx(maxvar)
759 c$$$      integer idum(1)
760 c$$$      integer i,nvar_restr
761 c$$$
762 c$$$
763 c$$$cmc      start_minim=.true.
764 c$$$      call deflt(2,iv,liv,lv,v)                                         
765 c$$$* 12 means fresh start, dont call deflt                                 
766 c$$$      iv(1)=12                                                          
767 c$$$* max num of fun calls                                                  
768 c$$$      if (maxfun.eq.0) maxfun=500
769 c$$$      iv(17)=maxfun
770 c$$$* max num of iterations                                                 
771 c$$$      if (maxmin.eq.0) maxmin=1000
772 c$$$      iv(18)=maxmin
773 c$$$* controls output                                                       
774 c$$$      iv(19)=1
775 c$$$* selects output unit                                                   
776 c$$$      iv(21)=0
777 c$$$c      iv(21)=iout               ! DEBUG
778 c$$$c      iv(21)=8                  ! DEBUG
779 c$$$* 1 means to print out result                                           
780 c$$$      iv(22)=0
781 c$$$c      iv(22)=1                  ! DEBUG
782 c$$$* 1 means to print out summary stats                                    
783 c$$$      iv(23)=0                                                          
784 c$$$c      iv(23)=1                  ! DEBUG
785 c$$$* 1 means to print initial x and d                                      
786 c$$$      iv(24)=0                                                          
787 c$$$c      iv(24)=1                  ! DEBUG
788 c$$$* min val for v(radfac) default is 0.1                                  
789 c$$$      v(24)=0.1D0                                                       
790 c$$$* max val for v(radfac) default is 4.0                                  
791 c$$$      v(25)=2.0D0                                                       
792 c$$$c     v(25)=4.0D0                                                       
793 c$$$* check false conv if (act fnctn decrease) .lt. v(26)*(exp decrease)    
794 c$$$* the sumsl default is 0.1                                              
795 c$$$      v(26)=0.1D0
796 c$$$* false conv if (act fnctn decrease) .lt. v(34)                         
797 c$$$* the sumsl default is 100*machep                                       
798 c$$$      v(34)=v(34)/100.0D0                                               
799 c$$$* absolute convergence                                                  
800 c$$$      if (tolf.eq.0.0D0) tolf=1.0D-4
801 c$$$      v(31)=tolf
802 c$$$* relative convergence                                                  
803 c$$$      if (rtolf.eq.0.0D0) rtolf=1.0D-1
804 c$$$      v(32)=rtolf
805 c$$$* controls initial step size                                            
806 c$$$       v(35)=1.0D-1                                                    
807 c$$$* large vals of d correspond to small components of step                
808 c$$$      do i=1,nphi
809 c$$$        d(i)=1.0D-1
810 c$$$      enddo
811 c$$$      do i=nphi+1,nvar
812 c$$$        d(i)=1.0D-1
813 c$$$      enddo
814 c$$$
815 c$$$      call geom_to_var(nvar,x)
816 c$$$      IF (mask_r) THEN
817 c$$$        do i=1,nres             ! Just in case...
818 c$$$          mask_phi(i)=0
819 c$$$          mask_theta(i)=0
820 c$$$        enddo
821 c$$$        call x2xx(x,xx,nvar_restr)
822 c$$$        call sumsl(nvar_restr,d,xx,func_sc,grad_sc,
823 c$$$     &       iv,liv,lv,v,idum,rdum,fdum)      
824 c$$$        call xx2x(x,xx)
825 c$$$      ELSE
826 c$$$c     When minimizing ALL side-chains, etotal_sc is a little
827 c$$$c     faster if we don't set mask_r
828 c$$$        do i=1,nres
829 c$$$          mask_phi(i)=0
830 c$$$          mask_theta(i)=0
831 c$$$          mask_side(i)=1
832 c$$$        enddo
833 c$$$        call x2xx(x,xx,nvar_restr)
834 c$$$        call sumsl(nvar_restr,d,xx,func_sc,grad_sc,
835 c$$$     &       iv,liv,lv,v,idum,rdum,fdum)      
836 c$$$        call xx2x(x,xx)
837 c$$$      ENDIF
838 c$$$      call var_to_geom(nvar,x)
839 c$$$      call chainbuild_sc
840 c$$$      etot_sc=v(10)                                                      
841 c$$$      iretcode=iv(1)
842 c$$$      nfun=iv(6)
843 c$$$      return  
844 c$$$      end  
845 c$$$
846 c$$$C--------------------------------------------------------------------------
847 c$$$
848 c$$$      subroutine chainbuild_sc
849 c$$$      implicit none
850 c$$$      include 'DIMENSIONS'
851 c$$$      include 'COMMON.VAR'
852 c$$$      include 'COMMON.INTERACT'
853 c$$$
854 c$$$c     Local variables
855 c$$$      integer i
856 c$$$
857 c$$$
858 c$$$      do i=nnt,nct
859 c$$$        if (.not.mask_r .or. mask_side(i).eq.1) then
860 c$$$          call locate_side_chain(i)
861 c$$$        endif
862 c$$$      enddo
863 c$$$
864 c$$$      return
865 c$$$      end
866 c$$$
867 c$$$C--------------------------------------------------------------------------
868 c$$$
869 c$$$      subroutine func_sc(n,x,nf,f,uiparm,urparm,ufparm)  
870 c$$$      implicit none
871 c$$$
872 c$$$c     Includes
873 c$$$      include 'DIMENSIONS'
874 c$$$      include 'COMMON.DERIV'
875 c$$$      include 'COMMON.VAR'
876 c$$$      include 'COMMON.MINIM'
877 c$$$      include 'COMMON.IOUNITS'
878 c$$$
879 c$$$c     Input arguments
880 c$$$      integer n
881 c$$$      double precision x(maxvar)
882 c$$$      double precision ufparm
883 c$$$      external ufparm
884 c$$$
885 c$$$c     Input/Output arguments
886 c$$$      integer nf
887 c$$$      integer uiparm(1)
888 c$$$      double precision urparm(1)
889 c$$$
890 c$$$c     Output arguments
891 c$$$      double precision f
892 c$$$
893 c$$$c     Local variables
894 c$$$      double precision energia(0:n_ene)
895 c$$$#ifdef OSF
896 c$$$c     Variables used to intercept NaNs
897 c$$$      double precision x_sum
898 c$$$      integer i_NAN
899 c$$$#endif
900 c$$$
901 c$$$
902 c$$$      nfl=nf
903 c$$$      icg=mod(nf,2)+1
904 c$$$
905 c$$$#ifdef OSF
906 c$$$c     Intercept NaNs in the coordinates, before calling etotal_sc
907 c$$$      x_sum=0.D0
908 c$$$      do i_NAN=1,n
909 c$$$        x_sum=x_sum+x(i_NAN)
910 c$$$      enddo
911 c$$$c     Calculate the energy only if the coordinates are ok
912 c$$$      if ((.not.(x_sum.lt.0.D0)) .and. (.not.(x_sum.ge.0.D0))) then
913 c$$$        write(iout,*)"   *** func_restr_sc : Found NaN in coordinates"
914 c$$$        f=1.0D+77
915 c$$$        nf=0
916 c$$$      else
917 c$$$#endif
918 c$$$
919 c$$$      call var_to_geom_restr(n,x)
920 c$$$      call zerograd
921 c$$$      call chainbuild_sc
922 c$$$      call etotal_sc(energia(0))
923 c$$$      f=energia(0)
924 c$$$      if (energia(1).eq.1.0D20 .or. energia(0).eq.1.0D99) nf=0
925 c$$$
926 c$$$#ifdef OSF
927 c$$$      endif
928 c$$$#endif
929 c$$$
930 c$$$      return                                                            
931 c$$$      end                                                               
932 c$$$
933 c$$$c-------------------------------------------------------
934 c$$$
935 c$$$      subroutine grad_sc(n,x,nf,g,uiparm,urparm,ufparm)
936 c$$$      implicit none
937 c$$$
938 c$$$c     Includes
939 c$$$      include 'DIMENSIONS'
940 c$$$      include 'COMMON.CHAIN'
941 c$$$      include 'COMMON.DERIV'
942 c$$$      include 'COMMON.VAR'
943 c$$$      include 'COMMON.INTERACT'
944 c$$$      include 'COMMON.MINIM'
945 c$$$
946 c$$$c     Input arguments
947 c$$$      integer n
948 c$$$      double precision x(maxvar)
949 c$$$      double precision ufparm
950 c$$$      external ufparm
951 c$$$
952 c$$$c     Input/Output arguments
953 c$$$      integer nf
954 c$$$      integer uiparm(1)
955 c$$$      double precision urparm(1)
956 c$$$
957 c$$$c     Output arguments
958 c$$$      double precision g(maxvar)
959 c$$$
960 c$$$c     Local variables
961 c$$$      double precision f,gphii,gthetai,galphai,gomegai
962 c$$$      integer ig,ind,i,j,k,igall,ij
963 c$$$
964 c$$$
965 c$$$      icg=mod(nf,2)+1
966 c$$$      if (nf-nfl+1) 20,30,40
967 c$$$   20 call func_sc(n,x,nf,f,uiparm,urparm,ufparm)
968 c$$$c     write (iout,*) 'grad 20'
969 c$$$      if (nf.eq.0) return
970 c$$$      goto 40
971 c$$$   30 call var_to_geom_restr(n,x)
972 c$$$      call chainbuild_sc
973 c$$$C
974 c$$$C Evaluate the derivatives of virtual bond lengths and SC vectors in variables.
975 c$$$C
976 c$$$   40 call cartder
977 c$$$C
978 c$$$C Convert the Cartesian gradient into internal-coordinate gradient.
979 c$$$C
980 c$$$
981 c$$$      ig=0
982 c$$$      ind=nres-2
983 c$$$      do i=2,nres-2
984 c$$$       IF (mask_phi(i+2).eq.1) THEN
985 c$$$        gphii=0.0D0
986 c$$$        do j=i+1,nres-1
987 c$$$          ind=ind+1
988 c$$$          do k=1,3
989 c$$$            gphii=gphii+dcdv(k+3,ind)*gradc(k,j,icg)
990 c$$$            gphii=gphii+dxdv(k+3,ind)*gradx(k,j,icg)
991 c$$$          enddo
992 c$$$        enddo
993 c$$$        ig=ig+1
994 c$$$        g(ig)=gphii
995 c$$$       ELSE
996 c$$$        ind=ind+nres-1-i
997 c$$$       ENDIF
998 c$$$      enddo                                        
999 c$$$
1000 c$$$
1001 c$$$      ind=0
1002 c$$$      do i=1,nres-2
1003 c$$$       IF (mask_theta(i+2).eq.1) THEN
1004 c$$$        ig=ig+1
1005 c$$$    gthetai=0.0D0
1006 c$$$    do j=i+1,nres-1
1007 c$$$          ind=ind+1
1008 c$$$      do k=1,3
1009 c$$$            gthetai=gthetai+dcdv(k,ind)*gradc(k,j,icg)
1010 c$$$            gthetai=gthetai+dxdv(k,ind)*gradx(k,j,icg)
1011 c$$$          enddo
1012 c$$$        enddo
1013 c$$$        g(ig)=gthetai
1014 c$$$       ELSE
1015 c$$$        ind=ind+nres-1-i
1016 c$$$       ENDIF
1017 c$$$      enddo
1018 c$$$
1019 c$$$      do i=2,nres-1
1020 c$$$    if (itype(i).ne.10) then
1021 c$$$         IF (mask_side(i).eq.1) THEN
1022 c$$$          ig=ig+1
1023 c$$$          galphai=0.0D0
1024 c$$$      do k=1,3
1025 c$$$        galphai=galphai+dxds(k,i)*gradx(k,i,icg)
1026 c$$$          enddo
1027 c$$$          g(ig)=galphai
1028 c$$$         ENDIF
1029 c$$$        endif
1030 c$$$      enddo
1031 c$$$
1032 c$$$      
1033 c$$$      do i=2,nres-1
1034 c$$$        if (itype(i).ne.10) then
1035 c$$$         IF (mask_side(i).eq.1) THEN
1036 c$$$          ig=ig+1
1037 c$$$      gomegai=0.0D0
1038 c$$$      do k=1,3
1039 c$$$        gomegai=gomegai+dxds(k+3,i)*gradx(k,i,icg)
1040 c$$$          enddo
1041 c$$$      g(ig)=gomegai
1042 c$$$         ENDIF
1043 c$$$        endif
1044 c$$$      enddo
1045 c$$$
1046 c$$$C
1047 c$$$C Add the components corresponding to local energy terms.
1048 c$$$C
1049 c$$$
1050 c$$$      ig=0
1051 c$$$      igall=0
1052 c$$$      do i=4,nres
1053 c$$$        igall=igall+1
1054 c$$$        if (mask_phi(i).eq.1) then
1055 c$$$          ig=ig+1
1056 c$$$          g(ig)=g(ig)+gloc(igall,icg)
1057 c$$$        endif
1058 c$$$      enddo
1059 c$$$
1060 c$$$      do i=3,nres
1061 c$$$        igall=igall+1
1062 c$$$        if (mask_theta(i).eq.1) then
1063 c$$$          ig=ig+1
1064 c$$$          g(ig)=g(ig)+gloc(igall,icg)
1065 c$$$        endif
1066 c$$$      enddo
1067 c$$$     
1068 c$$$      do ij=1,2
1069 c$$$      do i=2,nres-1
1070 c$$$        if (itype(i).ne.10) then
1071 c$$$          igall=igall+1
1072 c$$$          if (mask_side(i).eq.1) then
1073 c$$$            ig=ig+1
1074 c$$$            g(ig)=g(ig)+gloc(igall,icg)
1075 c$$$          endif
1076 c$$$        endif
1077 c$$$      enddo
1078 c$$$      enddo
1079 c$$$
1080 c$$$cd      do i=1,ig
1081 c$$$cd        write (iout,'(a2,i5,a3,f25.8)') 'i=',i,' g=',g(i)
1082 c$$$cd      enddo
1083 c$$$
1084 c$$$      return
1085 c$$$      end
1086 c$$$
1087 c$$$C-----------------------------------------------------------------------------
1088 c$$$
1089 c$$$      subroutine etotal_sc(energy_sc)
1090 c$$$      implicit none
1091 c$$$
1092 c$$$c     Includes
1093 c$$$      include 'DIMENSIONS'
1094 c$$$      include 'COMMON.VAR'
1095 c$$$      include 'COMMON.INTERACT'
1096 c$$$      include 'COMMON.DERIV'
1097 c$$$      include 'COMMON.FFIELD'
1098 c$$$
1099 c$$$c     Output arguments
1100 c$$$      double precision energy_sc(0:n_ene)
1101 c$$$
1102 c$$$c     Local variables
1103 c$$$      double precision evdw,escloc
1104 c$$$      integer i,j
1105 c$$$
1106 c$$$
1107 c$$$      do i=1,n_ene
1108 c$$$        energy_sc(i)=0.0D0
1109 c$$$      enddo
1110 c$$$
1111 c$$$      if (mask_r) then
1112 c$$$        call egb_sc(evdw)
1113 c$$$        call esc_sc(escloc)
1114 c$$$      else
1115 c$$$        call egb(evdw)
1116 c$$$        call esc(escloc)
1117 c$$$      endif
1118 c$$$
1119 c$$$      if (evdw.eq.1.0D20) then
1120 c$$$        energy_sc(0)=evdw
1121 c$$$      else
1122 c$$$        energy_sc(0)=wsc*evdw+wscloc*escloc
1123 c$$$      endif
1124 c$$$      energy_sc(1)=evdw
1125 c$$$      energy_sc(12)=escloc
1126 c$$$
1127 c$$$C
1128 c$$$C Sum up the components of the Cartesian gradient.
1129 c$$$C
1130 c$$$      do i=1,nct
1131 c$$$        do j=1,3
1132 c$$$          gradx(j,i,icg)=wsc*gvdwx(j,i)
1133 c$$$        enddo
1134 c$$$      enddo
1135 c$$$
1136 c$$$      return
1137 c$$$      end
1138 c$$$
1139 c$$$C-----------------------------------------------------------------------------
1140 c$$$
1141 c$$$      subroutine egb_sc(evdw)
1142 c$$$C
1143 c$$$C This subroutine calculates the interaction energy of nonbonded side chains
1144 c$$$C assuming the Gay-Berne potential of interaction.
1145 c$$$C
1146 c$$$      implicit real*8 (a-h,o-z)
1147 c$$$      include 'DIMENSIONS'
1148 c$$$      include 'COMMON.GEO'
1149 c$$$      include 'COMMON.VAR'
1150 c$$$      include 'COMMON.LOCAL'
1151 c$$$      include 'COMMON.CHAIN'
1152 c$$$      include 'COMMON.DERIV'
1153 c$$$      include 'COMMON.NAMES'
1154 c$$$      include 'COMMON.INTERACT'
1155 c$$$      include 'COMMON.IOUNITS'
1156 c$$$      include 'COMMON.CALC'
1157 c$$$      include 'COMMON.CONTROL'
1158 c$$$      logical lprn
1159 c$$$      evdw=0.0D0
1160 c$$$      energy_dec=.false.
1161 c$$$c     print *,'Entering EGB nnt=',nnt,' nct=',nct,' expon=',expon
1162 c$$$      evdw=0.0D0
1163 c$$$      lprn=.false.
1164 c$$$c     if (icall.eq.0) lprn=.false.
1165 c$$$      ind=0
1166 c$$$      do i=iatsc_s,iatsc_e
1167 c$$$        itypi=itype(i)
1168 c$$$        itypi1=itype(i+1)
1169 c$$$        xi=c(1,nres+i)
1170 c$$$        yi=c(2,nres+i)
1171 c$$$        zi=c(3,nres+i)
1172 c$$$        dxi=dc_norm(1,nres+i)
1173 c$$$        dyi=dc_norm(2,nres+i)
1174 c$$$        dzi=dc_norm(3,nres+i)
1175 c$$$c        dsci_inv=dsc_inv(itypi)
1176 c$$$        dsci_inv=vbld_inv(i+nres)
1177 c$$$c        write (iout,*) "i",i,dsc_inv(itypi),dsci_inv,1.0d0/vbld(i+nres)
1178 c$$$c        write (iout,*) "dcnori",dxi*dxi+dyi*dyi+dzi*dzi
1179 c$$$C
1180 c$$$C Calculate SC interaction energy.
1181 c$$$C
1182 c$$$        do iint=1,nint_gr(i)
1183 c$$$          do j=istart(i,iint),iend(i,iint)
1184 c$$$          IF (mask_side(j).eq.1.or.mask_side(i).eq.1) THEN
1185 c$$$            ind=ind+1
1186 c$$$            itypj=itype(j)
1187 c$$$c            dscj_inv=dsc_inv(itypj)
1188 c$$$            dscj_inv=vbld_inv(j+nres)
1189 c$$$c            write (iout,*) "j",j,dsc_inv(itypj),dscj_inv,
1190 c$$$c     &       1.0d0/vbld(j+nres)
1191 c$$$c            write (iout,*) "i",i," j", j," itype",itype(i),itype(j)
1192 c$$$            sig0ij=sigma(itypi,itypj)
1193 c$$$            chi1=chi(itypi,itypj)
1194 c$$$            chi2=chi(itypj,itypi)
1195 c$$$            chi12=chi1*chi2
1196 c$$$            chip1=chip(itypi)
1197 c$$$            chip2=chip(itypj)
1198 c$$$            chip12=chip1*chip2
1199 c$$$            alf1=alp(itypi)
1200 c$$$            alf2=alp(itypj)
1201 c$$$            alf12=0.5D0*(alf1+alf2)
1202 c$$$C For diagnostics only!!!
1203 c$$$c           chi1=0.0D0
1204 c$$$c           chi2=0.0D0
1205 c$$$c           chi12=0.0D0
1206 c$$$c           chip1=0.0D0
1207 c$$$c           chip2=0.0D0
1208 c$$$c           chip12=0.0D0
1209 c$$$c           alf1=0.0D0
1210 c$$$c           alf2=0.0D0
1211 c$$$c           alf12=0.0D0
1212 c$$$            xj=c(1,nres+j)-xi
1213 c$$$            yj=c(2,nres+j)-yi
1214 c$$$            zj=c(3,nres+j)-zi
1215 c$$$            dxj=dc_norm(1,nres+j)
1216 c$$$            dyj=dc_norm(2,nres+j)
1217 c$$$            dzj=dc_norm(3,nres+j)
1218 c$$$c            write (iout,*) "dcnorj",dxi*dxi+dyi*dyi+dzi*dzi
1219 c$$$c            write (iout,*) "j",j," dc_norm",
1220 c$$$c     &       dc_norm(1,nres+j),dc_norm(2,nres+j),dc_norm(3,nres+j)
1221 c$$$            rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
1222 c$$$            rij=dsqrt(rrij)
1223 c$$$C Calculate angle-dependent terms of energy and contributions to their
1224 c$$$C derivatives.
1225 c$$$            call sc_angular
1226 c$$$            sigsq=1.0D0/sigsq
1227 c$$$            sig=sig0ij*dsqrt(sigsq)
1228 c$$$            rij_shift=1.0D0/rij-sig+sig0ij
1229 c$$$c for diagnostics; uncomment
1230 c$$$c            rij_shift=1.2*sig0ij
1231 c$$$C I hate to put IF's in the loops, but here don't have another choice!!!!
1232 c$$$            if (rij_shift.le.0.0D0) then
1233 c$$$              evdw=1.0D20
1234 c$$$cd              write (iout,'(2(a3,i3,2x),17(0pf7.3))')
1235 c$$$cd     &        restyp(itypi),i,restyp(itypj),j,
1236 c$$$cd     &        rij_shift,1.0D0/rij,sig,sig0ij,sigsq,1-dsqrt(sigsq) 
1237 c$$$              return
1238 c$$$            endif
1239 c$$$            sigder=-sig*sigsq
1240 c$$$c---------------------------------------------------------------
1241 c$$$            rij_shift=1.0D0/rij_shift 
1242 c$$$            fac=rij_shift**expon
1243 c$$$            e1=fac*fac*aa(itypi,itypj)
1244 c$$$            e2=fac*bb(itypi,itypj)
1245 c$$$            evdwij=eps1*eps2rt*eps3rt*(e1+e2)
1246 c$$$            eps2der=evdwij*eps3rt
1247 c$$$            eps3der=evdwij*eps2rt
1248 c$$$c            write (iout,*) "sigsq",sigsq," sig",sig," eps2rt",eps2rt,
1249 c$$$c     &        " eps3rt",eps3rt," eps1",eps1," e1",e1," e2",e2
1250 c$$$            evdwij=evdwij*eps2rt*eps3rt
1251 c$$$            evdw=evdw+evdwij
1252 c$$$            if (lprn) then
1253 c$$$            sigm=dabs(aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
1254 c$$$            epsi=bb(itypi,itypj)**2/aa(itypi,itypj)
1255 c$$$            write (iout,'(2(a3,i3,2x),17(0pf7.3))')
1256 c$$$     &        restyp(itypi),i,restyp(itypj),j,
1257 c$$$     &        epsi,sigm,chi1,chi2,chip1,chip2,
1258 c$$$     &        eps1,eps2rt**2,eps3rt**2,sig,sig0ij,
1259 c$$$     &        om1,om2,om12,1.0D0/rij,1.0D0/rij_shift,
1260 c$$$     &        evdwij
1261 c$$$            endif
1262 c$$$
1263 c$$$            if (energy_dec) write (iout,'(a6,2i,0pf7.3)') 
1264 c$$$     &                        'evdw',i,j,evdwij
1265 c$$$
1266 c$$$C Calculate gradient components.
1267 c$$$            e1=e1*eps1*eps2rt**2*eps3rt**2
1268 c$$$            fac=-expon*(e1+evdwij)*rij_shift
1269 c$$$            sigder=fac*sigder
1270 c$$$            fac=rij*fac
1271 c$$$c            fac=0.0d0
1272 c$$$C Calculate the radial part of the gradient
1273 c$$$            gg(1)=xj*fac
1274 c$$$            gg(2)=yj*fac
1275 c$$$            gg(3)=zj*fac
1276 c$$$C Calculate angular part of the gradient.
1277 c$$$            call sc_grad
1278 c$$$          ENDIF
1279 c$$$          enddo      ! j
1280 c$$$        enddo        ! iint
1281 c$$$      enddo          ! i
1282 c$$$      energy_dec=.false.
1283 c$$$      return
1284 c$$$      end
1285 c$$$
1286 c$$$c-----------------------------------------------------------------------------
1287 c$$$
1288 c$$$      subroutine esc_sc(escloc)
1289 c$$$C Calculate the local energy of a side chain and its derivatives in the
1290 c$$$C corresponding virtual-bond valence angles THETA and the spherical angles 
1291 c$$$C ALPHA and OMEGA.
1292 c$$$      implicit real*8 (a-h,o-z)
1293 c$$$      include 'DIMENSIONS'
1294 c$$$      include 'COMMON.GEO'
1295 c$$$      include 'COMMON.LOCAL'
1296 c$$$      include 'COMMON.VAR'
1297 c$$$      include 'COMMON.INTERACT'
1298 c$$$      include 'COMMON.DERIV'
1299 c$$$      include 'COMMON.CHAIN'
1300 c$$$      include 'COMMON.IOUNITS'
1301 c$$$      include 'COMMON.NAMES'
1302 c$$$      include 'COMMON.FFIELD'
1303 c$$$      include 'COMMON.CONTROL'
1304 c$$$      double precision x(3),dersc(3),xemp(3),dersc0(3),dersc1(3),
1305 c$$$     &     ddersc0(3),ddummy(3),xtemp(3),temp(3)
1306 c$$$      common /sccalc/ time11,time12,time112,theti,it,nlobit
1307 c$$$      delta=0.02d0*pi
1308 c$$$      escloc=0.0D0
1309 c$$$c     write (iout,'(a)') 'ESC'
1310 c$$$      do i=loc_start,loc_end
1311 c$$$      IF (mask_side(i).eq.1) THEN
1312 c$$$        it=itype(i)
1313 c$$$        if (it.eq.10) goto 1
1314 c$$$        nlobit=nlob(it)
1315 c$$$c       print *,'i=',i,' it=',it,' nlobit=',nlobit
1316 c$$$c       write (iout,*) 'i=',i,' ssa=',ssa,' ssad=',ssad
1317 c$$$        theti=theta(i+1)-pipol
1318 c$$$        x(1)=dtan(theti)
1319 c$$$        x(2)=alph(i)
1320 c$$$        x(3)=omeg(i)
1321 c$$$
1322 c$$$        if (x(2).gt.pi-delta) then
1323 c$$$          xtemp(1)=x(1)
1324 c$$$          xtemp(2)=pi-delta
1325 c$$$          xtemp(3)=x(3)
1326 c$$$          call enesc(xtemp,escloci0,dersc0,ddersc0,.true.)
1327 c$$$          xtemp(2)=pi
1328 c$$$          call enesc(xtemp,escloci1,dersc1,ddummy,.false.)
1329 c$$$          call spline1(x(2),pi-delta,delta,escloci0,escloci1,dersc0(2),
1330 c$$$     &        escloci,dersc(2))
1331 c$$$          call spline2(x(2),pi-delta,delta,dersc0(1),dersc1(1),
1332 c$$$     &        ddersc0(1),dersc(1))
1333 c$$$          call spline2(x(2),pi-delta,delta,dersc0(3),dersc1(3),
1334 c$$$     &        ddersc0(3),dersc(3))
1335 c$$$          xtemp(2)=pi-delta
1336 c$$$          call enesc_bound(xtemp,esclocbi0,dersc0,dersc12,.true.)
1337 c$$$          xtemp(2)=pi
1338 c$$$          call enesc_bound(xtemp,esclocbi1,dersc1,chuju,.false.)
1339 c$$$          call spline1(x(2),pi-delta,delta,esclocbi0,esclocbi1,
1340 c$$$     &            dersc0(2),esclocbi,dersc02)
1341 c$$$          call spline2(x(2),pi-delta,delta,dersc0(1),dersc1(1),
1342 c$$$     &            dersc12,dersc01)
1343 c$$$          call splinthet(x(2),0.5d0*delta,ss,ssd)
1344 c$$$          dersc0(1)=dersc01
1345 c$$$          dersc0(2)=dersc02
1346 c$$$          dersc0(3)=0.0d0
1347 c$$$          do k=1,3
1348 c$$$            dersc(k)=ss*dersc(k)+(1.0d0-ss)*dersc0(k)
1349 c$$$          enddo
1350 c$$$          dersc(2)=dersc(2)+ssd*(escloci-esclocbi)
1351 c$$$c         write (iout,*) 'i=',i,x(2)*rad2deg,escloci0,escloci,
1352 c$$$c    &             esclocbi,ss,ssd
1353 c$$$          escloci=ss*escloci+(1.0d0-ss)*esclocbi
1354 c$$$c         escloci=esclocbi
1355 c$$$c         write (iout,*) escloci
1356 c$$$        else if (x(2).lt.delta) then
1357 c$$$          xtemp(1)=x(1)
1358 c$$$          xtemp(2)=delta
1359 c$$$          xtemp(3)=x(3)
1360 c$$$          call enesc(xtemp,escloci0,dersc0,ddersc0,.true.)
1361 c$$$          xtemp(2)=0.0d0
1362 c$$$          call enesc(xtemp,escloci1,dersc1,ddummy,.false.)
1363 c$$$          call spline1(x(2),delta,-delta,escloci0,escloci1,dersc0(2),
1364 c$$$     &        escloci,dersc(2))
1365 c$$$          call spline2(x(2),delta,-delta,dersc0(1),dersc1(1),
1366 c$$$     &        ddersc0(1),dersc(1))
1367 c$$$          call spline2(x(2),delta,-delta,dersc0(3),dersc1(3),
1368 c$$$     &        ddersc0(3),dersc(3))
1369 c$$$          xtemp(2)=delta
1370 c$$$          call enesc_bound(xtemp,esclocbi0,dersc0,dersc12,.true.)
1371 c$$$          xtemp(2)=0.0d0
1372 c$$$          call enesc_bound(xtemp,esclocbi1,dersc1,chuju,.false.)
1373 c$$$          call spline1(x(2),delta,-delta,esclocbi0,esclocbi1,
1374 c$$$     &            dersc0(2),esclocbi,dersc02)
1375 c$$$          call spline2(x(2),delta,-delta,dersc0(1),dersc1(1),
1376 c$$$     &            dersc12,dersc01)
1377 c$$$          dersc0(1)=dersc01
1378 c$$$          dersc0(2)=dersc02
1379 c$$$          dersc0(3)=0.0d0
1380 c$$$          call splinthet(x(2),0.5d0*delta,ss,ssd)
1381 c$$$          do k=1,3
1382 c$$$            dersc(k)=ss*dersc(k)+(1.0d0-ss)*dersc0(k)
1383 c$$$          enddo
1384 c$$$          dersc(2)=dersc(2)+ssd*(escloci-esclocbi)
1385 c$$$c         write (iout,*) 'i=',i,x(2)*rad2deg,escloci0,escloci,
1386 c$$$c    &             esclocbi,ss,ssd
1387 c$$$          escloci=ss*escloci+(1.0d0-ss)*esclocbi
1388 c$$$c         write (iout,*) escloci
1389 c$$$        else
1390 c$$$          call enesc(x,escloci,dersc,ddummy,.false.)
1391 c$$$        endif
1392 c$$$
1393 c$$$        escloc=escloc+escloci
1394 c$$$        if (energy_dec) write (iout,'(a6,i,0pf7.3)')
1395 c$$$     &     'escloc',i,escloci
1396 c$$$c       write (iout,*) 'i=',i,' escloci=',escloci,' dersc=',dersc
1397 c$$$
1398 c$$$        gloc(nphi+i-1,icg)=gloc(nphi+i-1,icg)+
1399 c$$$     &   wscloc*dersc(1)
1400 c$$$        gloc(ialph(i,1),icg)=wscloc*dersc(2)
1401 c$$$        gloc(ialph(i,1)+nside,icg)=wscloc*dersc(3)
1402 c$$$    1   continue
1403 c$$$      ENDIF
1404 c$$$      enddo
1405 c$$$      return
1406 c$$$      end
1407 c$$$
1408 c$$$C-----------------------------------------------------------------------------
1409 c$$$
1410 c$$$      subroutine egb_ij(i_sc,j_sc,evdw)
1411 c$$$C
1412 c$$$C This subroutine calculates the interaction energy of nonbonded side chains
1413 c$$$C assuming the Gay-Berne potential of interaction.
1414 c$$$C
1415 c$$$      implicit real*8 (a-h,o-z)
1416 c$$$      include 'DIMENSIONS'
1417 c$$$      include 'COMMON.GEO'
1418 c$$$      include 'COMMON.VAR'
1419 c$$$      include 'COMMON.LOCAL'
1420 c$$$      include 'COMMON.CHAIN'
1421 c$$$      include 'COMMON.DERIV'
1422 c$$$      include 'COMMON.NAMES'
1423 c$$$      include 'COMMON.INTERACT'
1424 c$$$      include 'COMMON.IOUNITS'
1425 c$$$      include 'COMMON.CALC'
1426 c$$$      include 'COMMON.CONTROL'
1427 c$$$      logical lprn
1428 c$$$      evdw=0.0D0
1429 c$$$      energy_dec=.false.
1430 c$$$c     print *,'Entering EGB nnt=',nnt,' nct=',nct,' expon=',expon
1431 c$$$      evdw=0.0D0
1432 c$$$      lprn=.false.
1433 c$$$      ind=0
1434 c$$$c$$$      do i=iatsc_s,iatsc_e
1435 c$$$      i=i_sc
1436 c$$$        itypi=itype(i)
1437 c$$$        itypi1=itype(i+1)
1438 c$$$        xi=c(1,nres+i)
1439 c$$$        yi=c(2,nres+i)
1440 c$$$        zi=c(3,nres+i)
1441 c$$$        dxi=dc_norm(1,nres+i)
1442 c$$$        dyi=dc_norm(2,nres+i)
1443 c$$$        dzi=dc_norm(3,nres+i)
1444 c$$$c        dsci_inv=dsc_inv(itypi)
1445 c$$$        dsci_inv=vbld_inv(i+nres)
1446 c$$$c        write (iout,*) "i",i,dsc_inv(itypi),dsci_inv,1.0d0/vbld(i+nres)
1447 c$$$c        write (iout,*) "dcnori",dxi*dxi+dyi*dyi+dzi*dzi
1448 c$$$C
1449 c$$$C Calculate SC interaction energy.
1450 c$$$C
1451 c$$$c$$$        do iint=1,nint_gr(i)
1452 c$$$c$$$          do j=istart(i,iint),iend(i,iint)
1453 c$$$        j=j_sc
1454 c$$$            ind=ind+1
1455 c$$$            itypj=itype(j)
1456 c$$$c            dscj_inv=dsc_inv(itypj)
1457 c$$$            dscj_inv=vbld_inv(j+nres)
1458 c$$$c            write (iout,*) "j",j,dsc_inv(itypj),dscj_inv,
1459 c$$$c     &       1.0d0/vbld(j+nres)
1460 c$$$c            write (iout,*) "i",i," j", j," itype",itype(i),itype(j)
1461 c$$$            sig0ij=sigma(itypi,itypj)
1462 c$$$            chi1=chi(itypi,itypj)
1463 c$$$            chi2=chi(itypj,itypi)
1464 c$$$            chi12=chi1*chi2
1465 c$$$            chip1=chip(itypi)
1466 c$$$            chip2=chip(itypj)
1467 c$$$            chip12=chip1*chip2
1468 c$$$            alf1=alp(itypi)
1469 c$$$            alf2=alp(itypj)
1470 c$$$            alf12=0.5D0*(alf1+alf2)
1471 c$$$C For diagnostics only!!!
1472 c$$$c           chi1=0.0D0
1473 c$$$c           chi2=0.0D0
1474 c$$$c           chi12=0.0D0
1475 c$$$c           chip1=0.0D0
1476 c$$$c           chip2=0.0D0
1477 c$$$c           chip12=0.0D0
1478 c$$$c           alf1=0.0D0
1479 c$$$c           alf2=0.0D0
1480 c$$$c           alf12=0.0D0
1481 c$$$            xj=c(1,nres+j)-xi
1482 c$$$            yj=c(2,nres+j)-yi
1483 c$$$            zj=c(3,nres+j)-zi
1484 c$$$            dxj=dc_norm(1,nres+j)
1485 c$$$            dyj=dc_norm(2,nres+j)
1486 c$$$            dzj=dc_norm(3,nres+j)
1487 c$$$c            write (iout,*) "dcnorj",dxi*dxi+dyi*dyi+dzi*dzi
1488 c$$$c            write (iout,*) "j",j," dc_norm",
1489 c$$$c     &       dc_norm(1,nres+j),dc_norm(2,nres+j),dc_norm(3,nres+j)
1490 c$$$            rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
1491 c$$$            rij=dsqrt(rrij)
1492 c$$$C Calculate angle-dependent terms of energy and contributions to their
1493 c$$$C derivatives.
1494 c$$$            call sc_angular
1495 c$$$            sigsq=1.0D0/sigsq
1496 c$$$            sig=sig0ij*dsqrt(sigsq)
1497 c$$$            rij_shift=1.0D0/rij-sig+sig0ij
1498 c$$$c for diagnostics; uncomment
1499 c$$$c            rij_shift=1.2*sig0ij
1500 c$$$C I hate to put IF's in the loops, but here don't have another choice!!!!
1501 c$$$            if (rij_shift.le.0.0D0) then
1502 c$$$              evdw=1.0D20
1503 c$$$cd              write (iout,'(2(a3,i3,2x),17(0pf7.3))')
1504 c$$$cd     &        restyp(itypi),i,restyp(itypj),j,
1505 c$$$cd     &        rij_shift,1.0D0/rij,sig,sig0ij,sigsq,1-dsqrt(sigsq) 
1506 c$$$              return
1507 c$$$            endif
1508 c$$$            sigder=-sig*sigsq
1509 c$$$c---------------------------------------------------------------
1510 c$$$            rij_shift=1.0D0/rij_shift 
1511 c$$$            fac=rij_shift**expon
1512 c$$$            e1=fac*fac*aa(itypi,itypj)
1513 c$$$            e2=fac*bb(itypi,itypj)
1514 c$$$            evdwij=eps1*eps2rt*eps3rt*(e1+e2)
1515 c$$$            eps2der=evdwij*eps3rt
1516 c$$$            eps3der=evdwij*eps2rt
1517 c$$$c            write (iout,*) "sigsq",sigsq," sig",sig," eps2rt",eps2rt,
1518 c$$$c     &        " eps3rt",eps3rt," eps1",eps1," e1",e1," e2",e2
1519 c$$$            evdwij=evdwij*eps2rt*eps3rt
1520 c$$$            evdw=evdw+evdwij
1521 c$$$            if (lprn) then
1522 c$$$            sigm=dabs(aa(itypi,itypj)/bb(itypi,itypj))**(1.0D0/6.0D0)
1523 c$$$            epsi=bb(itypi,itypj)**2/aa(itypi,itypj)
1524 c$$$            write (iout,'(2(a3,i3,2x),17(0pf7.3))')
1525 c$$$     &        restyp(itypi),i,restyp(itypj),j,
1526 c$$$     &        epsi,sigm,chi1,chi2,chip1,chip2,
1527 c$$$     &        eps1,eps2rt**2,eps3rt**2,sig,sig0ij,
1528 c$$$     &        om1,om2,om12,1.0D0/rij,1.0D0/rij_shift,
1529 c$$$     &        evdwij
1530 c$$$            endif
1531 c$$$
1532 c$$$            if (energy_dec) write (iout,'(a6,2i,0pf7.3)') 
1533 c$$$     &                        'evdw',i,j,evdwij
1534 c$$$
1535 c$$$C Calculate gradient components.
1536 c$$$            e1=e1*eps1*eps2rt**2*eps3rt**2
1537 c$$$            fac=-expon*(e1+evdwij)*rij_shift
1538 c$$$            sigder=fac*sigder
1539 c$$$            fac=rij*fac
1540 c$$$c            fac=0.0d0
1541 c$$$C Calculate the radial part of the gradient
1542 c$$$            gg(1)=xj*fac
1543 c$$$            gg(2)=yj*fac
1544 c$$$            gg(3)=zj*fac
1545 c$$$C Calculate angular part of the gradient.
1546 c$$$            call sc_grad
1547 c$$$c$$$          enddo      ! j
1548 c$$$c$$$        enddo        ! iint
1549 c$$$c$$$      enddo          ! i
1550 c$$$      energy_dec=.false.
1551 c$$$      return
1552 c$$$      end
1553 c$$$
1554 c$$$C-----------------------------------------------------------------------------
1555 c$$$
1556 c$$$      subroutine perturb_side_chain(i,angle)
1557 c$$$      implicit none
1558 c$$$
1559 c$$$c     Includes
1560 c$$$      include 'DIMENSIONS'
1561 c$$$      include 'COMMON.CHAIN'
1562 c$$$      include 'COMMON.GEO'
1563 c$$$      include 'COMMON.VAR'
1564 c$$$      include 'COMMON.LOCAL'
1565 c$$$      include 'COMMON.IOUNITS'
1566 c$$$
1567 c$$$c     External functions
1568 c$$$      external ran_number
1569 c$$$      double precision ran_number
1570 c$$$
1571 c$$$c     Input arguments
1572 c$$$      integer i
1573 c$$$      double precision angle    ! In degrees
1574 c$$$
1575 c$$$c     Local variables
1576 c$$$      integer i_sc
1577 c$$$      double precision rad_ang,rand_v(3),length,cost,sint
1578 c$$$
1579 c$$$
1580 c$$$      i_sc=i+nres
1581 c$$$      rad_ang=angle*deg2rad
1582 c$$$
1583 c$$$      length=0.0
1584 c$$$      do while (length.lt.0.01)
1585 c$$$        rand_v(1)=ran_number(0.01D0,1.0D0)
1586 c$$$        rand_v(2)=ran_number(0.01D0,1.0D0)
1587 c$$$        rand_v(3)=ran_number(0.01D0,1.0D0)
1588 c$$$        length=rand_v(1)*rand_v(1)+rand_v(2)*rand_v(2)+
1589 c$$$     +       rand_v(3)*rand_v(3)
1590 c$$$        length=sqrt(length)
1591 c$$$        rand_v(1)=rand_v(1)/length
1592 c$$$        rand_v(2)=rand_v(2)/length
1593 c$$$        rand_v(3)=rand_v(3)/length
1594 c$$$        cost=rand_v(1)*dc_norm(1,i_sc)+rand_v(2)*dc_norm(2,i_sc)+
1595 c$$$     +       rand_v(3)*dc_norm(3,i_sc)
1596 c$$$        length=1.0D0-cost*cost
1597 c$$$        if (length.lt.0.0D0) length=0.0D0
1598 c$$$        length=sqrt(length)
1599 c$$$        rand_v(1)=rand_v(1)-cost*dc_norm(1,i_sc)
1600 c$$$        rand_v(2)=rand_v(2)-cost*dc_norm(2,i_sc)
1601 c$$$        rand_v(3)=rand_v(3)-cost*dc_norm(3,i_sc)
1602 c$$$      enddo
1603 c$$$      rand_v(1)=rand_v(1)/length
1604 c$$$      rand_v(2)=rand_v(2)/length
1605 c$$$      rand_v(3)=rand_v(3)/length
1606 c$$$
1607 c$$$      cost=dcos(rad_ang)
1608 c$$$      sint=dsin(rad_ang)
1609 c$$$      dc(1,i_sc)=vbld(i_sc)*(dc_norm(1,i_sc)*cost+rand_v(1)*sint)
1610 c$$$      dc(2,i_sc)=vbld(i_sc)*(dc_norm(2,i_sc)*cost+rand_v(2)*sint)
1611 c$$$      dc(3,i_sc)=vbld(i_sc)*(dc_norm(3,i_sc)*cost+rand_v(3)*sint)
1612 c$$$      dc_norm(1,i_sc)=dc(1,i_sc)*vbld_inv(i_sc)
1613 c$$$      dc_norm(2,i_sc)=dc(2,i_sc)*vbld_inv(i_sc)
1614 c$$$      dc_norm(3,i_sc)=dc(3,i_sc)*vbld_inv(i_sc)
1615 c$$$      c(1,i_sc)=c(1,i)+dc(1,i_sc)
1616 c$$$      c(2,i_sc)=c(2,i)+dc(2,i_sc)
1617 c$$$      c(3,i_sc)=c(3,i)+dc(3,i_sc)
1618 c$$$
1619 c$$$      call chainbuild_cart
1620 c$$$
1621 c$$$      return
1622 c$$$      end
1623 c$$$
1624 c$$$c----------------------------------------------------------------------------
1625 c$$$
1626 c$$$      subroutine ss_relax3(i_in,j_in)
1627 c$$$      implicit none
1628 c$$$
1629 c$$$c     Includes
1630 c$$$      include 'DIMENSIONS'
1631 c$$$      include 'COMMON.VAR'
1632 c$$$      include 'COMMON.CHAIN'
1633 c$$$      include 'COMMON.IOUNITS'
1634 c$$$      include 'COMMON.INTERACT'
1635 c$$$
1636 c$$$c     External functions
1637 c$$$      external ran_number
1638 c$$$      double precision ran_number
1639 c$$$
1640 c$$$c     Input arguments
1641 c$$$      integer i_in,j_in
1642 c$$$
1643 c$$$c     Local variables
1644 c$$$      double precision energy_sc(0:n_ene),etot
1645 c$$$      double precision org_dc(3),org_dc_norm(3),org_c(3)
1646 c$$$      double precision ang_pert,rand_fact,exp_fact,beta
1647 c$$$      integer n,i_pert,i
1648 c$$$      logical notdone
1649 c$$$
1650 c$$$
1651 c$$$      beta=1.0D0
1652 c$$$
1653 c$$$      mask_r=.true.
1654 c$$$      do i=nnt,nct
1655 c$$$        mask_side(i)=0
1656 c$$$      enddo
1657 c$$$      mask_side(i_in)=1
1658 c$$$      mask_side(j_in)=1
1659 c$$$
1660 c$$$      call etotal_sc(energy_sc)
1661 c$$$      etot=energy_sc(0)
1662 c$$$c      write(iout,'(a,3d15.5)')"     SS_MC_START ",energy_sc(0),
1663 c$$$c     +     energy_sc(1),energy_sc(12)
1664 c$$$
1665 c$$$      notdone=.true.
1666 c$$$      n=0
1667 c$$$      do while (notdone)
1668 c$$$        if (mod(n,2).eq.0) then
1669 c$$$          i_pert=i_in
1670 c$$$        else
1671 c$$$          i_pert=j_in
1672 c$$$        endif
1673 c$$$        n=n+1
1674 c$$$
1675 c$$$        do i=1,3
1676 c$$$          org_dc(i)=dc(i,i_pert+nres)
1677 c$$$          org_dc_norm(i)=dc_norm(i,i_pert+nres)
1678 c$$$          org_c(i)=c(i,i_pert+nres)
1679 c$$$        enddo
1680 c$$$        ang_pert=ran_number(0.0D0,3.0D0)
1681 c$$$        call perturb_side_chain(i_pert,ang_pert)
1682 c$$$        call etotal_sc(energy_sc)
1683 c$$$        exp_fact=exp(beta*(etot-energy_sc(0)))
1684 c$$$        rand_fact=ran_number(0.0D0,1.0D0)
1685 c$$$        if (rand_fact.lt.exp_fact) then
1686 c$$$c          write(iout,'(a,3d15.5)')"     SS_MC_ACCEPT ",energy_sc(0),
1687 c$$$c     +     energy_sc(1),energy_sc(12)
1688 c$$$          etot=energy_sc(0)
1689 c$$$        else
1690 c$$$c          write(iout,'(a,3d15.5)')"     SS_MC_REJECT ",energy_sc(0),
1691 c$$$c     +     energy_sc(1),energy_sc(12)
1692 c$$$          do i=1,3
1693 c$$$            dc(i,i_pert+nres)=org_dc(i)
1694 c$$$            dc_norm(i,i_pert+nres)=org_dc_norm(i)
1695 c$$$            c(i,i_pert+nres)=org_c(i)
1696 c$$$          enddo
1697 c$$$        endif
1698 c$$$
1699 c$$$        if (n.eq.10000.or.etot.lt.30.0D0) notdone=.false.
1700 c$$$      enddo
1701 c$$$
1702 c$$$      mask_r=.false.
1703 c$$$
1704 c$$$      return
1705 c$$$      end
1706 c$$$
1707 c$$$c----------------------------------------------------------------------------
1708 c$$$
1709 c$$$      subroutine ss_relax2(etot,iretcode,nfun,i_in,j_in)
1710 c$$$      implicit none
1711 c$$$      include 'DIMENSIONS'
1712 c$$$      integer liv,lv
1713 c$$$      parameter (liv=60,lv=(77+maxres6*(maxres6+17)/2)) 
1714 c$$$*********************************************************************
1715 c$$$* OPTIMIZE sets up SUMSL or DFP and provides a simple interface for *
1716 c$$$* the calling subprogram.                                           *     
1717 c$$$* when d(i)=1.0, then v(35) is the length of the initial step,      *     
1718 c$$$* calculated in the usual pythagorean way.                          *     
1719 c$$$* absolute convergence occurs when the function is within v(31) of  *     
1720 c$$$* zero. unless you know the minimum value in advance, abs convg     *     
1721 c$$$* is probably not useful.                                           *     
1722 c$$$* relative convergence is when the model predicts that the function *   
1723 c$$$* will decrease by less than v(32)*abs(fun).                        *   
1724 c$$$*********************************************************************
1725 c$$$      include 'COMMON.IOUNITS'
1726 c$$$      include 'COMMON.VAR'
1727 c$$$      include 'COMMON.GEO'
1728 c$$$      include 'COMMON.MINIM'
1729 c$$$      include 'COMMON.CHAIN'
1730 c$$$
1731 c$$$      double precision orig_ss_dc,orig_ss_var,orig_ss_dist
1732 c$$$      common /orig_ss/ orig_ss_dc(3,0:maxres2),orig_ss_var(maxvar),
1733 c$$$     +     orig_ss_dist(maxres2,maxres2)
1734 c$$$
1735 c$$$      double precision etot
1736 c$$$      integer iretcode,nfun,i_in,j_in
1737 c$$$
1738 c$$$      external dist
1739 c$$$      double precision dist
1740 c$$$      external ss_func,fdum
1741 c$$$      double precision ss_func,fdum
1742 c$$$
1743 c$$$      integer iv(liv),uiparm(2)
1744 c$$$      double precision v(lv),x(maxres6),d(maxres6),rdum
1745 c$$$      integer i,j,k
1746 c$$$
1747 c$$$
1748 c$$$      call deflt(2,iv,liv,lv,v)                                         
1749 c$$$* 12 means fresh start, dont call deflt                                 
1750 c$$$      iv(1)=12                                                          
1751 c$$$* max num of fun calls                                                  
1752 c$$$      if (maxfun.eq.0) maxfun=500
1753 c$$$      iv(17)=maxfun
1754 c$$$* max num of iterations                                                 
1755 c$$$      if (maxmin.eq.0) maxmin=1000
1756 c$$$      iv(18)=maxmin
1757 c$$$* controls output                                                       
1758 c$$$      iv(19)=2                                                          
1759 c$$$* selects output unit                                                   
1760 c$$$c      iv(21)=iout                                                       
1761 c$$$      iv(21)=0
1762 c$$$* 1 means to print out result                                           
1763 c$$$      iv(22)=0                                                          
1764 c$$$* 1 means to print out summary stats                                    
1765 c$$$      iv(23)=0                                                          
1766 c$$$* 1 means to print initial x and d                                      
1767 c$$$      iv(24)=0                                                          
1768 c$$$* min val for v(radfac) default is 0.1                                  
1769 c$$$      v(24)=0.1D0                                                       
1770 c$$$* max val for v(radfac) default is 4.0                                  
1771 c$$$      v(25)=2.0D0                                                       
1772 c$$$c     v(25)=4.0D0                                                       
1773 c$$$* check false conv if (act fnctn decrease) .lt. v(26)*(exp decrease)    
1774 c$$$* the sumsl default is 0.1                                              
1775 c$$$      v(26)=0.1D0
1776 c$$$* false conv if (act fnctn decrease) .lt. v(34)                         
1777 c$$$* the sumsl default is 100*machep                                       
1778 c$$$      v(34)=v(34)/100.0D0                                               
1779 c$$$* absolute convergence                                                  
1780 c$$$      if (tolf.eq.0.0D0) tolf=1.0D-4
1781 c$$$      v(31)=tolf
1782 c$$$      v(31)=1.0D-1
1783 c$$$* relative convergence                                                  
1784 c$$$      if (rtolf.eq.0.0D0) rtolf=1.0D-4
1785 c$$$      v(32)=rtolf
1786 c$$$      v(32)=1.0D-1
1787 c$$$* controls initial step size                                            
1788 c$$$      v(35)=1.0D-1
1789 c$$$* large vals of d correspond to small components of step                
1790 c$$$      do i=1,6*nres
1791 c$$$        d(i)=1.0D0
1792 c$$$      enddo
1793 c$$$
1794 c$$$      do i=0,2*nres
1795 c$$$        do j=1,3
1796 c$$$          orig_ss_dc(j,i)=dc(j,i)
1797 c$$$        enddo
1798 c$$$      enddo
1799 c$$$      call geom_to_var(nvar,orig_ss_var)
1800 c$$$
1801 c$$$      do i=1,nres
1802 c$$$        do j=i,nres
1803 c$$$          orig_ss_dist(j,i)=dist(j,i)
1804 c$$$          orig_ss_dist(j+nres,i)=dist(j+nres,i)
1805 c$$$          orig_ss_dist(j,i+nres)=dist(j,i+nres)
1806 c$$$          orig_ss_dist(j+nres,i+nres)=dist(j+nres,i+nres)
1807 c$$$        enddo
1808 c$$$      enddo
1809 c$$$
1810 c$$$      k=0
1811 c$$$      do i=1,nres-1
1812 c$$$        do j=1,3
1813 c$$$          k=k+1
1814 c$$$          x(k)=dc(j,i)
1815 c$$$        enddo
1816 c$$$      enddo
1817 c$$$      do i=2,nres-1
1818 c$$$        if (ialph(i,1).gt.0) then
1819 c$$$        do j=1,3
1820 c$$$          k=k+1
1821 c$$$          x(k)=dc(j,i+nres)
1822 c$$$        enddo
1823 c$$$        endif
1824 c$$$      enddo
1825 c$$$
1826 c$$$      uiparm(1)=i_in
1827 c$$$      uiparm(2)=j_in
1828 c$$$      call smsno(k,d,x,ss_func,iv,liv,lv,v,uiparm,rdum,fdum)
1829 c$$$      etot=v(10)
1830 c$$$      iretcode=iv(1)
1831 c$$$      nfun=iv(6)+iv(30)
1832 c$$$
1833 c$$$      k=0
1834 c$$$      do i=1,nres-1
1835 c$$$        do j=1,3
1836 c$$$          k=k+1
1837 c$$$          dc(j,i)=x(k)
1838 c$$$        enddo
1839 c$$$      enddo
1840 c$$$      do i=2,nres-1
1841 c$$$        if (ialph(i,1).gt.0) then
1842 c$$$        do j=1,3
1843 c$$$          k=k+1
1844 c$$$          dc(j,i+nres)=x(k)
1845 c$$$        enddo
1846 c$$$        endif
1847 c$$$      enddo
1848 c$$$      call chainbuild_cart
1849 c$$$
1850 c$$$      return  
1851 c$$$      end  
1852 c$$$
1853 c$$$C-----------------------------------------------------------------------------
1854 c$$$
1855 c$$$      subroutine ss_func(n,x,nf,f,uiparm,urparm,ufparm)  
1856 c$$$      implicit none
1857 c$$$      include 'DIMENSIONS'
1858 c$$$      include 'COMMON.DERIV'
1859 c$$$      include 'COMMON.IOUNITS'
1860 c$$$      include 'COMMON.VAR'
1861 c$$$      include 'COMMON.CHAIN'
1862 c$$$      include 'COMMON.INTERACT'
1863 c$$$      include 'COMMON.SBRIDGE'
1864 c$$$
1865 c$$$      double precision orig_ss_dc,orig_ss_var,orig_ss_dist
1866 c$$$      common /orig_ss/ orig_ss_dc(3,0:maxres2),orig_ss_var(maxvar),
1867 c$$$     +     orig_ss_dist(maxres2,maxres2)
1868 c$$$
1869 c$$$      integer n
1870 c$$$      double precision x(maxres6)
1871 c$$$      integer nf
1872 c$$$      double precision f
1873 c$$$      integer uiparm(2)
1874 c$$$      real*8 urparm(1)
1875 c$$$      external ufparm
1876 c$$$      double precision ufparm
1877 c$$$
1878 c$$$      external dist
1879 c$$$      double precision dist
1880 c$$$
1881 c$$$      integer i,j,k,ss_i,ss_j
1882 c$$$      double precision tempf,var(maxvar)
1883 c$$$
1884 c$$$
1885 c$$$      ss_i=uiparm(1)
1886 c$$$      ss_j=uiparm(2)
1887 c$$$      f=0.0D0
1888 c$$$
1889 c$$$      k=0
1890 c$$$      do i=1,nres-1
1891 c$$$        do j=1,3
1892 c$$$          k=k+1
1893 c$$$          dc(j,i)=x(k)
1894 c$$$        enddo
1895 c$$$      enddo
1896 c$$$      do i=2,nres-1
1897 c$$$        if (ialph(i,1).gt.0) then
1898 c$$$        do j=1,3
1899 c$$$          k=k+1
1900 c$$$          dc(j,i+nres)=x(k)
1901 c$$$        enddo
1902 c$$$        endif
1903 c$$$      enddo
1904 c$$$      call chainbuild_cart
1905 c$$$
1906 c$$$      call geom_to_var(nvar,var)
1907 c$$$
1908 c$$$c     Constraints on all angles
1909 c$$$      do i=1,nvar
1910 c$$$        tempf=var(i)-orig_ss_var(i)
1911 c$$$        f=f+tempf*tempf
1912 c$$$      enddo
1913 c$$$
1914 c$$$c     Constraints on all distances
1915 c$$$      do i=1,nres-1
1916 c$$$        if (i.gt.1) then
1917 c$$$          tempf=dist(i+nres,i)-orig_ss_dist(i+nres,i)
1918 c$$$          f=f+tempf*tempf
1919 c$$$        endif
1920 c$$$        do j=i+1,nres
1921 c$$$          tempf=dist(j,i)-orig_ss_dist(j,i)
1922 c$$$          if (tempf.lt.0.0D0 .or. j.eq.i+1) f=f+tempf*tempf
1923 c$$$          tempf=dist(j+nres,i)-orig_ss_dist(j+nres,i)
1924 c$$$          if (tempf.lt.0.0D0) f=f+tempf*tempf
1925 c$$$          tempf=dist(j,i+nres)-orig_ss_dist(j,i+nres)
1926 c$$$          if (tempf.lt.0.0D0) f=f+tempf*tempf
1927 c$$$          tempf=dist(j+nres,i+nres)-orig_ss_dist(j+nres,i+nres)
1928 c$$$          if (tempf.lt.0.0D0) f=f+tempf*tempf
1929 c$$$        enddo
1930 c$$$      enddo
1931 c$$$
1932 c$$$c     Constraints for the relevant CYS-CYS
1933 c$$$      tempf=dist(nres+ss_i,nres+ss_j)-8.0D0
1934 c$$$      f=f+tempf*tempf
1935 c$$$CCCCCCCCCCCCCCCCC      ADD SOME ANGULAR STUFF
1936 c$$$
1937 c$$$c$$$      if (nf.ne.nfl) then
1938 c$$$c$$$        write(iout,'(a,i10,2d15.5)')"IN DIST_FUNC (NF,F,DIST)",nf,
1939 c$$$c$$$     +       f,dist(5+nres,14+nres)
1940 c$$$c$$$      endif
1941 c$$$
1942 c$$$      nfl=nf
1943 c$$$
1944 c$$$      return                                                            
1945 c$$$      end                                                               
1946 c$$$
1947 c$$$C-----------------------------------------------------------------------------