numerical gradient delta 10e-4
[unres.git] / source / unres / src_MD-M / energy_p_new_barrier.F
index dee7c2f..eef70b5 100644 (file)
@@ -85,7 +85,7 @@ C FG slaves receive the WEIGHTS array
         time_Bcastw=time_Bcastw+MPI_Wtime()-time00
 c        call chainbuild_cart
       endif
-c      print *,'Processor',myrank,' calling etotal ipot=',ipot
+C      print *,'Processor',myrank,' calling etotal ipot=',ipot
 c      print *,'Processor',myrank,' nnt=',nnt,' nct=',nct
 #else
 c      if (modecalc.eq.12.or.modecalc.eq.14) then
@@ -98,6 +98,7 @@ c      endif
 C 
 C Compute the side-chain and electrostatic interaction energy
 C
+C      write(iout,*) "zaczynam liczyc energie"
       goto (101,102,103,104,105,106) ipot
 C Lennard-Jones potential.
   101 call elj(evdw)
@@ -117,10 +118,18 @@ C Gay-Berne-Vorobjev potential (shifted LJ, angular dependence).
       goto 107
 C Soft-sphere potential
   106 call e_softsphere(evdw)
+C      write(iout,*) "skonczylem ipoty"
+
 C
 C Calculate electrostatic (H-bonding) energy of the main chain.
 C
   107 continue
+C           write(iout,*) "skonczylem ipoty"
+cmc
+cmc Sep-06: egb takes care of dynamic ss bonds too
+cmc
+c      if (dyn_ss) call dyn_set_nss
+
 c      print *,"Processor",myrank," computed USCSC"
 #ifdef TIMING
       time01=MPI_Wtime() 
@@ -298,6 +307,7 @@ C
       energia(21)=esccor
 c      print *," Processor",myrank," calls SUM_ENERGY"
       call sum_energy(energia,.true.)
+      if (dyn_ss) call dyn_set_nss
 c      print *," Processor",myrank," left SUM_ENERGY"
 #ifdef TIMING
       time_sumene=time_sumene+MPI_Wtime()-time00
@@ -433,9 +443,9 @@ cMS$ATTRIBUTES C ::  proc_proc
 #endif
 #ifdef MPI
       include 'mpif.h'
-      double precision gradbufc(3,maxres),gradbufx(3,maxres),
-     &  glocbuf(4*maxres),gradbufc_sum(3,maxres)
 #endif
+      double precision gradbufc(3,maxres),gradbufx(3,maxres),
+     &  glocbuf(4*maxres),gradbufc_sum(3,maxres),gloc_scbuf(3,maxres)
       include 'COMMON.SETUP'
       include 'COMMON.IOUNITS'
       include 'COMMON.FFIELD'
@@ -447,6 +457,7 @@ cMS$ATTRIBUTES C ::  proc_proc
       include 'COMMON.CONTROL'
       include 'COMMON.TIME1'
       include 'COMMON.MAXGRAD'
+      include 'COMMON.SCCOR'
 #ifdef TIMING
       time01=MPI_Wtime()
 #endif
@@ -689,7 +700,6 @@ c      enddo
      &   +wturn3*gel_loc_turn3(i)
      &   +wturn6*gel_loc_turn6(i)
      &   +wel_loc*gel_loc_loc(i)
-     &   +wsccor*gsccor_loc(i)
       enddo
 #ifdef DEBUG
       write (iout,*) "gloc after adding corr"
@@ -708,6 +718,21 @@ c      enddo
         do i=1,4*nres
           glocbuf(i)=gloc(i,icg)
         enddo
+#undef DEBUG
+#ifdef DEBUG
+      write (iout,*) "gloc_sc before reduce"
+      do i=1,nres
+       do j=1,1
+        write (iout,*) i,j,gloc_sc(j,i,icg)
+       enddo
+      enddo
+#endif
+#undef DEBUG
+        do i=1,nres
+         do j=1,3
+          gloc_scbuf(j,i)=gloc_sc(j,i,icg)
+         enddo
+        enddo
         time00=MPI_Wtime()
         call MPI_Barrier(FG_COMM,IERR)
         time_barrier_g=time_barrier_g+MPI_Wtime()-time00
@@ -719,6 +744,19 @@ c      enddo
         call MPI_Reduce(glocbuf(1),gloc(1,icg),4*nres,
      &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
         time_reduce=time_reduce+MPI_Wtime()-time00
+        call MPI_Reduce(gloc_scbuf(1,1),gloc_sc(1,1,icg),3*nres,
+     &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
+        time_reduce=time_reduce+MPI_Wtime()-time00
+#undef DEBUG
+#ifdef DEBUG
+      write (iout,*) "gloc_sc after reduce"
+      do i=1,nres
+       do j=1,1
+        write (iout,*) i,j,gloc_sc(j,i,icg)
+       enddo
+      enddo
+#endif
+#undef DEBUG
 #ifdef DEBUG
       write (iout,*) "gloc after reduce"
       do i=1,4*nres
@@ -1382,6 +1420,7 @@ C
       include 'COMMON.IOUNITS'
       include 'COMMON.CALC'
       include 'COMMON.CONTROL'
+      include 'COMMON.SBRIDGE'
       logical lprn
       evdw=0.0D0
 ccccc      energy_dec=.false.
@@ -1409,6 +1448,12 @@ C Calculate SC interaction energy.
 C
         do iint=1,nint_gr(i)
           do j=istart(i,iint),iend(i,iint)
+            IF (dyn_ss_mask(i).and.dyn_ss_mask(j)) THEN
+              call dyn_ssbond_ene(i,j,evdwij)
+              evdw=evdw+evdwij
+              if (energy_dec) write (iout,'(a6,2i5,0pf7.3,a3)') 
+     &                        'evdw',i,j,evdwij,' ss'
+            ELSE
             ind=ind+1
             itypj=itype(j)
             if (itypj.eq.21) cycle
@@ -1503,6 +1548,7 @@ C Calculate the radial part of the gradient
             gg(3)=zj*fac
 C Calculate angular part of the gradient.
             call sc_grad
+            ENDIF    ! dyn_ss            
           enddo      ! j
         enddo        ! iint
       enddo          ! i
@@ -2295,6 +2341,7 @@ C
         endif
 c        if (i.gt. iatel_s+2 .and. i.lt.iatel_e+5) then
         if (i.gt. nnt+2 .and. i.lt.nct+2) then
+c        write(iout,*) (itype(i-2))
           iti = itortyp(itype(i-2))
         else
           iti=ntortyp+1
@@ -2348,7 +2395,11 @@ c        if (i.gt. iatel_s+1 .and. i.lt.iatel_e+4) then
         do k=1,2
           mu(k,i-2)=Ub2(k,i-2)+b1(k,iti1)
         enddo
-cd        write (iout,*) 'mu ',mu(:,i-2)
+cd        write (iout,*) 'mu  ',mu(:,i-2),i-2
+cd        write (iout,*) 'b1  ',b1(:,iti1),i-2
+cd        write (iout,*) 'Ub2 ',Ub2(:,i-2),i-2
+cd        write (iout,*) 'Ug  ',Ug(:,:,i-2),i-2
+cd        write (iout,*) 'b2  ',b2(:,itortyp(itype(i))),i-2
 cd        write (iout,*) 'mu1',mu1(:,i-2)
 cd        write (iout,*) 'mu2',mu2(:,i-2)
         if (wcorr4.gt.0.0d0 .or. wcorr5.gt.0.0d0 .or.wcorr6.gt.0.0d0)
@@ -2675,7 +2726,7 @@ C
       dimension ggg(3),gggp(3),gggm(3),erij(3),dcosb(3),dcosg(3),
      &          erder(3,3),uryg(3,3),urzg(3,3),vryg(3,3),vrzg(3,3)
       double precision acipa(2,2),agg(3,4),aggi(3,4),aggi1(3,4),
-     &    aggj(3,4),aggj1(3,4),a_temp(2,2),muij(4)
+     &    aggj(3,4),aggj1(3,4),a_temp(2,2),muij(4),eel_loc_ij
       common /locel/ a_temp,agg,aggi,aggi1,aggj,aggj1,a22,a23,a32,a33,
      &    dxi,dyi,dzi,dx_normi,dy_normi,dz_normi,xmedi,ymedi,zmedi,
      &    num_conti,j1,j2
@@ -2718,6 +2769,7 @@ c        call vec_and_deriv
         time01=MPI_Wtime()
 #endif
         call set_matrices
+c        write (iout,*) "after set matrices"
 #ifdef TIMING
         time_mat=time_mat+MPI_Wtime()-time01
 #endif
@@ -2754,6 +2806,7 @@ c 9/27/08 AL Split the interaction loop to ensure load balancing of turn terms
 C
 C Loop over i,i+2 and i,i+3 pairs of the peptide groups
 C
+c      write(iout,*) "przed turnem3 loop"
       do i=iturn3_start,iturn3_end
         if (itype(i).eq.21 .or. itype(i+1).eq.21 
      &  .or. itype(i+2).eq.21 .or. itype(i+3).eq.21) cycle
@@ -2846,7 +2899,7 @@ C-------------------------------------------------------------------------------
       dimension ggg(3),gggp(3),gggm(3),erij(3),dcosb(3),dcosg(3),
      &          erder(3,3),uryg(3,3),urzg(3,3),vryg(3,3),vrzg(3,3)
       double precision acipa(2,2),agg(3,4),aggi(3,4),aggi1(3,4),
-     &    aggj(3,4),aggj1(3,4),a_temp(2,2),muij(4)
+     &    aggj(3,4),aggj1(3,4),a_temp(2,2),muij(4),a22,a23,a32,a33
       common /locel/ a_temp,agg,aggi,aggi1,aggj,aggj1,a22,a23,a32,a33,
      &    dxi,dyi,dzi,dx_normi,dy_normi,dz_normi,xmedi,ymedi,zmedi,
      &    num_conti,j1,j2
@@ -3616,7 +3669,7 @@ c        write (iout,*) "eturn4 i",i," j",j," j1",j1," j2",j2
         iti1=itortyp(itype(i+1))
         iti2=itortyp(itype(i+2))
         iti3=itortyp(itype(i+3))
-c        write(iout,*) "iti1",iti1," iti2",iti2," iti3",iti3
+C        write(iout,*) i,"iti1",iti1," iti2",iti2," iti3",iti3,itype(i+3)
         call transpose2(EUg(1,1,i+1),e1t(1,1))
         call transpose2(Eug(1,1,i+2),e2t(1,1))
         call transpose2(Eug(1,1,i+3),e3t(1,1))
@@ -4021,49 +4074,56 @@ C iii and jjj point to the residues for which the distance is assigned.
           iii=ii
           jjj=jj
         endif
-cd        write (iout,*) "i",i," ii",ii," iii",iii," jj",jj," jjj",jjj
+c        write (iout,*) "i",i," ii",ii," iii",iii," jj",jj," jjj",jjj,
+c     &    dhpb(i),dhpb1(i),forcon(i)
 C 24/11/03 AL: SS bridges handled separately because of introducing a specific
 C    distance and angle dependent SS bond potential.
-        if (ii.gt.nres .and. itype(iii).eq.1 .and. itype(jjj).eq.1) then
+cmc        if (ii.gt.nres .and. itype(iii).eq.1 .and. itype(jjj).eq.1) then
+C 18/07/06 MC: Use the convention that the first nss pairs are SS bonds
+        if (.not.dyn_ss .and. i.le.nss) then
+C 15/02/13 CC dynamic SSbond - additional check
+         if (ii.gt.nres 
+     &       .and. itype(iii).eq.1 .and. itype(jjj).eq.1) then 
           call ssbond_ene(iii,jjj,eij)
           ehpb=ehpb+2*eij
+         endif
 cd          write (iout,*) "eij",eij
         else
 C Calculate the distance between the two points and its difference from the
 C target distance.
-        dd=dist(ii,jj)
-        rdis=dd-dhpb(i)
+          dd=dist(ii,jj)
+            rdis=dd-dhpb(i)
 C Get the force constant corresponding to this distance.
-        waga=forcon(i)
+            waga=forcon(i)
 C Calculate the contribution to energy.
-        ehpb=ehpb+waga*rdis*rdis
+            ehpb=ehpb+waga*rdis*rdis
 C
 C Evaluate gradient.
 C
-        fac=waga*rdis/dd
+            fac=waga*rdis/dd
 cd      print *,'i=',i,' ii=',ii,' jj=',jj,' dhpb=',dhpb(i),' dd=',dd,
 cd   &   ' waga=',waga,' fac=',fac
-        do j=1,3
-          ggg(j)=fac*(c(j,jj)-c(j,ii))
-        enddo
+            do j=1,3
+              ggg(j)=fac*(c(j,jj)-c(j,ii))
+            enddo
 cd      print '(i3,3(1pe14.5))',i,(ggg(j),j=1,3)
 C If this is a SC-SC distance, we need to calculate the contributions to the
 C Cartesian gradient in the SC vectors (ghpbx).
-        if (iii.lt.ii) then
+          if (iii.lt.ii) then
           do j=1,3
             ghpbx(j,iii)=ghpbx(j,iii)-ggg(j)
             ghpbx(j,jjj)=ghpbx(j,jjj)+ggg(j)
           enddo
-        endif
+          endif
 cgrad        do j=iii,jjj-1
 cgrad          do k=1,3
 cgrad            ghpbc(k,j)=ghpbc(k,j)+ggg(k)
 cgrad          enddo
 cgrad        enddo
-        do k=1,3
-          ghpbc(k,jjj)=ghpbc(k,jjj)+ggg(k)
-          ghpbc(k,iii)=ghpbc(k,iii)-ggg(k)
-        enddo
+          do k=1,3
+            ghpbc(k,jjj)=ghpbc(k,jjj)+ggg(k)
+            ghpbc(k,iii)=ghpbc(k,iii)-ggg(k)
+          enddo
         endif
       enddo
       ehpb=0.5D0*ehpb
@@ -4189,7 +4249,7 @@ c
      &       "estr1",i,gnmr1(vbld(i),-1.0d0,distchainmax)
         else
         diff = vbld(i)-vbldp0
-        if (energy_dec) write (iout,*) 
+        if (energy_dec) write (iout,'(a7,i5,4f7.3)') 
      &     "estr bb",i,vbld(i),vbldp0,diff,AKP*diff*diff
         estr=estr+diff*diff
         do j=1,3
@@ -4490,7 +4550,8 @@ C
       logical lprn /.false./, lprn1 /.false./
       etheta=0.0D0
       do i=ithet_start,ithet_end
-        if (itype(i-1).eq.21) cycle
+        if ((itype(i-1).eq.ntyp1).or.(itype(i-2).eq.ntyp1).or.
+     &(itype(i).eq.ntyp1)) cycle
         dethetai=0.0d0
         dephii=0.0d0
         dephii1=0.0d0
@@ -4500,7 +4561,8 @@ C
           coskt(k)=dcos(k*theti2)
           sinkt(k)=dsin(k*theti2)
         enddo
-        if (i.gt.3 .and. itype(i-2).ne.21) then
+C        if (i.gt.3) then
+        if (i.gt.3 .and. itype(max0(i-3,1)).ne.ntyp1) then
 #ifdef OSF
           phii=phi(i)
           if (phii.ne.phii) phii=150.0
           enddo
         else
           phii=0.0d0
-          ityp1=nthetyp+1
+          ityp1=ithetyp(itype(i-2))
           do k=1,nsingle
             cosph1(k)=0.0d0
             sinph1(k)=0.0d0
           enddo 
         endif
-        if (i.lt.nres .and. itype(i).ne.21) then
+        if ((i.lt.nres).and. itype(i+1).ne.ntyp1) then
 #ifdef OSF
           phii1=phi(i+1)
           if (phii1.ne.phii1) phii1=150.0
@@ -4535,7 +4597,7 @@ C
           enddo
         else
           phii1=0.0d0
-          ityp3=nthetyp+1
+          ityp3=ithetyp(itype(i))
           do k=1,nsingle
             cosph2(k)=0.0d0
             sinph2(k)=0.0d0
@@ -4645,6 +4707,8 @@ C
      &   i,theta(i)*rad2deg,phii*rad2deg,
      &   phii1*rad2deg,ethetai
         etheta=etheta+ethetai
+        if (energy_dec) write (iout,'(a6,i5,0pf7.3)')
+     &      'ebend',i,ethetai
         if (i.gt.3) gloc(i-3,icg)=gloc(i-3,icg)+wang*dephii
         if (i.lt.nres) gloc(i-2,icg)=gloc(i-2,icg)+wang*dephii1
         gloc(nphi+i-2,icg)=wang*dethetai
@@ -5086,6 +5150,8 @@ c     &   sumene4,
 c     &   dscp1,dscp2,sumene
 c        sumene = enesc(x,xx,yy,zz,cost2tab(i+1),sint2tab(i+1))
         escloc = escloc + sumene
+        if (energy_dec) write (iout,'(a6,i5,0pf7.3)')
+     &     'escloc',i,sumene
 c        write (2,*) "i",i," escloc",sumene,escloc
 #ifdef DEBUG
 C
@@ -5503,8 +5569,9 @@ c     lprn=.true.
       etors=0.0D0
       do i=iphi_start,iphi_end
         if (itype(i-2).eq.21 .or. itype(i-1).eq.21 
-     &       .or. itype(i).eq.21) cycle
-        etors_ii=0.0D0
+     &       .or. itype(i).eq.21
+     &       .or. itype(i-3).eq.ntyp1) cycle
+      etors_ii=0.0D0
         itori=itortyp(itype(i-2))
         itori1=itortyp(itype(i-1))
         phii=phi(i)
@@ -5591,14 +5658,18 @@ C 6/23/01 Compute double torsional energy
       include 'COMMON.IOUNITS'
       include 'COMMON.FFIELD'
       include 'COMMON.TORCNSTR'
+      include 'COMMON.CONTROL'
       logical lprn
 C Set lprn=.true. for debugging
       lprn=.false.
 c     lprn=.true.
       etors_d=0.0D0
+C      write(iout,*) "a tu??"
       do i=iphid_start,iphid_end
         if (itype(i-2).eq.21 .or. itype(i-1).eq.21
-     &      .or. itype(i).eq.21 .or. itype(i+1).eq.21) cycle
+     &      .or. itype(i).eq.21 .or. itype(i+1).eq.21
+     &       .or. itype(i-3).eq.ntyp1) cycle
+        etors_d_ii=0.0D0
         itori=itortyp(itype(i-2))
         itori1=itortyp(itype(i-1))
         itori2=itortyp(itype(i))
@@ -5618,6 +5689,8 @@ C Regular cosine and sine terms
           sinphi2=dsin(j*phii1)
           etors_d=etors_d+v1cij*cosphi1+v1sij*sinphi1+
      &     v2cij*cosphi2+v2sij*sinphi2
+          if (energy_dec) etors_d_ii=etors_d_ii+
+     &     v1cij*cosphi1+v1sij*sinphi1+v2cij*cosphi2+v2sij*sinphi2
           gloci1=gloci1+j*(v1sij*cosphi1-v1cij*sinphi1)
           gloci2=gloci2+j*(v2sij*cosphi2-v2cij*sinphi2)
         enddo
@@ -5633,12 +5706,17 @@ C Regular cosine and sine terms
             sinphi1m2=dsin(l*phii-(k-l)*phii1)
             etors_d=etors_d+v1cdij*cosphi1p2+v2cdij*cosphi1m2+
      &        v1sdij*sinphi1p2+v2sdij*sinphi1m2
+            if (energy_dec) etors_d_ii=etors_d_ii+
+     &        v1cdij*cosphi1p2+v2cdij*cosphi1m2+
+     &        v1sdij*sinphi1p2+v2sdij*sinphi1m2
             gloci1=gloci1+l*(v1sdij*cosphi1p2+v2sdij*cosphi1m2
      &        -v1cdij*sinphi1p2-v2cdij*sinphi1m2)
             gloci2=gloci2+(k-l)*(v1sdij*cosphi1p2-v2sdij*cosphi1m2
      &        -v1cdij*sinphi1p2+v2cdij*sinphi1m2) 
           enddo
         enddo
+          if (energy_dec) write (iout,'(a6,i5,0pf7.3)')
+     &         'etor_d',i,etors_d_ii
         gloc(i-3,icg)=gloc(i-3,icg)+wtor_d*gloci1
         gloc(i-2,icg)=gloc(i-2,icg)+wtor_d*gloci2
       enddo
@@ -5671,29 +5749,58 @@ c        amino-acid residues.
 C Set lprn=.true. for debugging
       lprn=.false.
 c      lprn=.true.
-c      write (iout,*) "EBACK_SC_COR",iphi_start,iphi_end,nterm_sccor
+c      write (iout,*) "EBACK_SC_COR",itau_start,itau_end
       esccor=0.0D0
-      do i=iphi_start,iphi_end
-        if (itype(i-2).eq.21 .or. itype(i-1).eq.21) cycle
-        esccor_ii=0.0D0
-        itori=itype(i-2)
-        itori1=itype(i-1)
+      do i=itau_start,itau_end
+        if ((itype(i-2).eq.ntyp1).or.(itype(i-1).eq.ntyp1)) cycle
+
+        isccori=isccortyp(itype(i-2))
+        isccori1=isccortyp(itype(i-1))
+c      write (iout,*) "EBACK_SC_COR",i,nterm_sccor(isccori,isccori1)
         phii=phi(i)
+        do intertyp=1,3 !intertyp
+         esccor_ii=0.0D0
+cc Added 09 May 2012 (Adasko)
+cc  Intertyp means interaction type of backbone mainchain correlation: 
+c   1 = SC...Ca...Ca...Ca
+c   2 = Ca...Ca...Ca...SC
+c   3 = SC...Ca...Ca...SCi
         gloci=0.0D0
-        do j=1,nterm_sccor
-          v1ij=v1sccor(j,itori,itori1)
-          v2ij=v2sccor(j,itori,itori1)
-          cosphi=dcos(j*phii)
-          sinphi=dsin(j*phii)
+        if (((intertyp.eq.3).and.((itype(i-2).eq.10).or.
+     &      (itype(i-1).eq.10).or.(itype(i-2).eq.ntyp1).or.
+     &      (itype(i-1).eq.ntyp1)))
+     &    .or. ((intertyp.eq.1).and.((itype(i-2).eq.10)
+     &     .or.(itype(i-2).eq.ntyp1).or.(itype(i-1).eq.ntyp1)
+     &     .or.(itype(i).eq.ntyp1)))
+     &    .or.((intertyp.eq.2).and.((itype(i-1).eq.10).or.
+     &      (itype(i-1).eq.ntyp1).or.(itype(i-2).eq.ntyp1).or.
+     &      (itype(i-3).eq.ntyp1)))) cycle
+        if ((intertyp.eq.2).and.(i.eq.4).and.(itype(1).eq.ntyp1)) cycle
+        if ((intertyp.eq.1).and.(i.eq.nres).and.(itype(nres).eq.ntyp1))
+     & cycle
+       do j=1,nterm_sccor(isccori,isccori1)
+          v1ij=v1sccor(j,intertyp,isccori,isccori1)
+          v2ij=v2sccor(j,intertyp,isccori,isccori1)
+          cosphi=dcos(j*tauangle(intertyp,i))
+          sinphi=dsin(j*tauangle(intertyp,i))
+          if (energy_dec) esccor_ii=esccor_ii+v1ij*cosphi+v2ij*sinphi
           esccor=esccor+v1ij*cosphi+v2ij*sinphi
           gloci=gloci+j*(v2ij*cosphi-v1ij*sinphi)
         enddo
+          if (energy_dec) write (iout,'(a6,i5,i2,0pf7.3)')
+     &         'esccor',i,intertyp,esccor_ii
+cd       write (iout,*) "tau ",i,intertyp,tauangle(intertyp,i)*RAD2DEG
+c      write (iout,*) "EBACK_SC_COR",i,v1ij*cosphi+v2ij*sinphi,intertyp
+        gloc_sc(intertyp,i-3,icg)=gloc_sc(intertyp,i-3,icg)+wsccor*gloci
         if (lprn)
      &  write (iout,'(2(a3,2x,i3,2x),2i3,6f8.3/26x,6f8.3/)')
-     &  restyp(itype(i-2)),i-2,restyp(itype(i-1)),i-1,itori,itori1,
-     &  (v1sccor(j,itori,itori1),j=1,6),(v2sccor(j,itori,itori1),j=1,6)
+     &  restyp(itype(i-2)),i-2,restyp(itype(i-1)),i-1,isccori,isccori1,
+     &  (v1sccor(j,intertyp,isccori,isccori1),j=1,6)
+     & ,(v2sccor(j,intertyp,isccori,isccori1),j=1,6)
         gsccor_loc(i-3)=gsccor_loc(i-3)+gloci
+       enddo !intertyp
       enddo
+
       return
       end
 c----------------------------------------------------------------------------
@@ -7907,12 +8014,12 @@ C                                                                              C
 C          o             o                                                     C
 C     \   /l\           /j\   /                                                C
 C      \ /   \         /   \ /                                                 C
-C       o| o |         | o |o                                                  C
+C       o| o |         | o |o                                                  C                
 C     \ j|/k\|      \  |/k\|l                                                  C
 C      \ /   \       \ /   \                                                   C
 C       o             o                                                        C
-C       i             i                                                        C
-C                                                                              C
+C       i             i                                                        C 
+C                                                                              C           
 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
 cd      write (2,*) 'eello6_graph2: i,',i,' j',j,' k',k,' l',l
 C AL 7/4/01 s1 would occur in the sixth-order moment, 
@@ -8083,10 +8190,10 @@ c----------------------------------------------------------------------------
       double precision vv(2),pizda(2,2),auxmat(2,2),auxvec(2)
       logical swap
 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
-C                                                                              C
+C                                                                              C 
 C      Parallel       Antiparallel                                             C
 C                                                                              C
-C          o             o                                                     C
+C          o             o                                                     C 
 C         /l\   /   \   /j\                                                    C 
 C        /   \ /     \ /   \                                                   C
 C       /| o |o       o| o |\                                                  C
@@ -8200,7 +8307,7 @@ c----------------------------------------------------------------------------
      & auxvec1(2),auxmat1(2,2)
       logical swap
 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
-C                                                                              C
+C                                                                              C                       
 C      Parallel       Antiparallel                                             C
 C                                                                              C
 C          o             o                                                     C
@@ -8208,10 +8315,10 @@ C         /l\   /   \   /j\                                                    C
 C        /   \ /     \ /   \                                                   C
 C       /| o |o       o| o |\                                                  C
 C     \ j|/k\|      \  |/k\|l                                                  C
-C      \ /   \       \ /   \                                                   C
+C      \ /   \       \ /   \                                                   C 
 C       o     \       o     \                                                  C
 C       i             i                                                        C
-C                                                                              C
+C                                                                              C 
 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
 C
 C 4/7/01 AL Component s1 was removed, because it pertains to the respective