Adam's corrections
authorCezary Czaplewski <czarek@chem.univ.gda.pl>
Fri, 22 May 2020 21:12:14 +0000 (23:12 +0200)
committerCezary Czaplewski <czarek@chem.univ.gda.pl>
Fri, 22 May 2020 21:12:14 +0000 (23:12 +0200)
source/cluster/wham/src-HCD/energy_p_new.F
source/unres/src-HCD-5D/MREMD.F
source/unres/src-HCD-5D/readrtns_CSA.F
source/wham/src-HCD/DIMENSIONS
source/wham/src-HCD/energy_p_new.F
source/wham/src-HCD/read_constr_homology.F
source/wham/src-HCD/readpdb.F

index 27e944b..95898b4 100644 (file)
          enddo
          
 c         min_odl=minval(distancek)
-         do kk=1,constr_homology
-          if(l_homo(kk,ii)) then 
-            min_odl=distancek(kk)
-            exit
-          endif
-         enddo
-         do kk=1,constr_homology
-          if(l_homo(kk,ii) .and. distancek(kk).lt.min_odl) 
+         if (nexl.gt.0) then
+           min_odl=0.0d0
+         else
+           do kk=1,constr_homology
+            if(l_homo(kk,ii)) then
+              min_odl=distancek(kk)
+              exit
+            endif
+           enddo
+           do kk=1,constr_homology
+            if(l_homo(kk,ii) .and. distancek(kk).lt.min_odl)
      &              min_odl=distancek(kk)
-         enddo
+           enddo
+         endif
 c        write (iout,* )"min_odl",min_odl
 #ifdef DEBUG
          write (iout,*) "ij dij",i,j,dij
index 38db8a8..9191402 100644 (file)
@@ -62,6 +62,7 @@ cold      integer nup(0:maxprocs),ndown(0:maxprocs)
      & ene_i_iex,xxx,tmp,econstr_temp_iex,econstr_temp_i
       integer iran_num
       double precision ran_number
+      integer i_econstr/20/
 
 cdeb      imin_itime_old=0
       ntwx_cache=0
@@ -186,8 +187,9 @@ cd           write (*,*) me," After broadcast: file_exist",file_exist
              do i=nnt,nct
               if (itype(i).ne.ntyp1) stdforcsc(i)=stdfsc(iabs(itype(i)))
      &                   *dsqrt(gamsc(iabs(itype(i))))
-             enddo
+              enddo
            endif
+
        endif
 
        if(me.eq.king) then
@@ -250,20 +252,20 @@ cd           write (*,*) me," After broadcast: file_exist",file_exist
               write (iout,'(a6,i4,a1,100i4)') "ndowna",il,":",
      &                    (ndowna(i,il),i=1,ndowna(0,il))
              enddo
+
              stdfp=dsqrt(2*Rb*t_bath/d_time)
              do i=1,ntyp
                stdfsc(i)=dsqrt(2*Rb*t_bath/d_time)
              enddo
              if (lang.gt.0 .and. .not.surfarea) then
-               do i=nnt,nct-1
-                 stdforcp(i)=stdfp*dsqrt(gamp)
-               enddo
+             do i=nnt,nct-1
+               stdforcp(i)=stdfp*dsqrt(gamp)
+             enddo
              do i=nnt,nct
               if (itype(i).ne.ntyp1) stdforcsc(i)=stdfsc(iabs(itype(i)))
      &                   *dsqrt(gamsc(iabs(itype(i))))
-               enddo
-             endif
-
+             enddo
+         endif
         ELSE IF (.not.(rest.and.file_exist)) THEN
          do il=1,remd_m(1)
           ifirst(il)=il
@@ -319,7 +321,7 @@ cd       print *,'ttt',me,remd_tlist,(remd_t(i),i=1,nrep)
           iset=i2set(me)
           if(me.eq.king.or..not.out1file) 
      &     write(iout,*) me,"iset=",iset,"t_bath=",t_bath
-          call flush(iout)
+c          call flush(iout)
        endif        
 c
        stdfp=dsqrt(2*Rb*t_bath/d_time)
@@ -721,7 +723,7 @@ ctime           call flush(iout)
                write(iout,*) "MIN ii_write=",ii_write
              endif
             endif
-           call flush(iout)
+c           call flush(iout)
         endif
         if (synflag) then
 c Update the time safety limiy
@@ -855,7 +857,7 @@ cd end
           time05=MPI_WTIME()
           if (me.eq.king .or. .not. out1file) then
             write(iout,*) 'REMD writing traj time=',time05-time04
-            call flush(iout)
+c            call flush(iout)
           endif
 
 
@@ -1089,12 +1091,12 @@ cd            write(iout,*) "i_dir=",i_dir
                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
 c 9/1/17 AL: Correction; otherwise the restraint energies are mis-assigned 
 c       on failed replica exchange attempt
-               econstr_temp_i=remd_ene(20,i)
-               econstr_temp_iex=remd_ene(20,iex)
+               econstr_temp_i=remd_ene(i_econstr,i)
+               econstr_temp_iex=remd_ene(i_econstr,iex)
 c 9/11/17 AL: Adaptive sampling (temperature dependent restraints potentials)
                if (adaptive) then
-                 remd_ene(20,i)=remd_ene(n_ene+5,i)
-                 remd_ene(20,iex)=remd_ene(n_ene+6,iex)
+                 remd_ene(i_econstr,i)=remd_ene(n_ene+5,i)
+                 remd_ene(i_econstr,iex)=remd_ene(n_ene+6,iex)
                endif
             elseif(i_dir.eq.2 .and. mset(i_iset+1).gt.0)then
 
@@ -1103,10 +1105,10 @@ c 9/11/17 AL: Adaptive sampling (temperature dependent restraints potentials)
                i_iset1=i_iset+1
                i_mset1=iran_num(1,mset(i_iset1))
                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
-               econstr_temp_i=remd_ene(20,i)
-               econstr_temp_iex=remd_ene(20,iex)
-               remd_ene(20,i)=remd_ene(n_ene+3,i)
-               remd_ene(20,iex)=remd_ene(n_ene+4,iex)
+               econstr_temp_i=remd_ene(i_econstr,i)
+               econstr_temp_iex=remd_ene(i_econstr,iex)
+               remd_ene(i_econstr,i)=remd_ene(n_ene+3,i)
+               remd_ene(i_econstr,iex)=remd_ene(n_ene+4,iex)
 
             elseif(remd_m(i_temp+1).gt.0.and.mset(i_iset+1).gt.0)then
 
@@ -1115,14 +1117,14 @@ c 9/11/17 AL: Adaptive sampling (temperature dependent restraints potentials)
                i_iset1=i_iset+1
                i_mset1=iran_num(1,mset(i_iset1))
                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
-               econstr_temp_i=remd_ene(20,i)
-               econstr_temp_iex=remd_ene(20,iex)
+               econstr_temp_i=remd_ene(i_econstr,i)
+               econstr_temp_iex=remd_ene(i_econstr,iex)
                if (adaptive) then
-                 remd_ene(20,i)=remd_ene(n_ene+7,i)
-                 remd_ene(20,iex)=remd_ene(n_ene+8,iex)
+                 remd_ene(i_econstr,i)=remd_ene(n_ene+7,i)
+                 remd_ene(i_econstr,iex)=remd_ene(n_ene+8,iex)
                else
-                 remd_ene(20,i)=remd_ene(n_ene+3,i)
-                 remd_ene(20,iex)=remd_ene(n_ene+4,iex)
+                 remd_ene(i_econstr,i)=remd_ene(n_ene+3,i)
+                 remd_ene(i_econstr,iex)=remd_ene(n_ene+4,iex)
                endif
             else
                goto 444 
@@ -1222,8 +1224,8 @@ cd              write(iout,'(2i4,a6,2f12.5)') i,iex,' delta',delta,xxx
               else
                remd_ene(0,iex)=ene_iex_iex
                remd_ene(0,i)=ene_i_i
-               remd_ene(20,iex)=econstr_temp_iex
-               remd_ene(20,i)=econstr_temp_i
+               remd_ene(i_econstr,iex)=econstr_temp_iex
+               remd_ene(i_econstr,i)=econstr_temp_i
               endif
 
 cd      do il=1,nset
@@ -1259,7 +1261,7 @@ c-------------------------------------
               enddo
              endif
 
-             call flush(iout)
+c             call flush(iout)
 
 cd              write (iout,'(a6,100i4)') "ifirst",
 cd     &                    (ifirst(i),i=1,remd_m(1))
@@ -1305,7 +1307,6 @@ co     &    " rescaling weights with temperature",t_bath
          do i=1,ntyp
            stdfsc(i)=dsqrt(2*Rb*t_bath/d_time)
          enddo
-
 c Compute the standard deviations of stochastic forces for Langevin dynamics
 c if the friction coefficients do not depend on surface area
          if (lang.gt.0 .and. .not.surfarea) then
@@ -1317,11 +1318,11 @@ c if the friction coefficients do not depend on surface area
      &                   *dsqrt(gamsc(iabs(itype(i))))
            enddo
          endif
-cde         write(iout,*) 'REMD after',me,t_bath
+cde           write(iout,*) 'REMD after',me,t_bath
            time08=MPI_WTIME()
            if (me.eq.king .or. .not. out1file) then
             write(iout,*) 'REMD exchange time=',time08-time00
-            call flush(iout)
+c            call flush(iout)
            endif
         endif
       enddo
index c51e45e..16c0f37 100644 (file)
@@ -2979,10 +2979,11 @@ c    &    sigma_odl_temp(maxres,maxres,max_template)
       character*24 model_ki_dist, model_ki_angle
       character*500 controlcard
       integer ki,i,ii,j,k,l,ii_in_use(maxdim),i_tmp,idomain_tmp,irec,
-     & ik,iistart,iishift
+     & ik,iistart
       integer ilen
       external ilen
-      logical liiflag
+      logical liiflag,lfirst
+      integer i01,i10
 c
 c     FP - Nov. 2014 Temporary specifications for new vars
 c
@@ -3314,6 +3315,7 @@ c      write (iout,*) "waga_dist",waga_dist," nnt",nnt," nct",nct
       if (waga_dist.ne.0.0d0) then
         ii=0
         liiflag=.true.
+        lfirst=.true.
         do i=nnt,nct-2 
          do j=i+2,nct 
           ii=ii+1
@@ -3321,34 +3323,35 @@ c          if (idomain(k,i).eq.idomain(k,j).and.idomain(k,i).ne.0
 c     &            .and. distal.le.dist2_cut ) then
 c          write (iout,*) "i",i," j",j," ii",ii
 c          call flush(iout)
-          if (ii_in_use(ii).eq.0.and.liiflag) then
+          if (ii_in_use(ii).eq.0.and.liiflag.or.
+     &     ii_in_use(ii).eq.1.and.liiflag.and.ii.eq.lim_odl) then
             liiflag=.false.
-            iistart=ii
+            i10=ii
+            if (lfirst) then
+              lfirst=.false.
+              iistart=ii
+            else
+              if(i10.eq.lim_odl) i10=i10+1
+              do ki=0,i10-i01-1
+               ires_homo(iistart+ki)=ires_homo(ki+i01)
+               jres_homo(iistart+ki)=jres_homo(ki+i01)
+               ii_in_use(iistart+ki)=ii_in_use(ki+i01)
+               do k=1,constr_homology
+                odl(k,iistart+ki)=odl(k,ki+i01)
+                sigma_odl(k,iistart+ki)=sigma_odl(k,ki+i01)
+                l_homo(k,iistart+ki)=l_homo(k,ki+i01)
+               enddo
+              enddo
+              iistart=iistart+i10-i01
+            endif
           endif
-          if (ii_in_use(ii).ne.0.and..not.liiflag.or.
-     &                   .not.liiflag.and.ii.eq.lim_odl) then
-             if (ii.eq.lim_odl) then
-              iishift=ii-iistart+1
-             else
-              iishift=ii-iistart
-             endif
+          if (ii_in_use(ii).ne.0.and..not.liiflag) then
+             i01=ii
              liiflag=.true.
-             do ki=iistart,lim_odl-iishift
-              ires_homo(ki)=ires_homo(ki+iishift)
-              jres_homo(ki)=jres_homo(ki+iishift)
-              ii_in_use(ki)=ii_in_use(ki+iishift)
-              do k=1,constr_homology
-               odl(k,ki)=odl(k,ki+iishift)
-               sigma_odl(k,ki)=sigma_odl(k,ki+iishift)
-               l_homo(k,ki)=l_homo(k,ki+iishift)
-              enddo
-             enddo
-             ii=ii-iishift
-             lim_odl=lim_odl-iishift
-c          endif
           endif
          enddo
         enddo
+        lim_odl=iistart-1
       endif
 c      write (iout,*) "Removing distances completed"
 c      call flush(iout)
index 4bc527f..4d690f3 100644 (file)
@@ -24,7 +24,7 @@ c Max. number of chains
       parameter (maxchain=6)
 C Max number of symetries
        integer maxsym,maxperm
-       parameter (maxsym=maxchain,maxperm=1200)
+       parameter (maxsym=maxchain,maxperm=720)
 C Max. number of variables
       integer maxvar
       parameter (maxvar=4*maxres)
index bd69774..efba869 100644 (file)
@@ -159,6 +159,7 @@ c         write (iout,*) "Calling multibody_hbond"
          call multibody_hb(ecorr,ecorr5,ecorr6,n_corr,n_corr1)
       endif
 #endif
+c      write (iout,*) "nsaxs",nsaxs
 c      write (iout,*) "From Esaxs: Esaxs_constr",Esaxs_constr
       if (nsaxs.gt.0 .and. saxs_mode.eq.0) then
         call e_saxs(Esaxs_constr)
          enddo
          
 c         min_odl=minval(distancek)
-         do kk=1,constr_homology
-          if(l_homo(kk,ii)) then 
-            min_odl=distancek(kk)
-            exit
-          endif
-         enddo
-         do kk=1,constr_homology
-          if(l_homo(kk,ii) .and. distancek(kk).lt.min_odl) 
+         if (nexl.gt.0) then
+           min_odl=0.0d0
+         else
+           do kk=1,constr_homology
+            if(l_homo(kk,ii)) then
+              min_odl=distancek(kk)
+              exit
+            endif
+           enddo
+           do kk=1,constr_homology
+            if(l_homo(kk,ii) .and. distancek(kk).lt.min_odl)
      &              min_odl=distancek(kk)
-         enddo
+           enddo
+         endif
 c        write (iout,* )"min_odl",min_odl
 #ifdef DEBUG
          write (iout,*) "ij dij",i,j,dij
index dece50f..0af2bdb 100644 (file)
@@ -165,6 +165,7 @@ c        call readpdb
         write (iout,*) "read_constr_homology: after reading pdb file"
         call flush(iout)
 #endif
+
 c
 c     Distance restraints
 c
@@ -553,14 +554,18 @@ c     Distance restraints
 c
 c          ... --> odl(k,ii)
 C Copy the coordinates from reference coordinates (?)
+c        write (iout,*) "k",k
         do i=1,2*nres
           do j=1,3
             c(j,i)=chomo(j,i,k)
-c           write (iout,*) "c(",j,i,") =",c(j,i)
+c            write (iout,*) "c(",j,i,") =",c(j,i)
           enddo
         enddo
+c        call cartprint
+c        write (iout,*) "read_klapaucjusz: calling int_from_cart"
         call int_from_cart(.true.,.false.)
         call sc_loc_geom(.false.)
+c        write (iout,*) "en" 
         do i=1,nres
           thetaref(i)=theta(i)
           phiref(i)=phi(i)
index 83c63ab..55dcc04 100644 (file)
@@ -27,10 +27,9 @@ C geometry.
       ibeg=1
       ishift1=0
       sccalc=.false.
-      ires=0
       do
         read (ipdbin,'(a80)',end=10) card
-c        write (iout,'(a)') card
+!       write (iout,'(a)') card
         if (card(:5).eq.'HELIX') then
           nhfrag=nhfrag+1
           lsecondary=.true.
@@ -59,7 +58,7 @@ c        write (iout,'(a)') card
           iterter(ires_old)=1
           ishift1=ishift1+1
           ibeg=2
-          write (iout,*) "Chain ended",ires,ishift,ires_old
+!          write (iout,*) "Chain ended",ires,ishift,ires_old
           if (unres_pdb) then
             do j=1,3
               dc(j,ires)=sccor(j,iii)
@@ -120,7 +119,7 @@ c              write (iout,*) "BEG ires",ires
 ! Start a new chain
               ishift=-ires_old+ires-1 !!!!!
               ishift1=ishift1-1    !!!!!
-              write (iout,*) "New chain started",ires,ishift,ishift1,"!"
+!              write (iout,*) "New chain started",ires,ishift,ishift1,"!"
               ires=ires-ishift+ishift1
               ires_old=ires
               ibeg=0
@@ -326,7 +325,7 @@ c        write (iout,*) i,i-1,(c(j,i),j=1,3),(c(j,i-1),j=1,3),dist(i,i-1)
         if (itype(i-1).ne.ntyp1 .and. itype(i).ne.ntyp1 .and.
      &    (dist(i,i-1).lt.1.0D0 .or. dist(i,i-1).gt.6.0D0)) then
           write (iout,'(a,i4)') 'Bad Cartesians for residue',i
-          stop
+c          stop
         endif
         vbld(i)=dist(i-1,i)
         vbld_inv(i)=1.0d0/vbld(i)