working nano
[unres4.git] / source / unres / energy.f90
index 77c46dd..c166edf 100644 (file)
         gshieldc_ec,gshieldc_loc_ec,gshieldx_t3, &
         gshieldc_t3,gshieldc_loc_t3,gshieldx_t4,gshieldc_t4, &
         gshieldc_loc_t4,gshieldx_ll,gshieldc_ll,gshieldc_loc_ll,&
-        grad_shield !(3,maxres)
+        grad_shield,gg_tube,gg_tube_sc !(3,maxres)
 !      real(kind=8),dimension(:,:),allocatable :: gloc,gloc_x !(maxvar,2)
       real(kind=8),dimension(:,:),allocatable :: gel_loc,gel_loc_long,&
         gcorr3_turn,gcorr4_turn,gcorr6_turn,gradb,gradbx !(3,maxres)
       integer :: n_corr,n_corr1,ierror
       real(kind=8) :: etors,edihcnstr,etors_d,esccor,ehpb
       real(kind=8) :: evdw,evdw1,evdw2,evdw2_14,escloc,ees,eel_loc
-      real(kind=8) :: eello_turn3,eello_turn4,estr,ebe,eliptran
+      real(kind=8) :: eello_turn3,eello_turn4,estr,ebe,eliptran,etube
       real(kind=8) :: ecorr,ecorr5,ecorr6,eturn6
 
 #ifdef MPI      
       else
        eliptran=0.0d0
       endif
+      if (tubemode.eq.1) then
+       call calctube(etube)
+      else if (tubemode.eq.2) then
+       call calctube2(etube)
+      elseif (tubemode.eq.3) then
+       call calcnano(etube)
+      else
+       etube=0.0d0
+      endif
 
 #ifdef TIMING
       time_enecalc=time_enecalc+MPI_Wtime()-time00
       energia(20)=Uconst+Uconst_back
       energia(21)=esccor
       energia(22)=eliptran
+      energia(25)=etube
 !    Here are the energies showed per procesor if the are more processors 
 !    per molecule then we sum it up in sum_energy subroutine 
 !      print *," Processor",myrank," calls SUM_ENERGY"
       real(kind=8) :: evdw,evdw2,evdw2_14,ees,evdw1,ecorr,ecorr5,ecorr6
       real(kind=8) :: eel_loc,eello_turn3,eello_turn4,eturn6,ebe,escloc
       real(kind=8) :: etors,etors_d,ehpb,edihcnstr,estr,esccor,etot,   &
-        eliptran
+        eliptran,etube
       integer :: i
 #ifdef MPI
       integer :: ierr
       Uconst=energia(20)
       esccor=energia(21)
       eliptran=energia(22)
+      etube=energia(25)
 #ifdef SPLITELE
       etot=wsc*evdw+wscp*evdw2+welec*ees+wvdwpp*evdw1 &
        +wang*ebe+wtor*etors+wscloc*escloc &
        +wstrain*ehpb+wcorr*ecorr+wcorr5*ecorr5 &
        +wcorr6*ecorr6+wturn4*eello_turn4+wturn3*eello_turn3 &
        +wturn6*eturn6+wel_loc*eel_loc+edihcnstr+wtor_d*etors_d &
-       +wbond*estr+Uconst+wsccor*esccor+wliptran*eliptran
+       +wbond*estr+Uconst+wsccor*esccor+wliptran*eliptran+wtube*etube
 #else
       etot=wsc*evdw+wscp*evdw2+welec*(ees+evdw1) &
        +wang*ebe+wtor*etors+wscloc*escloc &
        +wstrain*ehpb+wcorr*ecorr+wcorr5*ecorr5 &
        +wcorr6*ecorr6+wturn4*eello_turn4+wturn3*eello_turn3 &
        +wturn6*eturn6+wel_loc*eel_loc+edihcnstr+wtor_d*etors_d &
-       +wbond*estr+Uconst+wsccor*esccor+wliptran*eliptran
+       +wbond*estr+Uconst+wsccor*esccor+wliptran*eliptran+wtube*etube
 #endif
       energia(0)=etot
 ! detecting NaNQ
 !el local variables
       real(kind=8) :: etot,evdw,evdw2,ees,evdw1,ecorr,ecorr5,ecorr6,eel_loc
       real(kind=8) :: eello_turn6,eello_turn3,eello_turn4,ebe,escloc
-      real(kind=8) :: etors,etors_d,ehpb,edihcnstr,estr,Uconst,esccor,eliptran
+      real(kind=8) :: etors,etors_d,ehpb,edihcnstr,estr,Uconst,esccor,eliptran,&
+       etube
 
       etot=energia(0)
       evdw=energia(1)
       Uconst=energia(20)
       esccor=energia(21)
       eliptran=energia(22)
-
+      etube=energia(25)
 #ifdef SPLITELE
       write (iout,10) evdw,wsc,evdw2,wscp,ees,welec,evdw1,wvdwpp,&
         estr,wbond,ebe,wang,&
         ecorr5,wcorr5,ecorr6,wcorr6,eel_loc,wel_loc,eello_turn3,wturn3,&
         eello_turn4,wturn4,eello_turn6,wturn6,esccor,wsccor,&
         edihcnstr,ebr*nss,&
-        Uconst,eliptran,wliptran,etot
+        Uconst,eliptran,wliptran,etube,wtube,etot
    10 format (/'Virtual-chain energies:'// &
        'EVDW=  ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-SC)'/ &
        'EVDW2= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-p)'/ &
        'ESS=   ',1pE16.6,' (disulfide-bridge intrinsic energy)'/ &
        'UCONST= ',1pE16.6,' (Constraint energy)'/ &
        'ELT=',1pE16.6, ' WEIGHT=',1pD16.6,' (Lipid transfer energy)'/&
+       'ETUBE=',1pE16.6, ' WEIGHT=',1pD16.6,' (cylindrical energy)'/ &
        'ETOT=  ',1pE16.6,' (total)')
 #else
       write (iout,10) evdw,wsc,evdw2,wscp,ees,welec,&
         ecorr,wcorr,&
         ecorr5,wcorr5,ecorr6,wcorr6,eel_loc,wel_loc,eello_turn3,wturn3,&
         eello_turn4,wturn4,eello_turn6,wturn6,esccor,wsccor,edihcnstr,&
-        ebr*nss,Uconst,eliptran,wliptran,etot
+        ebr*nss,Uconst,eliptran,wliptran,etube,wtube,etot
    10 format (/'Virtual-chain energies:'// &
        'EVDW=  ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-SC)'/ &
        'EVDW2= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-p)'/ &
        'ESS=   ',1pE16.6,' (disulfide-bridge intrinsic energy)'/ &
        'UCONST=',1pE16.6,' (Constraint energy)'/ &
        'ELT=',1pE16.6, ' WEIGHT=',1pD16.6,' (Lipid transfer energy)'/ &
+       'ETUBE=',1pE16.6, ' WEIGHT=',1pD16.6,' (cylindrical energy)'/ &
        'ETOT=  ',1pE16.6,' (total)')
 #endif
       return
                      +wcorr*gshieldc_ec(j,i) &
                      +wturn3*gshieldc_t3(j,i)&
                      +wturn4*gshieldc_t4(j,i)&
-                     +wel_loc*gshieldc_ll(j,i) 
+                     +wel_loc*gshieldc_ll(j,i)&
+                     +wtube*gg_tube(j,i)
 
 
         enddo
                      +welec*gshieldc(j,i)&
                      +wcorr*gshieldc_ec(j,i) &
                      +wturn4*gshieldc_t4(j,i) &
-                     +wel_loc*gshieldc_ll(j,i)
+                     +wel_loc*gshieldc_ll(j,i)&
+                     +wtube*gg_tube(j,i)
+
 
 
         enddo
                      +wturn4*gshieldc_t4(j,i) &
                      +wturn4*gshieldc_loc_t4(j,i) &
                      +wel_loc*gshieldc_ll(j,i) &
-                     +wel_loc*gshieldc_loc_ll(j,i) 
+                     +wel_loc*gshieldc_loc_ll(j,i) &
+                     +wtube*gg_tube(j,i)
+
 
 #else
           gradc(j,i,icg)=gradbufc(j,i)+welec*gelc(j,i)+ &
                      +wturn4*gshieldc_t4(j,i) &
                      +wturn4*gshieldc_loc_t4(j,i) &
                      +wel_loc*gshieldc_ll(j,i) &
-                     +wel_loc*gshieldc_loc_ll(j,i) 
+                     +wel_loc*gshieldc_loc_ll(j,i) &
+                     +wtube*gg_tube(j,i)
+
 
 
 #endif
                        +wcorr*gshieldx_ec(j,i)  &
                        +wturn3*gshieldx_t3(j,i) &
                        +wturn4*gshieldx_t4(j,i) &
-                       +wel_loc*gshieldx_ll(j,i)
+                       +wel_loc*gshieldx_ll(j,i)&
+                       +wtube*gg_tube_sc(j,i)
+
 
         enddo
       enddo 
@@ -15860,7 +15884,8 @@ write(iout,*) 'Calling CHECK_ECARTIN else.'
           gshieldx_ll(j,i)=0.0d0
           gshieldc_ll(j,i)=0.0d0
           gshieldc_loc_ll(j,i)=0.0d0
-
+          gg_tube(j,i)=0.0d0
+          gg_tube_sc(j,i)=0.0d0
           do intertyp=1,3
            gloc_sc(intertyp,i,icg)=0.0d0
           enddo
@@ -17904,6 +17929,634 @@ write(iout,*) 'Calling CHECK_ECARTIN else.'
        enddo
        return
        end  subroutine Eliptransfer
+!----------------------------------NANO FUNCTIONS
+!C-----------------------------------------------------------------------
+!C-----------------------------------------------------------
+!C This subroutine is to mimic the histone like structure but as well can be
+!C utilizet to nanostructures (infinit) small modification has to be used to 
+!C make it finite (z gradient at the ends has to be changes as well as the x,y
+!C gradient has to be modified at the ends 
+!C The energy function is Kihara potential 
+!C E=4esp*((sigma/(r-r0))^12 - (sigma/(r-r0))^6)
+!C 4eps is depth of well sigma is r_minimum r is distance from center of tube 
+!C and r0 is the excluded size of nanotube (can be set to 0 if we want just a 
+!C simple Kihara potential
+      subroutine calctube(Etube)
+      real(kind=8) :: vectube(3),enetube(nres*2)
+      real(kind=8) :: Etube,xtemp,xminact,yminact,& 
+       ytemp,xmin,ymin,tub_r,rdiff,rdiff6,fac,positi, &
+       sc_aa_tube,sc_bb_tube
+      integer :: i,j,iti
+      Etube=0.0d0
+      do i=itube_start,itube_end
+        enetube(i)=0.0d0
+        enetube(i+nres)=0.0d0
+      enddo
+!C first we calculate the distance from tube center
+!C for UNRES
+       do i=itube_start,itube_end
+!C lets ommit dummy atoms for now
+       if ((itype(i).eq.ntyp1).or.(itype(i+1).eq.ntyp1)) cycle
+!C now calculate distance from center of tube and direction vectors
+      xmin=boxxsize
+      ymin=boxysize
+! Find minimum distance in periodic box
+        do j=-1,1
+         vectube(1)=mod((c(1,i)+c(1,i+1))/2.0d0,boxxsize)
+         vectube(1)=vectube(1)+boxxsize*j
+         vectube(2)=mod((c(2,i)+c(2,i+1))/2.0d0,boxysize)
+         vectube(2)=vectube(2)+boxysize*j
+         xminact=abs(vectube(1)-tubecenter(1))
+         yminact=abs(vectube(2)-tubecenter(2))
+           if (xmin.gt.xminact) then
+            xmin=xminact
+            xtemp=vectube(1)
+           endif
+           if (ymin.gt.yminact) then
+             ymin=yminact
+             ytemp=vectube(2)
+            endif
+         enddo
+      vectube(1)=xtemp
+      vectube(2)=ytemp
+      vectube(1)=vectube(1)-tubecenter(1)
+      vectube(2)=vectube(2)-tubecenter(2)
+
+!C      print *,"x",(c(1,i)+c(1,i+1))/2.0d0,tubecenter(1)
+!C      print *,"y",(c(2,i)+c(2,i+1))/2.0d0,tubecenter(2)
+
+!C as the tube is infinity we do not calculate the Z-vector use of Z
+!C as chosen axis
+      vectube(3)=0.0d0
+!C now calculte the distance
+       tub_r=dsqrt(vectube(1)**2+vectube(2)**2+vectube(3)**2)
+!C now normalize vector
+      vectube(1)=vectube(1)/tub_r
+      vectube(2)=vectube(2)/tub_r
+!C calculte rdiffrence between r and r0
+      rdiff=tub_r-tubeR0
+!C and its 6 power
+      rdiff6=rdiff**6.0d0
+!C for vectorization reasons we will sumup at the end to avoid depenence of previous
+       enetube(i)=pep_aa_tube/rdiff6**2.0d0+pep_bb_tube/rdiff6
+!C       write(iout,*) "TU13",i,rdiff6,enetube(i)
+!C       print *,rdiff,rdiff6,pep_aa_tube
+!C pep_aa_tube and pep_bb_tube are precomputed values A=4eps*sigma^12 B=4eps*sigma^6
+!C now we calculate gradient
+       fac=(-12.0d0*pep_aa_tube/rdiff6- &
+            6.0d0*pep_bb_tube)/rdiff6/rdiff
+!C       write(iout,'(a5,i4,f12.1,3f12.5)') "TU13",i,rdiff6,enetube(i),
+!C     &rdiff,fac
+!C now direction of gg_tube vector
+        do j=1,3
+        gg_tube(j,i-1)=gg_tube(j,i-1)+vectube(j)*fac/2.0d0
+        gg_tube(j,i)=gg_tube(j,i)+vectube(j)*fac/2.0d0
+        enddo
+        enddo
+!C basically thats all code now we split for side-chains (REMEMBER to sum up at the END)
+!C        print *,gg_tube(1,0),"TU"
+
+
+       do i=itube_start,itube_end
+!C Lets not jump over memory as we use many times iti
+         iti=itype(i)
+!C lets ommit dummy atoms for now
+         if ((iti.eq.ntyp1)  &
+!C in UNRES uncomment the line below as GLY has no side-chain...
+!C      .or.(iti.eq.10)
+        ) cycle
+      xmin=boxxsize
+      ymin=boxysize
+        do j=-1,1
+         vectube(1)=mod((c(1,i+nres)),boxxsize)
+         vectube(1)=vectube(1)+boxxsize*j
+         vectube(2)=mod((c(2,i+nres)),boxysize)
+         vectube(2)=vectube(2)+boxysize*j
+
+         xminact=abs(vectube(1)-tubecenter(1))
+         yminact=abs(vectube(2)-tubecenter(2))
+           if (xmin.gt.xminact) then
+            xmin=xminact
+            xtemp=vectube(1)
+           endif
+           if (ymin.gt.yminact) then
+             ymin=yminact
+             ytemp=vectube(2)
+            endif
+         enddo
+      vectube(1)=xtemp
+      vectube(2)=ytemp
+!C          write(iout,*), "tututu", vectube(1),tubecenter(1),vectube(2),
+!C     &     tubecenter(2)
+      vectube(1)=vectube(1)-tubecenter(1)
+      vectube(2)=vectube(2)-tubecenter(2)
+
+!C as the tube is infinity we do not calculate the Z-vector use of Z
+!C as chosen axis
+      vectube(3)=0.0d0
+!C now calculte the distance
+       tub_r=dsqrt(vectube(1)**2+vectube(2)**2+vectube(3)**2)
+!C now normalize vector
+      vectube(1)=vectube(1)/tub_r
+      vectube(2)=vectube(2)/tub_r
+
+!C calculte rdiffrence between r and r0
+      rdiff=tub_r-tubeR0
+!C and its 6 power
+      rdiff6=rdiff**6.0d0
+!C for vectorization reasons we will sumup at the end to avoid depenence of previous
+       sc_aa_tube=sc_aa_tube_par(iti)
+       sc_bb_tube=sc_bb_tube_par(iti)
+       enetube(i+nres)=sc_aa_tube/rdiff6**2.0d0+sc_bb_tube/rdiff6
+       fac=-12.0d0*sc_aa_tube/rdiff6**2.0d0/rdiff-  &
+             6.0d0*sc_bb_tube/rdiff6/rdiff
+!C now direction of gg_tube vector
+         do j=1,3
+          gg_tube_SC(j,i)=gg_tube_SC(j,i)+vectube(j)*fac
+          gg_tube(j,i-1)=gg_tube(j,i-1)+vectube(j)*fac
+         enddo
+        enddo
+        do i=itube_start,itube_end
+          Etube=Etube+enetube(i)+enetube(i+nres)
+        enddo
+!C        print *,"ETUBE", etube
+        return
+        end subroutine calctube
+!C TO DO 1) add to total energy
+!C       2) add to gradient summation
+!C       3) add reading parameters (AND of course oppening of PARAM file)
+!C       4) add reading the center of tube
+!C       5) add COMMONs
+!C       6) add to zerograd
+!C       7) allocate matrices
+
+
+!C-----------------------------------------------------------------------
+!C-----------------------------------------------------------
+!C This subroutine is to mimic the histone like structure but as well can be
+!C utilizet to nanostructures (infinit) small modification has to be used to 
+!C make it finite (z gradient at the ends has to be changes as well as the x,y
+!C gradient has to be modified at the ends 
+!C The energy function is Kihara potential 
+!C E=4esp*((sigma/(r-r0))^12 - (sigma/(r-r0))^6)
+!C 4eps is depth of well sigma is r_minimum r is distance from center of tube 
+!C and r0 is the excluded size of nanotube (can be set to 0 if we want just a 
+!C simple Kihara potential
+      subroutine calctube2(Etube)
+      real(kind=8) :: vectube(3),enetube(nres*2)
+      real(kind=8) :: Etube,xtemp,xminact,yminact,&
+       ytemp,xmin,ymin,tub_r,rdiff,rdiff6,fac,positi,fracinbuf,&
+       sstube,ssgradtube,sc_aa_tube,sc_bb_tube
+      integer:: i,j,iti
+      Etube=0.0d0
+      do i=itube_start,itube_end
+        enetube(i)=0.0d0
+        enetube(i+nres)=0.0d0
+      enddo
+!C first we calculate the distance from tube center
+!C first sugare-phosphate group for NARES this would be peptide group 
+!C for UNRES
+       do i=itube_start,itube_end
+!C lets ommit dummy atoms for now
+
+       if ((itype(i).eq.ntyp1).or.(itype(i+1).eq.ntyp1)) cycle
+!C now calculate distance from center of tube and direction vectors
+!C      vectube(1)=mod((c(1,i)+c(1,i+1))/2.0d0,boxxsize)
+!C          if (vectube(1).lt.0) vectube(1)=vectube(1)+boxxsize
+!C      vectube(2)=mod((c(2,i)+c(2,i+1))/2.0d0,boxysize)
+!C          if (vectube(2).lt.0) vectube(2)=vectube(2)+boxysize
+      xmin=boxxsize
+      ymin=boxysize
+        do j=-1,1
+         vectube(1)=mod((c(1,i)+c(1,i+1))/2.0d0,boxxsize)
+         vectube(1)=vectube(1)+boxxsize*j
+         vectube(2)=mod((c(2,i)+c(2,i+1))/2.0d0,boxysize)
+         vectube(2)=vectube(2)+boxysize*j
+
+         xminact=abs(vectube(1)-tubecenter(1))
+         yminact=abs(vectube(2)-tubecenter(2))
+           if (xmin.gt.xminact) then
+            xmin=xminact
+            xtemp=vectube(1)
+           endif
+           if (ymin.gt.yminact) then
+             ymin=yminact
+             ytemp=vectube(2)
+            endif
+         enddo
+      vectube(1)=xtemp
+      vectube(2)=ytemp
+      vectube(1)=vectube(1)-tubecenter(1)
+      vectube(2)=vectube(2)-tubecenter(2)
+
+!C      print *,"x",(c(1,i)+c(1,i+1))/2.0d0,tubecenter(1)
+!C      print *,"y",(c(2,i)+c(2,i+1))/2.0d0,tubecenter(2)
+
+!C as the tube is infinity we do not calculate the Z-vector use of Z
+!C as chosen axis
+      vectube(3)=0.0d0
+!C now calculte the distance
+       tub_r=dsqrt(vectube(1)**2+vectube(2)**2+vectube(3)**2)
+!C now normalize vector
+      vectube(1)=vectube(1)/tub_r
+      vectube(2)=vectube(2)/tub_r
+!C calculte rdiffrence between r and r0
+      rdiff=tub_r-tubeR0
+!C and its 6 power
+      rdiff6=rdiff**6.0d0
+!C THIS FRAGMENT MAKES TUBE FINITE
+        positi=mod((c(3,i)+c(3,i+1))/2.0d0,boxzsize)
+        if (positi.le.0) positi=positi+boxzsize
+!C       print *,mod(c(3,i+nres),boxzsize),bordlipbot,bordliptop
+!c for each residue check if it is in lipid or lipid water border area
+!C       respos=mod(c(3,i+nres),boxzsize)
+!C       print *,positi,bordtubebot,buftubebot,bordtubetop
+       if ((positi.gt.bordtubebot)  &
+        .and.(positi.lt.bordtubetop)) then
+!C the energy transfer exist
+        if (positi.lt.buftubebot) then
+         fracinbuf=1.0d0-  &
+           ((positi-bordtubebot)/tubebufthick)
+!C lipbufthick is thickenes of lipid buffore
+         sstube=sscalelip(fracinbuf)
+         ssgradtube=-sscagradlip(fracinbuf)/tubebufthick
+!C         print *,ssgradtube, sstube,tubetranene(itype(i))
+         enetube(i)=enetube(i)+sstube*tubetranenepep
+!C         gg_tube_SC(3,i)=gg_tube_SC(3,i)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C         gg_tube(3,i-1)= gg_tube(3,i-1)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C         print *,"doing sccale for lower part"
+        elseif (positi.gt.buftubetop) then
+         fracinbuf=1.0d0-  &
+        ((bordtubetop-positi)/tubebufthick)
+         sstube=sscalelip(fracinbuf)
+         ssgradtube=sscagradlip(fracinbuf)/tubebufthick
+         enetube(i)=enetube(i)+sstube*tubetranenepep
+!C         gg_tube_SC(3,i)=gg_tube_SC(3,i)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C         gg_tube(3,i-1)= gg_tube(3,i-1)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C          print *, "doing sscalefor top part",sslip,fracinbuf
+        else
+         sstube=1.0d0
+         ssgradtube=0.0d0
+         enetube(i)=enetube(i)+sstube*tubetranenepep
+!C         print *,"I am in true lipid"
+        endif
+        else
+!C          sstube=0.0d0
+!C          ssgradtube=0.0d0
+        cycle
+        endif ! if in lipid or buffor
+
+!C for vectorization reasons we will sumup at the end to avoid depenence of previous
+       enetube(i)=enetube(i)+sstube* &
+        (pep_aa_tube/rdiff6**2.0d0+pep_bb_tube/rdiff6)
+!C       write(iout,*) "TU13",i,rdiff6,enetube(i)
+!C       print *,rdiff,rdiff6,pep_aa_tube
+!C pep_aa_tube and pep_bb_tube are precomputed values A=4eps*sigma^12 B=4eps*sigma^6
+!C now we calculate gradient
+       fac=(-12.0d0*pep_aa_tube/rdiff6-  &
+             6.0d0*pep_bb_tube)/rdiff6/rdiff*sstube
+!C       write(iout,'(a5,i4,f12.1,3f12.5)') "TU13",i,rdiff6,enetube(i),
+!C     &rdiff,fac
+
+!C now direction of gg_tube vector
+       do j=1,3
+        gg_tube(j,i-1)=gg_tube(j,i-1)+vectube(j)*fac/2.0d0
+        gg_tube(j,i)=gg_tube(j,i)+vectube(j)*fac/2.0d0
+        enddo
+         gg_tube(3,i)=gg_tube(3,i)  &
+       +ssgradtube*enetube(i)/sstube/2.0d0
+         gg_tube(3,i-1)= gg_tube(3,i-1)  &
+       +ssgradtube*enetube(i)/sstube/2.0d0
+
+        enddo
+!C basically thats all code now we split for side-chains (REMEMBER to sum up at the END)
+!C        print *,gg_tube(1,0),"TU"
+        do i=itube_start,itube_end
+!C Lets not jump over memory as we use many times iti
+         iti=itype(i)
+!C lets ommit dummy atoms for now
+         if ((iti.eq.ntyp1) &
+!!C in UNRES uncomment the line below as GLY has no side-chain...
+           .or.(iti.eq.10) &
+          ) cycle
+          vectube(1)=c(1,i+nres)
+          vectube(1)=mod(vectube(1),boxxsize)
+          if (vectube(1).lt.0) vectube(1)=vectube(1)+boxxsize
+          vectube(2)=c(2,i+nres)
+          vectube(2)=mod(vectube(2),boxysize)
+          if (vectube(2).lt.0) vectube(2)=vectube(2)+boxysize
+
+      vectube(1)=vectube(1)-tubecenter(1)
+      vectube(2)=vectube(2)-tubecenter(2)
+!C THIS FRAGMENT MAKES TUBE FINITE
+        positi=(mod(c(3,i+nres),boxzsize))
+        if (positi.le.0) positi=positi+boxzsize
+!C       print *,mod(c(3,i+nres),boxzsize),bordlipbot,bordliptop
+!c for each residue check if it is in lipid or lipid water border area
+!C       respos=mod(c(3,i+nres),boxzsize)
+!C       print *,positi,bordtubebot,buftubebot,bordtubetop
+
+       if ((positi.gt.bordtubebot)  &
+        .and.(positi.lt.bordtubetop)) then
+!C the energy transfer exist
+        if (positi.lt.buftubebot) then
+         fracinbuf=1.0d0- &
+            ((positi-bordtubebot)/tubebufthick)
+!C lipbufthick is thickenes of lipid buffore
+         sstube=sscalelip(fracinbuf)
+         ssgradtube=-sscagradlip(fracinbuf)/tubebufthick
+!C         print *,ssgradtube, sstube,tubetranene(itype(i))
+         enetube(i+nres)=enetube(i+nres)+sstube*tubetranene(itype(i))
+!C         gg_tube_SC(3,i)=gg_tube_SC(3,i)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C         gg_tube(3,i-1)= gg_tube(3,i-1)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C         print *,"doing sccale for lower part"
+        elseif (positi.gt.buftubetop) then
+         fracinbuf=1.0d0- &
+        ((bordtubetop-positi)/tubebufthick)
+
+         sstube=sscalelip(fracinbuf)
+         ssgradtube=sscagradlip(fracinbuf)/tubebufthick
+         enetube(i+nres)=enetube(i+nres)+sstube*tubetranene(itype(i))
+!C         gg_tube_SC(3,i)=gg_tube_SC(3,i)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C         gg_tube(3,i-1)= gg_tube(3,i-1)
+!C     &+ssgradtube*tubetranene(itype(i))
+!C          print *, "doing sscalefor top part",sslip,fracinbuf
+        else
+         sstube=1.0d0
+         ssgradtube=0.0d0
+         enetube(i+nres)=enetube(i+nres)+sstube*tubetranene(itype(i))
+!C         print *,"I am in true lipid"
+        endif
+        else
+!C          sstube=0.0d0
+!C          ssgradtube=0.0d0
+        cycle
+        endif ! if in lipid or buffor
+!CEND OF FINITE FRAGMENT
+!C as the tube is infinity we do not calculate the Z-vector use of Z
+!C as chosen axis
+      vectube(3)=0.0d0
+!C now calculte the distance
+       tub_r=dsqrt(vectube(1)**2+vectube(2)**2+vectube(3)**2)
+!C now normalize vector
+      vectube(1)=vectube(1)/tub_r
+      vectube(2)=vectube(2)/tub_r
+!C calculte rdiffrence between r and r0
+      rdiff=tub_r-tubeR0
+!C and its 6 power
+      rdiff6=rdiff**6.0d0
+!C for vectorization reasons we will sumup at the end to avoid depenence of previous
+       sc_aa_tube=sc_aa_tube_par(iti)
+       sc_bb_tube=sc_bb_tube_par(iti)
+       enetube(i+nres)=(sc_aa_tube/rdiff6**2.0d0+sc_bb_tube/rdiff6)&
+                       *sstube+enetube(i+nres)
+!C pep_aa_tube and pep_bb_tube are precomputed values A=4eps*sigma^12 B=4eps*sigma^6
+!C now we calculate gradient
+       fac=(-12.0d0*sc_aa_tube/rdiff6**2.0d0/rdiff-&
+            6.0d0*sc_bb_tube/rdiff6/rdiff)*sstube
+!C now direction of gg_tube vector
+         do j=1,3
+          gg_tube_SC(j,i)=gg_tube_SC(j,i)+vectube(j)*fac
+          gg_tube(j,i-1)=gg_tube(j,i-1)+vectube(j)*fac
+         enddo
+         gg_tube_SC(3,i)=gg_tube_SC(3,i) &
+       +ssgradtube*enetube(i+nres)/sstube
+         gg_tube(3,i-1)= gg_tube(3,i-1) &
+       +ssgradtube*enetube(i+nres)/sstube
+
+        enddo
+        do i=itube_start,itube_end
+          Etube=Etube+enetube(i)+enetube(i+nres)
+        enddo
+!C        print *,"ETUBE", etube
+        return
+        end subroutine calctube2
+!=====================================================================================================================================
+      subroutine calcnano(Etube)
+      real(kind=8) :: vectube(3),enetube(nres*2), &
+      enecavtube(nres*2)
+      real(kind=8) :: Etube,xtemp,xminact,yminact,&
+       ytemp,xmin,ymin,tub_r,rdiff,rdiff6,fac,denominator,faccav,&
+       sc_aa_tube,sc_bb_tube,zmin,ztemp,zminact
+       integer:: i,j,iti
+
+      Etube=0.0d0
+      print *,itube_start,itube_end,"poczatek"
+      do i=itube_start,itube_end
+        enetube(i)=0.0d0
+        enetube(i+nres)=0.0d0
+      enddo
+!C first we calculate the distance from tube center
+!C first sugare-phosphate group for NARES this would be peptide group 
+!C for UNRES
+       do i=itube_start,itube_end
+!C lets ommit dummy atoms for now
+       if ((itype(i).eq.ntyp1).or.(itype(i+1).eq.ntyp1)) cycle
+!C now calculate distance from center of tube and direction vectors
+      xmin=boxxsize
+      ymin=boxysize
+      zmin=boxzsize
+
+        do j=-1,1
+         vectube(1)=dmod((c(1,i)+c(1,i+1))/2.0d0,boxxsize)
+         vectube(1)=vectube(1)+boxxsize*j
+         vectube(2)=dmod((c(2,i)+c(2,i+1))/2.0d0,boxysize)
+         vectube(2)=vectube(2)+boxysize*j
+         vectube(3)=dmod((c(3,i)+c(3,i+1))/2.0d0,boxzsize)
+         vectube(3)=vectube(3)+boxzsize*j
+
+
+         xminact=dabs(vectube(1)-tubecenter(1))
+         yminact=dabs(vectube(2)-tubecenter(2))
+         zminact=dabs(vectube(3)-tubecenter(3))
+
+           if (xmin.gt.xminact) then
+            xmin=xminact
+            xtemp=vectube(1)
+           endif
+           if (ymin.gt.yminact) then
+             ymin=yminact
+             ytemp=vectube(2)
+            endif
+           if (zmin.gt.zminact) then
+             zmin=zminact
+             ztemp=vectube(3)
+            endif
+         enddo
+      vectube(1)=xtemp
+      vectube(2)=ytemp
+      vectube(3)=ztemp
+
+      vectube(1)=vectube(1)-tubecenter(1)
+      vectube(2)=vectube(2)-tubecenter(2)
+      vectube(3)=vectube(3)-tubecenter(3)
+
+!C      print *,"x",(c(1,i)+c(1,i+1))/2.0d0,tubecenter(1)
+!C      print *,"y",(c(2,i)+c(2,i+1))/2.0d0,tubecenter(2)
+!C as the tube is infinity we do not calculate the Z-vector use of Z
+!C as chosen axis
+!C      vectube(3)=0.0d0
+!C now calculte the distance
+       tub_r=dsqrt(vectube(1)**2+vectube(2)**2+vectube(3)**2)
+!C now normalize vector
+      vectube(1)=vectube(1)/tub_r
+      vectube(2)=vectube(2)/tub_r
+      vectube(3)=vectube(3)/tub_r
+!C calculte rdiffrence between r and r0
+      rdiff=tub_r-tubeR0
+!C and its 6 power
+      rdiff6=rdiff**6.0d0
+!C for vectorization reasons we will sumup at the end to avoid depenence of previous
+       enetube(i)=pep_aa_tube/rdiff6**2.0d0+pep_bb_tube/rdiff6
+!C       write(iout,*) "TU13",i,rdiff6,enetube(i)
+!C       print *,rdiff,rdiff6,pep_aa_tube
+!C pep_aa_tube and pep_bb_tube are precomputed values A=4eps*sigma^12 B=4eps*sigma^6
+!C now we calculate gradient
+       fac=(-12.0d0*pep_aa_tube/rdiff6-   &
+            6.0d0*pep_bb_tube)/rdiff6/rdiff
+!C       write(iout,'(a5,i4,f12.1,3f12.5)') "TU13",i,rdiff6,enetube(i),
+!C     &rdiff,fac
+         if (acavtubpep.eq.0.0d0) then
+!C go to 667
+         enecavtube(i)=0.0
+         faccav=0.0
+         else
+         denominator=(1.0d0+dcavtubpep*rdiff6*rdiff6)
+         enecavtube(i)=  &
+        (bcavtubpep*rdiff+acavtubpep*dsqrt(rdiff)+ccavtubpep) &
+        /denominator
+         enecavtube(i)=0.0
+         faccav=((bcavtubpep*1.0d0+acavtubpep/2.0d0/dsqrt(rdiff)) &
+        *denominator-(bcavtubpep*rdiff+acavtubpep*dsqrt(rdiff)   &
+        +ccavtubpep)*rdiff6**2.0d0/rdiff*dcavtubpep*12.0d0)      &
+        /denominator**2.0d0
+!C         faccav=0.0
+!C         fac=fac+faccav
+!C 667     continue
+         endif
+
+        do j=1,3
+        gg_tube(j,i-1)=gg_tube(j,i-1)+vectube(j)*fac/2.0d0
+        gg_tube(j,i)=gg_tube(j,i)+vectube(j)*fac/2.0d0
+        enddo
+        enddo
+
+       do i=itube_start,itube_end
+        enecavtube(i)=0.0d0
+!C Lets not jump over memory as we use many times iti
+         iti=itype(i)
+!C lets ommit dummy atoms for now
+         if ((iti.eq.ntyp1) &
+!C in UNRES uncomment the line below as GLY has no side-chain...
+!C      .or.(iti.eq.10)
+         ) cycle
+      xmin=boxxsize
+      ymin=boxysize
+      zmin=boxzsize
+        do j=-1,1
+         vectube(1)=dmod((c(1,i+nres)),boxxsize)
+         vectube(1)=vectube(1)+boxxsize*j
+         vectube(2)=dmod((c(2,i+nres)),boxysize)
+         vectube(2)=vectube(2)+boxysize*j
+         vectube(3)=dmod((c(3,i+nres)),boxzsize)
+         vectube(3)=vectube(3)+boxzsize*j
+
+
+         xminact=dabs(vectube(1)-tubecenter(1))
+         yminact=dabs(vectube(2)-tubecenter(2))
+         zminact=dabs(vectube(3)-tubecenter(3))
+
+           if (xmin.gt.xminact) then
+            xmin=xminact
+            xtemp=vectube(1)
+           endif
+           if (ymin.gt.yminact) then
+             ymin=yminact
+             ytemp=vectube(2)
+            endif
+           if (zmin.gt.zminact) then
+             zmin=zminact
+             ztemp=vectube(3)
+            endif
+         enddo
+      vectube(1)=xtemp
+      vectube(2)=ytemp
+      vectube(3)=ztemp
+
+!C          write(iout,*), "tututu", vectube(1),tubecenter(1),vectube(2),
+!C     &     tubecenter(2)
+      vectube(1)=vectube(1)-tubecenter(1)
+      vectube(2)=vectube(2)-tubecenter(2)
+      vectube(3)=vectube(3)-tubecenter(3)
+!C now calculte the distance
+       tub_r=dsqrt(vectube(1)**2+vectube(2)**2+vectube(3)**2)
+!C now normalize vector
+      vectube(1)=vectube(1)/tub_r
+      vectube(2)=vectube(2)/tub_r
+      vectube(3)=vectube(3)/tub_r
+
+!C calculte rdiffrence between r and r0
+      rdiff=tub_r-tubeR0
+!C and its 6 power
+      rdiff6=rdiff**6.0d0
+       sc_aa_tube=sc_aa_tube_par(iti)
+       sc_bb_tube=sc_bb_tube_par(iti)
+       enetube(i+nres)=sc_aa_tube/rdiff6**2.0d0+sc_bb_tube/rdiff6
+!C       enetube(i+nres)=0.0d0
+!C pep_aa_tube and pep_bb_tube are precomputed values A=4eps*sigma^12 B=4eps*sigma^6
+!C now we calculate gradient
+       fac=-12.0d0*sc_aa_tube/rdiff6**2.0d0/rdiff- &
+            6.0d0*sc_bb_tube/rdiff6/rdiff
+!C       fac=0.0
+!C now direction of gg_tube vector
+!C Now cavity term E=a(x+bsqrt(x)+c)/(1+dx^12)
+         if (acavtub(iti).eq.0.0d0) then
+!C go to 667
+         enecavtube(i+nres)=0.0d0
+         faccav=0.0d0
+         else
+         denominator=(1.0d0+dcavtub(iti)*rdiff6*rdiff6)
+         enecavtube(i+nres)=   &
+        (bcavtub(iti)*rdiff+acavtub(iti)*dsqrt(rdiff)+ccavtub(iti)) &
+        /denominator
+!C         enecavtube(i)=0.0
+         faccav=((bcavtub(iti)*1.0d0+acavtub(iti)/2.0d0/dsqrt(rdiff)) &
+        *denominator-(bcavtub(iti)*rdiff+acavtub(iti)*dsqrt(rdiff)   &
+        +ccavtub(iti))*rdiff6**2.0d0/rdiff*dcavtub(iti)*12.0d0)      &
+        /denominator**2.0d0
+!C         faccav=0.0
+         fac=fac+faccav
+!C 667     continue
+         endif
+!C         print *,"TUT",i,iti,rdiff,rdiff6,acavtub(iti),denominator,
+!C     &   enecavtube(i),faccav
+!C         print *,"licz=",
+!C     & (bcavtub(iti)*rdiff+acavtub(iti)*sqrt(rdiff)+ccavtub(iti))
+!C         print *,"finene=",enetube(i+nres)+enecavtube(i)
+         do j=1,3
+          gg_tube_SC(j,i)=gg_tube_SC(j,i)+vectube(j)*fac
+          gg_tube(j,i-1)=gg_tube(j,i-1)+vectube(j)*fac
+         enddo
+        enddo
+
+
+
+        do i=itube_start,itube_end
+          Etube=Etube+enetube(i)+enetube(i+nres)+enecavtube(i) &
+         +enecavtube(i+nres)
+        enddo
+!C        print *,"ETUBE", etube
+        return
+        end subroutine calcnano
+
+!===============================================
 !--------------------------------------------------------------------------------
 !C first for shielding is setting of function of side-chains
 
@@ -18327,6 +18980,8 @@ write(iout,*) 'Calling CHECK_ECARTIN else.'
       allocate(gshieldc_ll(3,-1:nres))
       allocate(gshieldc_loc_ll(3,-1:nres))
       allocate(grad_shield(3,-1:nres))
+      allocate(gg_tube_sc(3,-1:nres))
+      allocate(gg_tube(3,-1:nres))
 !(3,maxres)
       allocate(grad_shield_side(3,50,nres))
       allocate(grad_shield_loc(3,50,nres))