Introduction of PBC to cluster energy and reading, turn off wham one of the debug
authorAdam Sieradzan <adasko@piasek4.chem.univ.gda.pl>
Tue, 10 Jun 2014 16:03:57 +0000 (18:03 +0200)
committerAdam Sieradzan <adasko@piasek4.chem.univ.gda.pl>
Tue, 10 Jun 2014 16:03:57 +0000 (18:03 +0200)
source/cluster/wham/src-M/energy_p_new.F
source/cluster/wham/src-M/readrtns.F
source/wham/src-M/energy_p_new.F

index c588bee..d5ccc6d 100644 (file)
@@ -754,6 +754,7 @@ C
       common /srutu/icall
       integer icant
       external icant
+      integer xshift,yshift,zshift
 c     print *,'Entering EGB nnt=',nnt,' nct=',nct,' expon=',expon
       evdw=0.0D0
       evdw_t=0.0d0
@@ -767,6 +768,12 @@ c      if (icall.gt.0) lprn=.true.
         xi=c(1,nres+i)
         yi=c(2,nres+i)
         zi=c(3,nres+i)
+          xi=mod(xi,boxxsize)
+          if (xi.lt.0) xi=xi+boxxsize
+          yi=mod(yi,boxysize)
+          if (yi.lt.0) yi=yi+boxysize
+          zi=mod(zi,boxzsize)
+          if (zi.lt.0) zi=zi+boxzsize
         dxi=dc_norm(1,nres+i)
         dyi=dc_norm(2,nres+i)
         dzi=dc_norm(3,nres+i)
@@ -800,15 +807,55 @@ c           chip12=0.0D0
 c           alf1=0.0D0
 c           alf2=0.0D0
 c           alf12=0.0D0
-            xj=c(1,nres+j)-xi
-            yj=c(2,nres+j)-yi
-            zj=c(3,nres+j)-zi
+            xj=c(1,nres+j)
+            yj=c(2,nres+j)
+            zj=c(3,nres+j)
+          xj=mod(xj,boxxsize)
+          if (xj.lt.0) xj=xj+boxxsize
+          yj=mod(yj,boxysize)
+          if (yj.lt.0) yj=yj+boxysize
+          zj=mod(zj,boxzsize)
+          if (zj.lt.0) zj=zj+boxzsize
+      dist_init=(xj-xi)**2+(yj-yi)**2+(zj-zi)**2
+      xj_safe=xj
+      yj_safe=yj
+      zj_safe=zj
+      subchap=0
+      do xshift=-1,1
+      do yshift=-1,1
+      do zshift=-1,1
+          xj=xj_safe+xshift*boxxsize
+          yj=yj_safe+yshift*boxysize
+          zj=zj_safe+zshift*boxzsize
+          dist_temp=(xj-xi)**2+(yj-yi)**2+(zj-zi)**2
+          if(dist_temp.lt.dist_init) then
+            dist_init=dist_temp
+            xj_temp=xj
+            yj_temp=yj
+            zj_temp=zj
+            subchap=1
+          endif
+       enddo
+       enddo
+       enddo
+       if (subchap.eq.1) then
+          xj=xj_temp-xi
+          yj=yj_temp-yi
+          zj=zj_temp-zi
+       else
+          xj=xj_safe-xi
+          yj=yj_safe-yi
+          zj=zj_safe-zi
+       endif
             dxj=dc_norm(1,nres+j)
             dyj=dc_norm(2,nres+j)
             dzj=dc_norm(3,nres+j)
 c            write (iout,*) i,j,xj,yj,zj
             rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
             rij=dsqrt(rrij)
+            sss=sscale((1.0d0/rij)/sigma(itypi,itypj))
+            sssgrad=sscagrad((1.0d0/rij)/sigma(itypi,itypj))
+            if (sss.le.0.0d0) cycle
 C Calculate angle-dependent terms of energy and contributions to their
 C derivatives.
             call sc_angular
@@ -831,9 +878,9 @@ c---------------------------------------------------------------
             eps3der=evdwij*eps2rt
             evdwij=evdwij*eps2rt*eps3rt
             if (bb(itypi,itypj).gt.0) then
-              evdw=evdw+evdwij
+              evdw=evdw+evdwij*sss
             else
-              evdw_t=evdw_t+evdwij
+              evdw_t=evdw_t+evdwij*sss
             endif
             ij=icant(itypi,itypj)
             aux=eps1*eps2rt**2*eps3rt**2
@@ -857,6 +904,7 @@ C Calculate gradient components.
             fac=-expon*(e1+evdwij)*rij_shift
             sigder=fac*sigder
             fac=rij*fac
+            fac=fac+evdwij/sss*sssgrad/sigma(itypi,itypj)*rij
 C Calculate the radial part of the gradient
             gg(1)=xj*fac
             gg(2)=yj*fac
@@ -1806,7 +1854,15 @@ cd      write (iout,*) 'iatel_s=',iatel_s,' iatel_e=',iatel_e
         gcorr_loc(i)=0.0d0
       enddo
       do i=iatel_s,iatel_e
-        if (itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1) cycle
+          if (i.eq.1) then
+           if (itype(i).eq.ntyp1.or. itype(i+1).eq.ntyp1
+     &  .or. itype(i+2).eq.ntyp1) cycle
+          else
+        if (itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1
+     &  .or. itype(i+2).eq.ntyp1
+     &  .or. itype(i-1).eq.ntyp1
+     &) cycle
+         endif
         if (itel(i).eq.0) goto 1215
         dxi=dc(1,i)
         dyi=dc(2,i)
@@ -1817,10 +1873,25 @@ cd      write (iout,*) 'iatel_s=',iatel_s,' iatel_e=',iatel_e
         xmedi=c(1,i)+0.5d0*dxi
         ymedi=c(2,i)+0.5d0*dyi
         zmedi=c(3,i)+0.5d0*dzi
+          xmedi=mod(xmedi,boxxsize)
+          if (xmedi.lt.0) xmedi=xmedi+boxxsize
+          ymedi=mod(ymedi,boxysize)
+          if (ymedi.lt.0) ymedi=ymedi+boxysize
+          zmedi=mod(zmedi,boxzsize)
+          if (zmedi.lt.0) zmedi=zmedi+boxzsize
         num_conti=0
 c        write (iout,*) 'i',i,' ielstart',ielstart(i),' ielend',ielend(i)
         do j=ielstart(i),ielend(i)
-          if (itype(j).eq.ntyp1 .or. itype(j+1).eq.ntyp1) cycle
+          if (j.eq.1) then
+           if (itype(j).eq.ntyp1 .or. itype(j+1).eq.ntyp1
+     & .or.itype(j+2).eq.ntyp1
+     &) cycle
+          else
+          if (itype(j).eq.ntyp1 .or. itype(j+1).eq.ntyp1
+     & .or.itype(j+2).eq.ntyp1
+     & .or.itype(j-1).eq.ntyp1
+     &) cycle
+         endif
           if (itel(j).eq.0) goto 1216
           ind=ind+1
           iteli=itel(i)
@@ -1842,10 +1913,50 @@ C End diagnostics
           dx_normj=dc_norm(1,j)
           dy_normj=dc_norm(2,j)
           dz_normj=dc_norm(3,j)
-          xj=c(1,j)+0.5D0*dxj-xmedi
-          yj=c(2,j)+0.5D0*dyj-ymedi
-          zj=c(3,j)+0.5D0*dzj-zmedi
+          xj=c(1,j)+0.5D0*dxj
+          yj=c(2,j)+0.5D0*dyj
+          zj=c(3,j)+0.5D0*dzj
+         xj=mod(xj,boxxsize)
+          if (xj.lt.0) xj=xj+boxxsize
+          yj=mod(yj,boxysize)
+          if (yj.lt.0) yj=yj+boxysize
+          zj=mod(zj,boxzsize)
+          if (zj.lt.0) zj=zj+boxzsize
+      dist_init=(xj-xi)**2+(yj-yi)**2+(zj-zi)**2
+      xj_safe=xj
+      yj_safe=yj
+      zj_safe=zj
+      isubchap=0
+      do xshift=-1,1
+      do yshift=-1,1
+      do zshift=-1,1
+          xj=xj_safe+xshift*boxxsize
+          yj=yj_safe+yshift*boxysize
+          zj=zj_safe+zshift*boxzsize
+          dist_temp=(xj-xi)**2+(yj-yi)**2+(zj-zi)**2
+          if(dist_temp.lt.dist_init) then
+            dist_init=dist_temp
+            xj_temp=xj
+            yj_temp=yj
+            zj_temp=zj
+            isubchap=1
+          endif
+       enddo
+       enddo
+       enddo
+       if (isubchap.eq.1) then
+          xj=xj_temp-xmedi
+          yj=yj_temp-ymedi
+          zj=zj_temp-zmedi
+       else
+          xj=xj_safe-xmedi
+          yj=yj_safe-ymedi
+          zj=zj_safe-zmedi
+       endif
+
           rij=xj*xj+yj*yj+zj*zj
+            sss=sscale(sqrt(rij))
+            sssgrad=sscagrad(sqrt(rij))
           rrmij=1.0D0/rij
           rij=dsqrt(rij)
           rmij=1.0D0/rij
@@ -1869,7 +1980,7 @@ c          write (iout,*) "i",i,iteli," j",j,itelj," eesij",eesij
 C 12/26/95 - for the evaluation of multi-body H-bonding interactions
           ees0ij=4.0D0+fac*fac-3.0D0*(cosb*cosb+cosg*cosg)
           ees=ees+eesij
-          evdw1=evdw1+evdwij
+          evdw1=evdw1+evdwij*sss
 cd          write(iout,'(2(2i3,2x),7(1pd12.4)/2(3(1pd12.4),5x)/)')
 cd     &      iteli,i,itelj,j,aaa,bbb,ael6i,ael3i,
 cd     &      1.0D0/dsqrt(rrmij),evdwij,eesij,
@@ -1878,7 +1989,7 @@ C
 C Calculate contributions to the Cartesian gradient.
 C
 #ifdef SPLITELE
-          facvdw=-6*rrmij*(ev1+evdwij) 
+          facvdw=-6*rrmij*(ev1+evdwij)*sss
           facel=-3*rrmij*(el1+eesij)
           fac1=fac
           erij(1)=xj*rmij
@@ -1904,9 +2015,18 @@ C
               gelc(l,k)=gelc(l,k)+ggg(l)
             enddo
           enddo
-          ggg(1)=facvdw*xj
-          ggg(2)=facvdw*yj
-          ggg(3)=facvdw*zj
+C          ggg(1)=facvdw*xj
+C          ggg(2)=facvdw*yj
+C          ggg(3)=facvdw*zj
+          if (sss.gt.0.0) then
+          ggg(1)=facvdw*xj+sssgrad*rmij*evdwij*xj
+          ggg(2)=facvdw*yj+sssgrad*rmij*evdwij*yj
+          ggg(3)=facvdw*zj+sssgrad*rmij*evdwij*zj
+          else
+          ggg(1)=0.0
+          ggg(2)=0.0
+          ggg(3)=0.0
+          endif
           do k=1,3
             ghalf=0.5D0*ggg(k)
             gvdwpp(k,i)=gvdwpp(k,i)+ghalf
@@ -1921,7 +2041,7 @@ C
             enddo
           enddo
 #else
-          facvdw=ev1+evdwij 
+          facvdw=(ev1+evdwij)*sss
           facel=el1+eesij  
           fac1=fac
           fac=-3*rrmij*(facvdw+facvdw+facel)
@@ -2765,6 +2885,13 @@ c     &   " iscp",(iscpstart(i,j),iscpend(i,j),j=1,nscp_gr(i))
         xi=0.5D0*(c(1,i)+c(1,i+1))
         yi=0.5D0*(c(2,i)+c(2,i+1))
         zi=0.5D0*(c(3,i)+c(3,i+1))
+C    Returning the ith atom to box
+          xi=mod(xi,boxxsize)
+          if (xi.lt.0) xi=xi+boxxsize
+          yi=mod(yi,boxysize)
+          if (yi.lt.0) yi=yi+boxysize
+          zi=mod(zi,boxzsize)
+          if (zi.lt.0) zi=zi+boxzsize
 
         do iint=1,nscp_gr(i)
 
@@ -2776,26 +2903,72 @@ c         xj=c(1,nres+j)-xi
 c         yj=c(2,nres+j)-yi
 c         zj=c(3,nres+j)-zi
 C Uncomment following three lines for Ca-p interactions
-          xj=c(1,j)-xi
-          yj=c(2,j)-yi
-          zj=c(3,j)-zi
+          xj=c(1,j)
+          yj=c(2,j)
+          zj=c(3,j)
+C returning the jth atom to box
+          xj=mod(xj,boxxsize)
+          if (xj.lt.0) xj=xj+boxxsize
+          yj=mod(yj,boxysize)
+          if (yj.lt.0) yj=yj+boxysize
+          zj=mod(zj,boxzsize)
+          if (zj.lt.0) zj=zj+boxzsize
+      dist_init=(xj-xi)**2+(yj-yi)**2+(zj-zi)**2
+      xj_safe=xj
+      yj_safe=yj
+      zj_safe=zj
+      subchap=0
+C Finding the closest jth atom
+      do xshift=-1,1
+      do yshift=-1,1
+      do zshift=-1,1
+          xj=xj_safe+xshift*boxxsize
+          yj=yj_safe+yshift*boxysize
+          zj=zj_safe+zshift*boxzsize
+          dist_temp=(xj-xi)**2+(yj-yi)**2+(zj-zi)**2
+          if(dist_temp.lt.dist_init) then
+            dist_init=dist_temp
+            xj_temp=xj
+            yj_temp=yj
+            zj_temp=zj
+            subchap=1
+          endif
+       enddo
+       enddo
+       enddo
+       if (subchap.eq.1) then
+          xj=xj_temp-xi
+          yj=yj_temp-yi
+          zj=zj_temp-zi
+       else
+          xj=xj_safe-xi
+          yj=yj_safe-yi
+          zj=zj_safe-zi
+       endif
+
           rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
+C sss is scaling function for smoothing the cutoff gradient otherwise
+C the gradient would not be continuouse
+          sss=sscale(1.0d0/(dsqrt(rrij)))
+          if (sss.le.0.0d0) cycle
+          sssgrad=sscagrad(1.0d0/(dsqrt(rrij)))
           fac=rrij**expon2
           e1=fac*fac*aad(itypj,iteli)
           e2=fac*bad(itypj,iteli)
           if (iabs(j-i) .le. 2) then
             e1=scal14*e1
             e2=scal14*e2
-            evdw2_14=evdw2_14+e1+e2
+            evdw2_14=evdw2_14+(e1+e2)*sss
           endif
           evdwij=e1+e2
 c          write (iout,*) i,j,evdwij
-          evdw2=evdw2+evdwij
+          evdw2=evdw2+evdwij*sss
           if (calc_grad) then
 C
 C Calculate contributions to the gradient in the virtual-bond and SC vectors.
 C
-          fac=-(evdwij+e1)*rrij
+           fac=-(evdwij+e1)*rrij*sss
+           fac=fac+(evdwij)*sssgrad*dsqrt(rrij)/expon
           ggg(1)=xj*fac
           ggg(2)=yj*fac
           ggg(3)=zj*fac
       estr=0.0d0
       estr1=0.0d0
       do i=nnt+1,nct
-        if (itype(i-1).eq.ntyp1 .or. itype(i).eq.ntyp1) then
-          estr1=estr1+gnmr1(vbld(i),-1.0d0,distchainmax)
-          do j=1,3
-          gradb(j,i-1)=gnmr1prim(vbld(i),-1.0d0,distchainmax)
-     &      *dc(j,i-1)/vbld(i)
-          enddo
-          if (energy_dec) write(iout,*)
-     &       "estr1",i,gnmr1(vbld(i),-1.0d0,distchainmax)
-        else
+        if (itype(i-1).eq.ntyp1 .and. itype(i).eq.ntyp1) cycle
+C          estr1=estr1+gnmr1(vbld(i),-1.0d0,distchainmax)
+C          do j=1,3
+C          gradb(j,i-1)=gnmr1prim(vbld(i),-1.0d0,distchainmax)
+C     &      *dc(j,i-1)/vbld(i)
+C          enddo
+C          if (energy_dec) write(iout,*)
+C     &       "estr1",i,vbld(i),distchainmax,
+C     &       gnmr1(vbld(i),-1.0d0,distchainmax)
+C        else
+         if (itype(i-1).eq.ntyp1 .or. itype(i).eq.ntyp1) then
+        diff = vbld(i)-vbldpDUM
+         else
           diff = vbld(i)-vbldp0
 c          write (iout,*) i,vbld(i),vbldp0,diff,AKP*diff*diff
+         endif
           estr=estr+diff*diff
           do j=1,3
             gradb(j,i-1)=AKP*diff*dc(j,i-1)/vbld(i)
           enddo
-        endif
-
+C        endif
+C        write (iout,'(a7,i5,4f7.3)')
+C     &     "estr bb",i,vbld(i),vbldp0,diff,AKP*diff*diff
       enddo
       estr=0.5d0*AKP*estr+estr1
 c
@@ -3129,7 +3308,9 @@ c      write (iout,*) "nres",nres
 c     write (*,'(a,i2)') 'EBEND ICG=',icg
 c      write (iout,*) ithet_start,ithet_end
       do i=ithet_start,ithet_end
-        if (itype(i-1).eq.ntyp1) cycle
+        if (i.le.2) cycle
+        if ((itype(i-1).eq.ntyp1).or.itype(i-2).eq.ntyp1
+     &  .or.itype(i).eq.ntyp1) cycle
 C Zero the energy function and its derivative at 0 or pi.
         call splinthet(theta(i),0.5d0*delta,ss,ssd)
         it=itype(i-1)
@@ -3145,7 +3326,11 @@ C Zero the energy function and its derivative at 0 or pi.
           ichir21=isign(1,itype(i))
           ichir22=isign(1,itype(i))
          endif
-        if (i.gt.3 .and. itype(i-2).ne.ntyp1) then
+         if (i.eq.3) then
+          y(1)=0.0D0
+          y(2)=0.0D0
+          else
+        if (i.gt.3 .and. itype(i-3).ne.ntyp1) then
 #ifdef OSF
           phii=phi(i)
           icrc=0
@@ -3160,7 +3345,8 @@ C Zero the energy function and its derivative at 0 or pi.
           y(1)=0.0D0
           y(2)=0.0D0
         endif
-        if (i.lt.nres .and. itype(i).ne.ntyp1) then
+        endif
+        if (i.lt.nres .and. itype(i+1).ne.ntyp1) then
 #ifdef OSF
           phii1=phi(i+1)
           icrc=0
       etheta=0.0D0
 c      write (iout,*) "ithetyp",(ithetyp(i),i=1,ntyp1)
       do i=ithet_start,ithet_end
-        if (itype(i-1).eq.ntyp1) cycle
+        if (i.le.2) cycle
+        if ((itype(i-1).eq.ntyp1).or.itype(i-2).eq.ntyp1
+     &  .or.itype(i).eq.ntyp1) cycle
         if (iabs(itype(i+1)).eq.20) iblock=2
         if (iabs(itype(i+1)).ne.20) iblock=1
         dethetai=0.0d0
         dephii=0.0d0
         dephii1=0.0d0
         theti2=0.5d0*theta(i)
-CC Ta zmina jest niewlasciwa
         ityp2=ithetyp((itype(i-1)))
         do k=1,nntheterm
           coskt(k)=dcos(k*theti2)
           sinkt(k)=dsin(k*theti2)
         enddo
-        if (i.gt.3 .and. itype(i-2).ne.ntyp1) then
+        if (i.eq.3) then
+          phii=0.0d0
+          ityp1=nthetyp+1
+          do k=1,nsingle
+            cosph1(k)=0.0d0
+            sinph1(k)=0.0d0
+          enddo
+        else
+        if (i.gt.3 .and. itype(i-3).ne.ntyp1) then
 #ifdef OSF
           phii=phi(i)
           if (phii.ne.phii) phii=150.0
@@ -3408,7 +3603,8 @@ CC Ta zmina jest niewlasciwa
             sinph1(k)=0.0d0
           enddo 
         endif
-        if (i.lt.nres .and. itype(i).ne.ntyp1) then
+        endif
+        if (i.lt.nres .and. itype(i+1).ne.ntyp1) then
 #ifdef OSF
           phii1=phi(i+1)
           if (phii1.ne.phii1) phii1=150.0
@@ -4342,8 +4538,9 @@ C Set lprn=.true. for debugging
 c      lprn=.true.
       etors=0.0D0
       do i=iphi_start,iphi_end
-        if (itype(i-2).eq.ntyp1 .or. itype(i-1).eq.ntyp1
-     &       .or. itype(i).eq.ntyp1) cycle
+        if (i.le.2) cycle
+        if (itype(i-2).eq.ntyp1.or. itype(i-1).eq.ntyp1
+     &      .or. itype(i).eq.ntyp1 .or. itype(i-3).eq.ntyp1) cycle
         if (itel(i-2).eq.0 .or. itel(i-1).eq.0) goto 1215
          if (iabs(itype(i)).eq.20) then
          iblock=2
@@ -4441,8 +4638,10 @@ C Set lprn=.true. for debugging
 c     lprn=.true.
       etors_d=0.0D0
       do i=iphi_start,iphi_end-1
-        if (itype(i-2).eq.ntyp1 .or. itype(i-1).eq.ntyp1
-     &      .or. itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1) cycle
+        if (i.le.3) cycle
+         if ((itype(i-2).eq.ntyp1).or.itype(i-3).eq.ntyp1.or.
+     &  (itype(i-1).eq.ntyp1).or.(itype(i).eq.ntyp1).or.
+     &  (itype(i+1).eq.ntyp1)) cycle
         if (itel(i-2).eq.0 .or. itel(i-1).eq.0 .or. itel(i).eq.0) 
      &     goto 1215
         itori=itortyp(itype(i-2))
@@ -4698,9 +4897,9 @@ c------------------------------------------------------------------------------
       integer dimen1,dimen2,atom,indx
       double precision buffer(dimen1,dimen2)
       double precision zapas 
-      common /contacts_hb/ zapas(3,20,maxres,7),
-     &         facont_hb(20,maxres),ees0p(20,maxres),ees0m(20,maxres),
-     &         num_cont_hb(maxres),jcont_hb(20,maxres)
+      common /contacts_hb/ zapas(3,ntyp,maxres,7),
+     &     facont_hb(ntyp,maxres),ees0p(ntyp,maxres),ees0m(ntyp,maxres),
+     &         num_cont_hb(maxres),jcont_hb(ntyp,maxres)
       num_kont=buffer(1,indx+26)
       num_kont_old=num_cont_hb(atom)
       num_cont_hb(atom)=num_kont+num_kont_old
@@ -7445,4 +7644,34 @@ C-----------------------------------------------------------------------------
       scalar=sc
       return
       end
+C-----------------------------------------------------------------------
+      double precision function sscale(r)
+      double precision r,gamm
+      include "COMMON.SPLITELE"
+      if(r.lt.r_cut-rlamb) then
+        sscale=1.0d0
+      else if(r.le.r_cut.and.r.ge.r_cut-rlamb) then
+        gamm=(r-(r_cut-rlamb))/rlamb
+        sscale=1.0d0+gamm*gamm*(2*gamm-3.0d0)
+      else
+        sscale=0d0
+      endif
+      return
+      end
+C-----------------------------------------------------------------------
+C-----------------------------------------------------------------------
+      double precision function sscagrad(r)
+      double precision r,gamm
+      include "COMMON.SPLITELE"
+      if(r.lt.r_cut-rlamb) then
+        sscagrad=0.0d0
+      else if(r.le.r_cut.and.r.ge.r_cut-rlamb) then
+        gamm=(r-(r_cut-rlamb))/rlamb
+        sscagrad=gamm*(6*gamm-6.0d0)/rlamb
+      else
+        sscagrad=0.0d0
+      endif
+      return
+      end
+C-----------------------------------------------------------------------
 
index 321e11e..f29d6f7 100644 (file)
@@ -28,6 +28,13 @@ C
       call readi(controlcard,'RESCALE',rescale_mode,2)
       call reada(controlcard,'DISTCHAINMAX',distchainmax,50.0d0)
       write (iout,*) "DISTCHAINMAX",distchainmax
+C Reading the dimensions of box in x,y,z coordinates
+      call reada(controlcard,'BOXX',boxxsize,100.0d0)
+      call reada(controlcard,'BOXY',boxysize,100.0d0)
+      call reada(controlcard,'BOXZ',boxzsize,100.0d0)
+c Cutoff range for interactions
+      call reada(controlcard,"R_CUT",r_cut,15.0d0)
+      call reada(controlcard,"LAMBDA",rlamb,0.3d0)
       call readi(controlcard,'PDBOUT',outpdb,0)
       call readi(controlcard,'MOL2OUT',outmol2,0)
       refstr=(index(controlcard,'REFSTR').gt.0)
index 7217542..752f02e 100644 (file)
@@ -3330,8 +3330,8 @@ C     &      AKSC(1,iti),AKSC(1,iti)*diff*diff
               usum=usum+uprod1
               usumsqder=usumsqder+ud(j)*uprod2
             enddo
-            write (iout,*) i,iti,vbld(i+nres),(vbldsc0(j,iti),
-     &      AKSC(j,iti),abond0(j,iti),u(j),j=1,nbi)
+c            write (iout,*) i,iti,vbld(i+nres),(vbldsc0(j,iti),
+c     &      AKSC(j,iti),abond0(j,iti),u(j),j=1,nbi)
             estr=estr+uprod/usum
             do j=1,3
              gradbx(j,i)=usumsqder/(usum*usum)*dc(j,i+nres)/vbld(i+nres)