X-Git-Url: http://mmka.chem.univ.gda.pl/gitweb/?a=blobdiff_plain;f=source%2Funres%2Fsrc_MD-M%2Fenergy_p_new_barrier.F;h=ff2d579c5db15cd9f2893b76c5bff389635153e8;hb=8999b2afa91f8ac4948fa5ea826253408e17e8d4;hp=91353c14e34e4cd66e44c2b1fbd37b679b714fd6;hpb=bda982c0a4335703f4b46cdfa25a7021ceb6dbdf;p=unres.git diff --git a/source/unres/src_MD-M/energy_p_new_barrier.F b/source/unres/src_MD-M/energy_p_new_barrier.F index 91353c1..ff2d579 100644 --- a/source/unres/src_MD-M/energy_p_new_barrier.F +++ b/source/unres/src_MD-M/energy_p_new_barrier.F @@ -55,6 +55,7 @@ C FG slaves as WEIGHTS array. weights_(17)=wbond weights_(18)=scal14 weights_(21)=wsccor + weights_(25)=wsaxs C FG Master broadcasts the WEIGHTS_ array call MPI_Bcast(weights_(1),n_ene, & MPI_DOUBLE_PRECISION,king,FG_COMM,IERROR) @@ -81,6 +82,7 @@ C FG slaves receive the WEIGHTS array wbond=weights(17) scal14=weights(18) wsccor=weights(21) + wsaxs=weights(25) endif time_Bcast=time_Bcast+MPI_Wtime()-time00 time_Bcastw=time_Bcastw+MPI_Wtime()-time00 @@ -214,6 +216,18 @@ cd print *,'nterm=',nterm etors=0 edihcnstr=0 endif + + if (constr_homology.ge.1) then + call e_modeller(ehomology_constr) +c print *,'iset=',iset,'me=',me,ehomology_constr, +c & 'Processor',fg_rank,' CG group',kolor, +c & ' absolute rank',MyRank + else + ehomology_constr=0.0d0 + endif + + +c write(iout,*) ehomology_constr c print *,"Processor",myrank," computed Utor" C C 6/23/01 Calculate double-torsional energy @@ -254,6 +268,16 @@ cd &" ecorr",ecorr," ecorr5",ecorr5," ecorr6",ecorr6," eturn6",eturn6 call multibody_hb(ecorr,ecorr5,ecorr6,n_corr,n_corr1) cd write (iout,*) "multibody_hb ecorr",ecorr endif +c write (iout,*) "nsaxs",nsaxs," saxs_mode",saxs_mode + if (nsaxs.gt.0 .and. saxs_mode.eq.0) then + call e_saxs(Esaxs_constr) +c write (iout,*) "From Esaxs: Esaxs_constr",Esaxs_constr + else if (nsaxs.gt.0 .and. saxs_mode.gt.0) then + call e_saxsC(Esaxs_constr) +c write (iout,*) "From EsaxsC: Esaxs_constr",Esaxs_constr + else + Esaxs_constr = 0.0d0 + endif c print *,"Processor",myrank," computed Ucorr" C C If performing constraint dynamics, call the constraint energy @@ -321,6 +345,8 @@ C energia(21)=esccor energia(22)=eliptran energia(23)=Eafmforce + energia(24)=ehomology_constr + energia(25)=Esaxs_constr c Here are the energies showed per procesor if the are more processors c per molecule then we sum it up in sum_energy subroutine c print *," Processor",myrank," calls SUM_ENERGY" @@ -414,20 +440,26 @@ cMS$ATTRIBUTES C :: proc_proc esccor=energia(21) eliptran=energia(22) Eafmforce=energia(23) + ehomology_constr=energia(24) + esaxs_constr=energia(25) +c write (iout,*) "sum_energy esaxs_constr",esaxs_constr, +c & " wsaxs",wsaxs #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+Eafmforce + & +wbond*estr+Uconst+wsccor*esccor+ehomology_constr + & +wsaxs*esaxs_constr+wliptran*eliptran+Eafmforce #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+ehomology_constr + & +wsaxs*esaxs_constr+wliptran*eliptran & +Eafmforce #endif energia(0)=etot @@ -480,17 +512,26 @@ cMS$ATTRIBUTES C :: proc_proc include 'COMMON.TIME1' include 'COMMON.MAXGRAD' include 'COMMON.SCCOR' + include 'COMMON.MD' #ifdef TIMING time01=MPI_Wtime() #endif #ifdef DEBUG write (iout,*) "sum_gradient gvdwc, gvdwx" - do i=1,nres - write (iout,'(i3,3f10.5,5x,3f10.5,5x,f10.5)') + do i=0,nres + write (iout,'(i3,3e15.5,5x,3e15.5)') & i,(gvdwx(j,i),j=1,3),(gvdwc(j,i),j=1,3) enddo call flush(iout) #endif +#ifdef DEBUG + write (iout,*) "sum_gradient gsaxsc, gsaxsx" + do i=0,nres + write (iout,'(i3,3e15.5,5x,3e15.5)') + & i,(gsaxsc(j,i),j=1,3),(gsaxsx(j,i),j=1,3) + enddo + call flush(iout) +#endif #ifdef MPI C FG slaves call the following matching MPI_Bcast in ERGASTULUM if (nfgtasks.gt.1 .and. fg_rank.eq.0) @@ -530,7 +571,8 @@ c enddo & wcorr5*gradcorr5_long(j,i)+ & wcorr6*gradcorr6_long(j,i)+ & wturn6*gcorr6_turn_long(j,i)+ - & wstrain*ghpbc(j,i) + & wstrain*ghpbc(j,i)+ + & wsaxs*gsaxsc(j,i) & +wliptran*gliptranc(j,i) & +gradafm(j,i) @@ -548,7 +590,8 @@ c enddo & wcorr5*gradcorr5_long(j,i)+ & wcorr6*gradcorr6_long(j,i)+ & wturn6*gcorr6_turn_long(j,i)+ - & wstrain*ghpbc(j,i) + & wstrain*ghpbc(j,i)+ + & wsaxs*gsaxsc(j,i) & +wliptran*gliptranc(j,i) & +gradafm(j,i) @@ -608,7 +651,8 @@ c enddo do j=1,3 gradbufc(j,nres-1)=gradbufc_sum(j,nres) enddo - do i=nres-2,-1,-1 + do i=nres-2,0,-1 +c do i=nres-2,nnt,-1 do j=1,3 gradbufc(j,i)=gradbufc(j,i+1)+gradbufc_sum(j,i+1) enddo @@ -616,7 +660,7 @@ c enddo #ifdef DEBUG write (iout,*) "gradbufc after summing" do i=1,nres - write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3) + write (iout,'(i3,3e15.5)') i,(gradbufc(j,i),j=1,3) enddo call flush(iout) #endif @@ -624,8 +668,8 @@ c enddo #endif #ifdef DEBUG write (iout,*) "gradbufc" - do i=1,nres - write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3) + do i=0,nres + write (iout,'(i3,3e15.5)') i,(gradbufc(j,i),j=1,3) enddo call flush(iout) #endif @@ -638,7 +682,8 @@ c enddo do j=1,3 gradbufc(j,nres-1)=gradbufc_sum(j,nres) enddo - do i=nres-2,-1,-1 + do i=nres-2,0,-1 +c do i=nres-2,nnt,-1 do j=1,3 gradbufc(j,i)=gradbufc(j,i+1)+gradbufc_sum(j,i+1) enddo @@ -655,8 +700,8 @@ c enddo c enddo #ifdef DEBUG write (iout,*) "gradbufc after summing" - do i=1,nres - write (iout,'(i3,3f10.5)') i,(gradbufc(j,i),j=1,3) + do i=0,nres + write (iout,'(i3,3e15.5)') i,(gradbufc(j,i),j=1,3) enddo call flush(iout) #endif @@ -714,12 +759,21 @@ c enddo #endif gradx(j,i,icg)=wsc*gvdwx(j,i)+wscp*gradx_scp(j,i)+ & wbond*gradbx(j,i)+ - & wstrain*ghpbx(j,i)+wcorr*gradxorr(j,i)+ - & wsccor*gsccorx(j,i) + & wstrain*ghpbx(j,i)+wcorr*gradxorr(j,i) + & +wsaxs*gsaxsx(j,i) + & +wsccor*gsccorx(j,i) & +wscloc*gsclocx(j,i) & +wliptran*gliptranx(j,i) enddo enddo + if (constr_homology.gt.0) then + do i=1,nct + do j=1,3 + gradc(j,i,icg)=gradc(j,i,icg)+duscdiff(j,i) + gradx(j,i,icg)=gradx(j,i,icg)+duscdiffx(j,i) + enddo + enddo + endif #ifdef DEBUG write (iout,*) "gloc before adding corr" do i=1,4*nres @@ -904,7 +958,7 @@ c #ifdef DEBUG write (iout,*) "gradc gradx gloc" do i=1,nres - write (iout,'(i5,3f10.5,5x,3f10.5,5x,f10.5)') + write (iout,'(i5,3e15.5,5x,3e15.5,5x,e15.5)') & i,(gradc(j,i,icg),j=1,3),(gradx(j,i,icg),j=1,3),gloc(i,icg) enddo #endif @@ -1006,6 +1060,8 @@ C------------------------------------------------------------------------ estr=energia(17) Uconst=energia(20) esccor=energia(21) + ehomology_constr=energia(24) + esaxs_constr=energia(25) eliptran=energia(22) Eafmforce=energia(23) #ifdef SPLITELE @@ -1015,7 +1071,7 @@ C------------------------------------------------------------------------ & ecorr,wcorr, & ecorr5,wcorr5,ecorr6,wcorr6,eel_loc,wel_loc,eello_turn3,wturn3, & eello_turn4,wturn4,eello_turn6,wturn6,esccor,wsccor, - & edihcnstr,ebr*nss, + & edihcnstr,ehomology_constr,esaxs_constr*wsaxs, ebr*nss, & Uconst,eliptran,wliptran,Eafmforce,etot 10 format (/'Virtual-chain energies:'// & 'EVDW= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-SC)'/ @@ -1027,7 +1083,7 @@ C------------------------------------------------------------------------ & 'ESC= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC local)'/ & 'ETORS= ',1pE16.6,' WEIGHT=',1pD16.6,' (torsional)'/ & 'ETORSD=',1pE16.6,' WEIGHT=',1pD16.6,' (double torsional)'/ - & 'EHBP= ',1pE16.6,' WEIGHT=',1pD16.6, + & 'EHPB= ',1pE16.6,' WEIGHT=',1pD16.6, & ' (SS bridges & dist. cnstr.)'/ & 'ECORR4=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/ & 'ECORR5=',1pE16.6,' WEIGHT=',1pD16.6,' (multi-body)'/ @@ -1038,6 +1094,8 @@ C------------------------------------------------------------------------ & 'ETURN6=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 6th order)'/ & 'ESCCOR=',1pE16.6,' WEIGHT=',1pD16.6,' (backbone-rotamer corr)'/ & 'EDIHC= ',1pE16.6,' (dihedral angle constraints)'/ + & 'H_CONS=',1pE16.6,' (Homology model constraints energy)'/ + & 'E_SAXS=',1pE16.6,' (SAXS restraints)'/ & 'ESS= ',1pE16.6,' (disulfide-bridge intrinsic energy)'/ & 'UCONST= ',1pE16.6,' (Constraint energy)'/ & 'ELT=',1pE16.6, ' WEIGHT=',1pD16.6,' (Lipid transfer energy)'/ @@ -1051,7 +1109,9 @@ C------------------------------------------------------------------------ & ecorr,wcorr, & ecorr5,wcorr5,ecorr6,wcorr6,eel_loc,wel_loc,eello_turn3,wturn3, & eello_turn4,wturn4,eello_turn6,wturn6,esccor,wsccro,edihcnstr, - & ebr*nss,Uconst,eliptran,wliptran,Eafmforc,etot + & ehomology_constr,esaxs_constr*wsaxs,ebr*nss,Uconst, + & eliptran,wliptran,Eafmforc, + & etot 10 format (/'Virtual-chain energies:'// & 'EVDW= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-SC)'/ & 'EVDW2= ',1pE16.6,' WEIGHT=',1pD16.6,' (SC-p)'/ @@ -1072,6 +1132,8 @@ C------------------------------------------------------------------------ & 'ETURN6=',1pE16.6,' WEIGHT=',1pD16.6,' (turns, 6th order)'/ & 'ESCCOR=',1pE16.6,' WEIGHT=',1pD16.6,' (backbone-rotamer corr)'/ & 'EDIHC= ',1pE16.6,' (dihedral angle constraints)'/ + & 'H_CONS=',1pE16.6,' (Homology model constraints energy)'/ + & 'E_SAXS=',1pE16.6,' (SAXS restraints)'/ & 'ESS= ',1pE16.6,' (disulfide-bridge intrinsic energy)'/ & 'UCONST=',1pE16.6,' (Constraint energy)'/ & 'ELT=',1pE16.6, ' WEIGHT=',1pD16.6,' (Lipid transfer energy)'/ @@ -2830,7 +2892,8 @@ c if (i.gt. iatel_s+1 .and. i.lt.iatel_e+4) then cd write (iout,*) '*******i',i,' iti1',iti cd write (iout,*) 'b1',b1(:,iti) cd write (iout,*) 'b2',b2(:,iti) -cd write (iout,*) 'Ug',Ug(:,:,i-2) +cd write (iout,*) "phi(",i,")=",phi(i)," sin1",sin1," cos1",cos1 +cd write (iout,*) 'Ug',Ug(:,:,i-2) c if (i .gt. iatel_s+2) then if (i .gt. nnt+2) then call matvec2(Ug(1,1,i-2),b2(1,i-2),Ub2(1,i-2)) @@ -2885,7 +2948,11 @@ c if (i.gt. iatel_s+1 .and. i.lt.iatel_e+4) then mu(k,i-2)=Ub2(k,i-2)+b1(k,i-1) enddo C write (iout,*) 'mumu',i,b1(1,i-1),Ub2(1,i-2) -c write (iout,*) 'mu ',mu(:,i-2),i-2 +cd write (iout,*) 'mu ',mu(:,i-2),i-2 +cd write (iout,*) 'b1 ',b1(:,i-1),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(:,i-2),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) @@ -3294,21 +3361,21 @@ C Loop over i,i+2 and i,i+3 pairs of the peptide groups C C 14/01/2014 TURN3,TUNR4 does no go under periodic boundry condition do i=iturn3_start,iturn3_end - if (i.le.1) cycle +CAna if (i.le.1) cycle C write(iout,*) "tu jest i",i if (itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1 C changes suggested by Ana to avoid out of bounds - & .or.((i+4).gt.nres) - & .or.((i-1).le.0) +CAna & .or.((i+4).gt.nres) +CAna & .or.((i-1).le.0) C end of changes by Ana & .or. itype(i+2).eq.ntyp1 & .or. itype(i+3).eq.ntyp1) cycle - if(i.gt.1)then - if(itype(i-1).eq.ntyp1)cycle - end if - if(i.LT.nres-3)then - if (itype(i+4).eq.ntyp1) cycle - end if +CAna if(i.gt.1)then +CAna if(itype(i-1).eq.ntyp1)cycle +CAna end if +CAna if(i.LT.nres-3)then +CAna if (itype(i+4).eq.ntyp1) cycle +CAna end if dxi=dc(1,i) dyi=dc(2,i) dzi=dc(3,i) @@ -3330,17 +3397,17 @@ C end of changes by Ana num_cont_hb(i)=num_conti enddo do i=iturn4_start,iturn4_end - if (i.le.1) cycle +cAna if (i.le.1) cycle if (itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1 C changes suggested by Ana to avoid out of bounds - & .or.((i+5).gt.nres) - & .or.((i-1).le.0) +cAna & .or.((i+5).gt.nres) +cAna & .or.((i-1).le.0) C end of changes suggested by Ana & .or. itype(i+3).eq.ntyp1 & .or. itype(i+4).eq.ntyp1 - & .or. itype(i+5).eq.ntyp1 - & .or. itype(i).eq.ntyp1 - & .or. itype(i-1).eq.ntyp1 +cAna & .or. itype(i+5).eq.ntyp1 +cAna & .or. itype(i).eq.ntyp1 +cAna & .or. itype(i-1).eq.ntyp1 & ) cycle dxi=dc(1,i) dyi=dc(2,i) @@ -3398,14 +3465,14 @@ c c Loop over all pairs of interacting peptide groups except i,i+2 and i,i+3 c do i=iatel_s,iatel_e - if (i.le.1) cycle +cAna if (i.le.1) cycle if (itype(i).eq.ntyp1 .or. itype(i+1).eq.ntyp1 C changes suggested by Ana to avoid out of bounds - & .or.((i+2).gt.nres) - & .or.((i-1).le.0) +cAna & .or.((i+2).gt.nres) +cAna & .or.((i-1).le.0) C end of changes by Ana - & .or. itype(i+2).eq.ntyp1 - & .or. itype(i-1).eq.ntyp1 +cAna & .or. itype(i+2).eq.ntyp1 +cAna & .or. itype(i-1).eq.ntyp1 & ) cycle dxi=dc(1,i) dyi=dc(2,i) @@ -3456,14 +3523,14 @@ c write (iout,*) 'i',i,' ielstart',ielstart(i),' ielend',ielend(i) num_conti=num_cont_hb(i) do j=ielstart(i),ielend(i) C write (iout,*) i,j - if (j.le.1) cycle +cAna if (j.le.1) cycle if (itype(j).eq.ntyp1.or. itype(j+1).eq.ntyp1 C changes suggested by Ana to avoid out of bounds - & .or.((j+2).gt.nres) - & .or.((j-1).le.0) +cAna & .or.((j+2).gt.nres) +cAna & .or.((j-1).le.0) C end of changes by Ana - & .or.itype(j+2).eq.ntyp1 - & .or.itype(j-1).eq.ntyp1 +cAna & .or.itype(j+2).eq.ntyp1 +cAna & .or.itype(j-1).eq.ntyp1 &) cycle call eelecij(i,j,ees,evdw1,eel_loc) enddo ! j @@ -3649,7 +3716,7 @@ cd & xmedi,ymedi,zmedi,xj,yj,zj if (energy_dec) then write (iout,'(a6,2i5,0pf7.3,2i5,2e11.3)') &'evdw1',i,j,evdwij - &,iteli,itelj,aaa,evdw1 +c &,iteli,itelj,aaa,evdw1 write (iout,'(a6,2i5,0pf7.3)') 'ees',i,j,eesij endif @@ -4358,6 +4425,7 @@ C Derivatives in gamma(i+1) gel_loc_turn3(i+1)=gel_loc_turn3(i+1) & +0.5d0*(pizda(1,1)+pizda(2,2)) C Cartesian derivatives +!DIR$ UNROLL(0) do l=1,3 c ghalf1=0.5d0*agg(l,1) c ghalf2=0.5d0*agg(l,2) @@ -4514,8 +4582,8 @@ c i+3 eello_turn4=eello_turn4-(s1+s2+s3) c write(iout,*)'chujOWO', auxvec(1),b1(1,iti2) - if (energy_dec) write (iout,'(a6,2i5,0pf7.3,3f7.3)') - & 'eturn4',i,j,-(s1+s2+s3),s1,s2,s3 +c if (energy_dec) write (iout,'(a6,2i5,0pf7.3,3f7.3)') +c & 'eturn4',i,j,-(s1+s2+s3),s1,s2,s3 cd write (2,*) 'i,',i,' j',j,'eello_turn4',-(s1+s2+s3), cd & ' eello_turn4_num',8*eello_turn4_num #ifdef NEWCORR @@ -5041,8 +5109,8 @@ c if (sss.eq.0) print *,'czasem jest OK' evdwij=e1+e2 evdw2=evdw2+evdwij*sss if (energy_dec) write (iout,'(a6,2i5,0pf7.3,2i3,3e11.3)') - & 'evdw2',i,j,evdwij,iteli,itypj,fac,aad(itypj,iteli), - & bad(itypj,iteli) + & 'evdw2',i,j,evdwij +c & ,iteli,itypj,fac,aad(itypj,iteli),bad(itypj,iteli) C C Calculate contributions to the gradient in the virtual-bond and SC vectors. C @@ -5654,7 +5722,7 @@ C logical lprn /.false./, lprn1 /.false./ etheta=0.0D0 do i=ithet_start,ithet_end - if (i.eq.2) cycle +c if (i.eq.2) cycle c print *,i,itype(i-1),itype(i),itype(i-2) if ((itype(i-1).eq.ntyp1).or.(itype(i-2).eq.ntyp1) & .or.(itype(i).eq.ntyp1)) cycle @@ -5671,7 +5739,7 @@ C In current verion the ALL DUMMY ATOM POTENTIALS ARE OFF coskt(k)=dcos(k*theti2) sinkt(k)=dsin(k*theti2) enddo - if (i.gt.3 .and. itype(i-3).ne.ntyp1) 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 @@ -5686,7 +5754,7 @@ C propagation of chirality for glycine type 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 @@ -5707,7 +5775,7 @@ C propagation of chirality for glycine type enddo else phii1=0.0d0 - ityp3=nthetyp+1 + ityp3=ithetyp(itype(i)) do k=1,nsingle cosph2(k)=0.0d0 sinph2(k)=0.0d0 @@ -6265,6 +6333,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,it,itype(i) c & ,zz,xx,yy c#define DEBUG @@ -6661,6 +6731,15 @@ c write (iout,*) 'i=',i,' gloc=',gloc(i-3,icg) return end c------------------------------------------------------------------------------ +c LICZENIE WIEZOW Z ROWNANIA ENERGII MODELLERA + subroutine e_modeller(ehomology_constr) + ehomology_constr=0.0d0 + write (iout,*) "!!!!!UWAGA, JESTEM W DZIWNEJ PETLI, TEST!!!!!" + return + end +C !!!!!!!! NIE CZYTANE !!!!!!!!!!! + +c------------------------------------------------------------------------------ subroutine etor_d(etors_d) etors_d=0.0d0 return @@ -6775,6 +6854,611 @@ cd write (iout,*) 'edihcnstr',edihcnstr return end c---------------------------------------------------------------------------- +c MODELLER restraint function + subroutine e_modeller(ehomology_constr) + implicit real*8 (a-h,o-z) + include 'DIMENSIONS' + + integer nnn, i, j, k, ki, irec, l + integer katy, odleglosci, test7 + real*8 odleg, odleg2, odleg3, kat, kat2, kat3, gdih(max_template) + real*8 Eval,Erot + real*8 distance(max_template),distancek(max_template), + & min_odl,godl(max_template),dih_diff(max_template) + +c +c FP - 30/10/2014 Temporary specifications for homology restraints +c + double precision utheta_i,gutheta_i,sum_gtheta,sum_sgtheta, + & sgtheta + double precision, dimension (maxres) :: guscdiff,usc_diff + double precision, dimension (max_template) :: + & gtheta,dscdiff,uscdiffk,guscdiff2,guscdiff3, + & theta_diff +c + + include 'COMMON.SBRIDGE' + include 'COMMON.CHAIN' + include 'COMMON.GEO' + include 'COMMON.DERIV' + include 'COMMON.LOCAL' + include 'COMMON.INTERACT' + include 'COMMON.VAR' + include 'COMMON.IOUNITS' + include 'COMMON.MD' + include 'COMMON.CONTROL' +c +c From subroutine Econstr_back +c + include 'COMMON.NAMES' + include 'COMMON.TIME1' +c + + + do i=1,19 + distancek(i)=9999999.9 + enddo + + + odleg=0.0d0 + +c Pseudo-energy and gradient from homology restraints (MODELLER-like +c function) +C AL 5/2/14 - Introduce list of restraints +c write(iout,*) "waga_theta",waga_theta,"waga_d",waga_d +#ifdef DEBUG + write(iout,*) "------- dist restrs start -------" +#endif + do ii = link_start_homo,link_end_homo + i = ires_homo(ii) + j = jres_homo(ii) + dij=dist(i,j) +c write (iout,*) "dij(",i,j,") =",dij + do k=1,constr_homology +c write(iout,*) ii,k,i,j,l_homo(k,ii),dij,odl(k,ii) + if(.not.l_homo(k,ii)) cycle + distance(k)=odl(k,ii)-dij +c write (iout,*) "distance(",k,") =",distance(k) +c +c For Gaussian-type Urestr +c + distancek(k)=0.5d0*distance(k)**2*sigma_odl(k,ii) ! waga_dist rmvd from Gaussian argument +c write (iout,*) "sigma_odl(",k,ii,") =",sigma_odl(k,ii) +c write (iout,*) "distancek(",k,") =",distancek(k) +c distancek(k)=0.5d0*waga_dist*distance(k)**2*sigma_odl(k,ii) +c +c For Lorentzian-type Urestr +c + if (waga_dist.lt.0.0d0) then + sigma_odlir(k,ii)=dsqrt(1/sigma_odl(k,ii)) + distancek(k)=distance(k)**2/(sigma_odlir(k,ii)* + & (distance(k)**2+sigma_odlir(k,ii)**2)) + endif + 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) + & min_odl=distancek(kk) + enddo + +c write (iout,* )"min_odl",min_odl +#ifdef DEBUG + write (iout,*) "ij dij",i,j,dij + write (iout,*) "distance",(distance(k),k=1,constr_homology) + write (iout,*) "distancek",(distancek(k),k=1,constr_homology) + write (iout,* )"min_odl",min_odl +#endif + odleg2=0.0d0 + do k=1,constr_homology +c Nie wiem po co to liczycie jeszcze raz! +c odleg3=-waga_dist(iset)*((distance(i,j,k)**2)/ +c & (2*(sigma_odl(i,j,k))**2)) + if(.not.l_homo(k,ii)) cycle + if (waga_dist.ge.0.0d0) then +c +c For Gaussian-type Urestr +c + godl(k)=dexp(-distancek(k)+min_odl) + odleg2=odleg2+godl(k) +c +c For Lorentzian-type Urestr +c + else + odleg2=odleg2+distancek(k) + endif + +ccc write(iout,779) i,j,k, "odleg2=",odleg2, "odleg3=", odleg3, +ccc & "dEXP(odleg3)=", dEXP(odleg3),"distance(i,j,k)^2=", +ccc & distance(i,j,k)**2, "dist(i+1,j+1)=", dist(i+1,j+1), +ccc & "sigma_odl(i,j,k)=", sigma_odl(i,j,k) + + enddo +c write (iout,*) "godl",(godl(k),k=1,constr_homology) ! exponents +c write (iout,*) "ii i j",ii,i,j," odleg2",odleg2 ! sum of exps +#ifdef DEBUG + write (iout,*) "godl",(godl(k),k=1,constr_homology) ! exponents + write (iout,*) "ii i j",ii,i,j," odleg2",odleg2 ! sum of exps +#endif + if (waga_dist.ge.0.0d0) then +c +c For Gaussian-type Urestr +c + odleg=odleg-dLOG(odleg2/constr_homology)+min_odl +c +c For Lorentzian-type Urestr +c + else + odleg=odleg+odleg2/constr_homology + endif +c +c write (iout,*) "odleg",odleg ! sum of -ln-s +c Gradient +c +c For Gaussian-type Urestr +c + if (waga_dist.ge.0.0d0) sum_godl=odleg2 + sum_sgodl=0.0d0 + do k=1,constr_homology +c godl=dexp(((-(distance(i,j,k)**2)/(2*(sigma_odl(i,j,k))**2)) +c & *waga_dist)+min_odl +c sgodl=-godl(k)*distance(k)*sigma_odl(k,ii)*waga_dist +c + if(.not.l_homo(k,ii)) cycle + if (waga_dist.ge.0.0d0) then +c For Gaussian-type Urestr +c + sgodl=-godl(k)*distance(k)*sigma_odl(k,ii) ! waga_dist rmvd +c +c For Lorentzian-type Urestr +c + else + sgodl=-2*sigma_odlir(k,ii)*(distance(k)/(distance(k)**2+ + & sigma_odlir(k,ii)**2)**2) + endif + sum_sgodl=sum_sgodl+sgodl + +c sgodl2=sgodl2+sgodl +c write(iout,*) i, j, k, distance(i,j,k), "W GRADIENCIE1" +c write(iout,*) "constr_homology=",constr_homology +c write(iout,*) i, j, k, "TEST K" + enddo + if (waga_dist.ge.0.0d0) then +c +c For Gaussian-type Urestr +c + grad_odl3=waga_homology(iset)*waga_dist + & *sum_sgodl/(sum_godl*dij) +c +c For Lorentzian-type Urestr +c + else +c Original grad expr modified by analogy w Gaussian-type Urestr grad +c grad_odl3=-waga_homology(iset)*waga_dist*sum_sgodl + grad_odl3=-waga_homology(iset)*waga_dist* + & sum_sgodl/(constr_homology*dij) + endif +c +c grad_odl3=sum_sgodl/(sum_godl*dij) + + +c write(iout,*) i, j, k, distance(i,j,k), "W GRADIENCIE2" +c write(iout,*) (distance(i,j,k)**2), (2*(sigma_odl(i,j,k))**2), +c & (-(distance(i,j,k)**2)/(2*(sigma_odl(i,j,k))**2)) + +ccc write(iout,*) godl, sgodl, grad_odl3 + +c grad_odl=grad_odl+grad_odl3 + + do jik=1,3 + ggodl=grad_odl3*(c(jik,i)-c(jik,j)) +ccc write(iout,*) c(jik,i+1), c(jik,j+1), (c(jik,i+1)-c(jik,j+1)) +ccc write(iout,746) "GRAD_ODL_1", i, j, jik, ggodl, +ccc & ghpbc(jik,i+1), ghpbc(jik,j+1) + ghpbc(jik,i)=ghpbc(jik,i)+ggodl + ghpbc(jik,j)=ghpbc(jik,j)-ggodl +ccc write(iout,746) "GRAD_ODL_2", i, j, jik, ggodl, +ccc & ghpbc(jik,i+1), ghpbc(jik,j+1) +c if (i.eq.25.and.j.eq.27) then +c write(iout,*) "jik",jik,"i",i,"j",j +c write(iout,*) "sum_sgodl",sum_sgodl,"sgodl",sgodl +c write(iout,*) "grad_odl3",grad_odl3 +c write(iout,*) "c(",jik,i,")",c(jik,i),"c(",jik,j,")",c(jik,j) +c write(iout,*) "ggodl",ggodl +c write(iout,*) "ghpbc(",jik,i,")", +c & ghpbc(jik,i),"ghpbc(",jik,j,")", +c & ghpbc(jik,j) +c endif + enddo +ccc write(iout,778)"TEST: odleg2=", odleg2, "DLOG(odleg2)=", +ccc & dLOG(odleg2),"-odleg=", -odleg + + enddo ! ii-loop for dist +#ifdef DEBUG + write(iout,*) "------- dist restrs end -------" +c if (waga_angle.eq.1.0d0 .or. waga_theta.eq.1.0d0 .or. +c & waga_d.eq.1.0d0) call sum_gradient +#endif +c Pseudo-energy and gradient from dihedral-angle restraints from +c homology templates +c write (iout,*) "End of distance loop" +c call flush(iout) + kat=0.0d0 +c write (iout,*) idihconstr_start_homo,idihconstr_end_homo +#ifdef DEBUG + write(iout,*) "------- dih restrs start -------" + do i=idihconstr_start_homo,idihconstr_end_homo + write (iout,*) "gloc_init(",i,icg,")",gloc(i,icg) + enddo +#endif + do i=idihconstr_start_homo,idihconstr_end_homo + kat2=0.0d0 +c betai=beta(i,i+1,i+2,i+3) + betai = phi(i) +c write (iout,*) "betai =",betai + do k=1,constr_homology + dih_diff(k)=pinorm(dih(k,i)-betai) +cd write (iout,'(a8,2i4,2f15.8)') "dih_diff",i,k,dih_diff(k) +cd & ,sigma_dih(k,i) +c if (dih_diff(i,k).gt.3.14159) dih_diff(i,k)= +c & -(6.28318-dih_diff(i,k)) +c if (dih_diff(i,k).lt.-3.14159) dih_diff(i,k)= +c & 6.28318+dih_diff(i,k) + + kat3=-0.5d0*dih_diff(k)**2*sigma_dih(k,i) ! waga_angle rmvd from Gaussian argument +c kat3=-0.5d0*waga_angle*dih_diff(k)**2*sigma_dih(k,i) + gdih(k)=dexp(kat3) + kat2=kat2+gdih(k) +c write(iout,*) "kat2=", kat2, "exp(kat3)=", exp(kat3) +c write(*,*)"" + enddo +c write (iout,*) "gdih",(gdih(k),k=1,constr_homology) ! exps +c write (iout,*) "i",i," betai",betai," kat2",kat2 ! sum of exps +#ifdef DEBUG + write (iout,*) "i",i," betai",betai," kat2",kat2 + write (iout,*) "gdih",(gdih(k),k=1,constr_homology) +#endif + if (kat2.le.1.0d-14) cycle + kat=kat-dLOG(kat2/constr_homology) +c write (iout,*) "kat",kat ! sum of -ln-s + +ccc write(iout,778)"TEST: kat2=", kat2, "DLOG(kat2)=", +ccc & dLOG(kat2), "-kat=", -kat + +c ---------------------------------------------------------------------- +c Gradient +c ---------------------------------------------------------------------- + + sum_gdih=kat2 + sum_sgdih=0.0d0 + do k=1,constr_homology + sgdih=-gdih(k)*dih_diff(k)*sigma_dih(k,i) ! waga_angle rmvd +c sgdih=-gdih(k)*dih_diff(k)*sigma_dih(k,i)*waga_angle + sum_sgdih=sum_sgdih+sgdih + enddo +c grad_dih3=sum_sgdih/sum_gdih + grad_dih3=waga_homology(iset)*waga_angle*sum_sgdih/sum_gdih + +c write(iout,*)i,k,gdih,sgdih,beta(i+1,i+2,i+3,i+4),grad_dih3 +ccc write(iout,747) "GRAD_KAT_1", i, nphi, icg, grad_dih3, +ccc & gloc(nphi+i-3,icg) + gloc(i-3,icg)=gloc(i-3,icg)+grad_dih3 +c if (i.eq.25) then +c write(iout,*) "i",i,"icg",icg,"gloc(",i,icg,")",gloc(i,icg) +c endif +ccc write(iout,747) "GRAD_KAT_2", i, nphi, icg, grad_dih3, +ccc & gloc(nphi+i-3,icg) + + enddo ! i-loop for dih +#ifdef DEBUG + write(iout,*) "------- dih restrs end -------" +#endif + +c Pseudo-energy and gradient for theta angle restraints from +c homology templates +c FP 01/15 - inserted from econstr_local_test.F, loop structure +c adapted + +c +c For constr_homology reference structures (FP) +c +c Uconst_back_tot=0.0d0 + Eval=0.0d0 + Erot=0.0d0 +c Econstr_back legacy + do i=1,nres +c do i=ithet_start,ithet_end + dutheta(i)=0.0d0 +c enddo +c do i=loc_start,loc_end + do j=1,3 + duscdiff(j,i)=0.0d0 + duscdiffx(j,i)=0.0d0 + enddo + enddo +c +c do iref=1,nref +c write (iout,*) "ithet_start =",ithet_start,"ithet_end =",ithet_end +c write (iout,*) "waga_theta",waga_theta + if (waga_theta.gt.0.0d0) then +#ifdef DEBUG + write (iout,*) "usampl",usampl + write(iout,*) "------- theta restrs start -------" +c do i=ithet_start,ithet_end +c write (iout,*) "gloc_init(",nphi+i,icg,")",gloc(nphi+i,icg) +c enddo +#endif +c write (iout,*) "maxres",maxres,"nres",nres + + do i=ithet_start,ithet_end +c +c do i=1,nfrag_back +c ii = ifrag_back(2,i,iset)-ifrag_back(1,i,iset) +c +c Deviation of theta angles wrt constr_homology ref structures +c + utheta_i=0.0d0 ! argument of Gaussian for single k + gutheta_i=0.0d0 ! Sum of Gaussians over constr_homology ref structures +c do j=ifrag_back(1,i,iset)+2,ifrag_back(2,i,iset) ! original loop +c over residues in a fragment +c write (iout,*) "theta(",i,")=",theta(i) + do k=1,constr_homology +c +c dtheta_i=theta(j)-thetaref(j,iref) +c dtheta_i=thetaref(k,i)-theta(i) ! original form without indexing + theta_diff(k)=thetatpl(k,i)-theta(i) +cd write (iout,'(a8,2i4,2f15.8)') "theta_diff",i,k,theta_diff(k) +cd & ,sigma_theta(k,i) + +c + utheta_i=-0.5d0*theta_diff(k)**2*sigma_theta(k,i) ! waga_theta rmvd from Gaussian argument +c utheta_i=-0.5d0*waga_theta*theta_diff(k)**2*sigma_theta(k,i) ! waga_theta? + gtheta(k)=dexp(utheta_i) ! + min_utheta_i? + gutheta_i=gutheta_i+dexp(utheta_i) ! Sum of Gaussians (pk) +c Gradient for single Gaussian restraint in subr Econstr_back +c dutheta(j-2)=dutheta(j-2)+wfrag_back(1,i,iset)*dtheta_i/(ii-1) +c + enddo +c write (iout,*) "gtheta",(gtheta(k),k=1,constr_homology) ! exps +c write (iout,*) "i",i," gutheta_i",gutheta_i ! sum of exps + +c +c Gradient for multiple Gaussian restraint + sum_gtheta=gutheta_i + sum_sgtheta=0.0d0 + do k=1,constr_homology +c New generalized expr for multiple Gaussian from Econstr_back + sgtheta=-gtheta(k)*theta_diff(k)*sigma_theta(k,i) ! waga_theta rmvd +c +c sgtheta=-gtheta(k)*theta_diff(k)*sigma_theta(k,i)*waga_theta ! right functional form? + sum_sgtheta=sum_sgtheta+sgtheta ! cum variable + enddo +c Final value of gradient using same var as in Econstr_back + gloc(nphi+i-2,icg)=gloc(nphi+i-2,icg) + & +sum_sgtheta/sum_gtheta*waga_theta + & *waga_homology(iset) +c dutheta(i-2)=sum_sgtheta/sum_gtheta*waga_theta +c & *waga_homology(iset) +c dutheta(i)=sum_sgtheta/sum_gtheta +c +c Uconst_back=Uconst_back+waga_theta*utheta(i) ! waga_theta added as weight + Eval=Eval-dLOG(gutheta_i/constr_homology) +c write (iout,*) "utheta(",i,")=",utheta(i) ! -ln of sum of exps +c write (iout,*) "Uconst_back",Uconst_back ! sum of -ln-s +c Uconst_back=Uconst_back+utheta(i) + enddo ! (i-loop for theta) +#ifdef DEBUG + write(iout,*) "------- theta restrs end -------" +#endif + endif +c +c Deviation of local SC geometry +c +c Separation of two i-loops (instructed by AL - 11/3/2014) +c +c write (iout,*) "loc_start =",loc_start,"loc_end =",loc_end +c write (iout,*) "waga_d",waga_d + +#ifdef DEBUG + write(iout,*) "------- SC restrs start -------" + write (iout,*) "Initial duscdiff,duscdiffx" + do i=loc_start,loc_end + write (iout,*) i,(duscdiff(jik,i),jik=1,3), + & (duscdiffx(jik,i),jik=1,3) + enddo +#endif + do i=loc_start,loc_end + usc_diff_i=0.0d0 ! argument of Gaussian for single k + guscdiff(i)=0.0d0 ! Sum of Gaussians over constr_homology ref structures +c do j=ifrag_back(1,i,iset)+1,ifrag_back(2,i,iset)-1 ! Econstr_back legacy +c write(iout,*) "xxtab, yytab, zztab" +c write(iout,'(i5,3f8.2)') i,xxtab(i),yytab(i),zztab(i) + do k=1,constr_homology +c + dxx=-xxtpl(k,i)+xxtab(i) ! Diff b/w x component of ith SC vector in model and kth ref str? +c Original sign inverted for calc of gradients (s. Econstr_back) + dyy=-yytpl(k,i)+yytab(i) ! ibid y + dzz=-zztpl(k,i)+zztab(i) ! ibid z +c write(iout,*) "dxx, dyy, dzz" +cd write(iout,'(2i5,4f8.2)') k,i,dxx,dyy,dzz,sigma_d(k,i) +c + usc_diff_i=-0.5d0*(dxx**2+dyy**2+dzz**2)*sigma_d(k,i) ! waga_d rmvd from Gaussian argument +c usc_diff(i)=-0.5d0*waga_d*(dxx**2+dyy**2+dzz**2)*sigma_d(k,i) ! waga_d? +c uscdiffk(k)=usc_diff(i) + guscdiff2(k)=dexp(usc_diff_i) ! without min_scdiff + guscdiff(i)=guscdiff(i)+dexp(usc_diff_i) !Sum of Gaussians (pk) +c write (iout,'(i5,6f10.5)') j,xxtab(j),yytab(j),zztab(j), +c & xxref(j),yyref(j),zzref(j) + enddo +c +c Gradient +c +c Generalized expression for multiple Gaussian acc to that for a single +c Gaussian in Econstr_back as instructed by AL (FP - 03/11/2014) +c +c Original implementation +c sum_guscdiff=guscdiff(i) +c +c sum_sguscdiff=0.0d0 +c do k=1,constr_homology +c sguscdiff=-guscdiff2(k)*dscdiff(k)*sigma_d(k,i)*waga_d !waga_d? +c sguscdiff=-guscdiff3(k)*dscdiff(k)*sigma_d(k,i)*waga_d ! w min_uscdiff +c sum_sguscdiff=sum_sguscdiff+sguscdiff +c enddo +c +c Implementation of new expressions for gradient (Jan. 2015) +c +c grad_uscdiff=sum_sguscdiff/(sum_guscdiff*dtab) !? + do k=1,constr_homology +c +c New calculation of dxx, dyy, and dzz corrected by AL (07/11), was missing and wrong +c before. Now the drivatives should be correct +c + dxx=-xxtpl(k,i)+xxtab(i) ! Diff b/w x component of ith SC vector in model and kth ref str? +c Original sign inverted for calc of gradients (s. Econstr_back) + dyy=-yytpl(k,i)+yytab(i) ! ibid y + dzz=-zztpl(k,i)+zztab(i) ! ibid z +c +c New implementation +c + sum_guscdiff=guscdiff2(k)*!(dsqrt(dxx*dxx+dyy*dyy+dzz*dzz))* -> wrong! + & sigma_d(k,i) ! for the grad wrt r' +c sum_sguscdiff=sum_sguscdiff+sum_guscdiff +c +c +c New implementation + sum_guscdiff = waga_homology(iset)*waga_d*sum_guscdiff + do jik=1,3 + duscdiff(jik,i-1)=duscdiff(jik,i-1)+ + & sum_guscdiff*(dXX_C1tab(jik,i)*dxx+ + & dYY_C1tab(jik,i)*dyy+dZZ_C1tab(jik,i)*dzz)/guscdiff(i) + duscdiff(jik,i)=duscdiff(jik,i)+ + & sum_guscdiff*(dXX_Ctab(jik,i)*dxx+ + & dYY_Ctab(jik,i)*dyy+dZZ_Ctab(jik,i)*dzz)/guscdiff(i) + duscdiffx(jik,i)=duscdiffx(jik,i)+ + & sum_guscdiff*(dXX_XYZtab(jik,i)*dxx+ + & dYY_XYZtab(jik,i)*dyy+dZZ_XYZtab(jik,i)*dzz)/guscdiff(i) +c +#ifdef DEBUG + write(iout,*) "jik",jik,"i",i + write(iout,*) "dxx, dyy, dzz" + write(iout,'(2i5,3f8.2)') k,i,dxx,dyy,dzz + write(iout,*) "guscdiff2(",k,")",guscdiff2(k) +c write(iout,*) "sum_sguscdiff",sum_sguscdiff +cc write(iout,*) "dXX_Ctab(",jik,i,")",dXX_Ctab(jik,i) +c write(iout,*) "dYY_Ctab(",jik,i,")",dYY_Ctab(jik,i) +c write(iout,*) "dZZ_Ctab(",jik,i,")",dZZ_Ctab(jik,i) +c write(iout,*) "dXX_C1tab(",jik,i,")",dXX_C1tab(jik,i) +c write(iout,*) "dYY_C1tab(",jik,i,")",dYY_C1tab(jik,i) +c write(iout,*) "dZZ_C1tab(",jik,i,")",dZZ_C1tab(jik,i) +c write(iout,*) "dXX_XYZtab(",jik,i,")",dXX_XYZtab(jik,i) +c write(iout,*) "dYY_XYZtab(",jik,i,")",dYY_XYZtab(jik,i) +c write(iout,*) "dZZ_XYZtab(",jik,i,")",dZZ_XYZtab(jik,i) +c write(iout,*) "duscdiff(",jik,i-1,")",duscdiff(jik,i-1) +c write(iout,*) "duscdiff(",jik,i,")",duscdiff(jik,i) +c write(iout,*) "duscdiffx(",jik,i,")",duscdiffx(jik,i) +c endif +#endif + enddo + enddo +c +c uscdiff(i)=-dLOG(guscdiff(i)/(ii-1)) ! Weighting by (ii-1) required? +c usc_diff(i)=-dLOG(guscdiff(i)/constr_homology) ! + min_uscdiff ? +c +c write (iout,*) i," uscdiff",uscdiff(i) +c +c Put together deviations from local geometry + +c Uconst_back=Uconst_back+wfrag_back(1,i,iset)*utheta(i)+ +c & wfrag_back(3,i,iset)*uscdiff(i) + Erot=Erot-dLOG(guscdiff(i)/constr_homology) +c write (iout,*) "usc_diff(",i,")=",usc_diff(i) ! -ln of sum of exps +c write (iout,*) "Uconst_back",Uconst_back ! cum sum of -ln-s +c Uconst_back=Uconst_back+usc_diff(i) +c +c Gradient of multiple Gaussian restraint (FP - 04/11/2014 - right?) +c +c New implment: multiplied by sum_sguscdiff +c + + enddo ! (i-loop for dscdiff) + +c endif + +#ifdef DEBUG + write(iout,*) "------- SC restrs end -------" + write (iout,*) "------ After SC loop in e_modeller ------" + do i=loc_start,loc_end + write (iout,*) "i",i," gradc",(gradc(j,i,icg),j=1,3) + write (iout,*) "i",i," gradx",(gradx(j,i,icg),j=1,3) + enddo + if (waga_theta.eq.1.0d0) then + write (iout,*) "in e_modeller after SC restr end: dutheta" + do i=ithet_start,ithet_end + write (iout,*) i,dutheta(i) + enddo + endif + if (waga_d.eq.1.0d0) then + write (iout,*) "e_modeller after SC loop: duscdiff/x" + do i=1,nres + write (iout,*) i,(duscdiff(j,i),j=1,3) + write (iout,*) i,(duscdiffx(j,i),j=1,3) + enddo + endif +#endif + +c Total energy from homology restraints +#ifdef DEBUG + write (iout,*) "odleg",odleg," kat",kat +#endif +c +c Addition of energy of theta angle and SC local geom over constr_homologs ref strs +c +c ehomology_constr=odleg+kat +c +c For Lorentzian-type Urestr +c + + if (waga_dist.ge.0.0d0) then +c +c For Gaussian-type Urestr +c + ehomology_constr=(waga_dist*odleg+waga_angle*kat+ + & waga_theta*Eval+waga_d*Erot)*waga_homology(iset) +c write (iout,*) "ehomology_constr=",ehomology_constr + else +c +c For Lorentzian-type Urestr +c + ehomology_constr=(-waga_dist*odleg+waga_angle*kat+ + & waga_theta*Eval+waga_d*Erot)*waga_homology(iset) +c write (iout,*) "ehomology_constr=",ehomology_constr + endif +#ifdef DEBUG + write (iout,*) "odleg",waga_dist,odleg," kat",waga_angle,kat, + & "Eval",waga_theta,eval, + & "Erot",waga_d,Erot + write (iout,*) "ehomology_constr",ehomology_constr +#endif + return +c +c FP 01/15 end +c + 748 format(a8,f12.3,a6,f12.3,a7,f12.3) + 747 format(a12,i4,i4,i4,f8.3,f8.3) + 746 format(a12,i4,i4,i4,f8.3,f8.3,f8.3) + 778 format(a7,1X,f10.3,1X,a4,1X,f10.3,1X,a5,1X,f10.3) + 779 format(i3,1X,i3,1X,i2,1X,a7,1X,f7.3,1X,a7,1X,f7.3,1X,a13,1X, + & f7.3,1X,a17,1X,f9.3,1X,a10,1X,f8.3,1X,a10,1X,f8.3) + end + +c------------------------------------------------------------------------------ subroutine etor_d(etors_d) C 6/23/01 Compute double torsional energy implicit real*8 (a-h,o-z) @@ -6906,12 +7590,12 @@ c write (iout,*) "EBACK_SC_COR",itau_start,itau_end esccor=0.0D0 do i=itau_start,itau_end if ((itype(i-2).eq.ntyp1).or.(itype(i-1).eq.ntyp1)) cycle - esccor_ii=0.0D0 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 @@ -6935,9 +7619,12 @@ c 3 = SC...Ca...Ca...SCi 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 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) @@ -10362,3 +11049,398 @@ C print *,'AFM',Eafmforce,totTafm*velAFMconst,dist return end +c---------------------------------------------------------------------------- + double precision function sscale2(r,r_cut,r0,rlamb) + implicit none + double precision r,gamm,r_cut,r0,rlamb,rr + rr = dabs(r-r0) +c write (2,*) "r",r," r_cut",r_cut," r0",r0," rlamb",rlamb +c write (2,*) "rr",rr + if(rr.lt.r_cut-rlamb) then + sscale2=1.0d0 + else if(rr.le.r_cut.and.rr.ge.r_cut-rlamb) then + gamm=(rr-(r_cut-rlamb))/rlamb + sscale2=1.0d0+gamm*gamm*(2*gamm-3.0d0) + else + sscale2=0d0 + endif + return + end +C----------------------------------------------------------------------- + double precision function sscalgrad2(r,r_cut,r0,rlamb) + implicit none + double precision r,gamm,r_cut,r0,rlamb,rr + rr = dabs(r-r0) + if(rr.lt.r_cut-rlamb) then + sscalgrad2=0.0d0 + else if(rr.le.r_cut.and.rr.ge.r_cut-rlamb) then + gamm=(rr-(r_cut-rlamb))/rlamb + sscalgrad2=gamm*(6*gamm-6.0d0)/rlamb + if (r.lt.r0) sscalgrad2=-sscalgrad2 + else + sscalgrad2=0.0d0 + endif + return + end +c---------------------------------------------------------------------------- + subroutine e_saxs(Esaxs_constr) + implicit none + include 'DIMENSIONS' +#ifdef MPI + include "mpif.h" + include "COMMON.SETUP" + integer IERR +#endif + include 'COMMON.SBRIDGE' + include 'COMMON.CHAIN' + include 'COMMON.GEO' + include 'COMMON.DERIV' + include 'COMMON.LOCAL' + include 'COMMON.INTERACT' + include 'COMMON.VAR' + include 'COMMON.IOUNITS' + include 'COMMON.MD' + include 'COMMON.CONTROL' + include 'COMMON.NAMES' + include 'COMMON.TIME1' + include 'COMMON.FFIELD' +c + double precision Esaxs_constr + integer i,iint,j,k,l + double precision PgradC(maxSAXS,3,maxres), + & PgradX(maxSAXS,3,maxres),Pcalc(maxSAXS) +#ifdef MPI + double precision PgradC_(maxSAXS,3,maxres), + & PgradX_(maxSAXS,3,maxres),Pcalc_(maxSAXS) +#endif + double precision dk,dijCACA,dijCASC,dijSCCA,dijSCSC, + & sigma2CACA,sigma2CASC,sigma2SCCA,sigma2SCSC,expCACA,expCASC, + & expSCCA,expSCSC,CASCgrad,SCCAgrad,SCSCgrad,aux,auxC,auxC1, + & auxX,auxX1,CACAgrad,Cnorm + double precision sss2,ssgrad2,rrr,sscalgrad2,sscale2 + double precision dist + external dist +c SAXS restraint penalty function +#ifdef DEBUG + write(iout,*) "------- SAXS penalty function start -------" + write (iout,*) "nsaxs",nsaxs + write (iout,*) "Esaxs: iatsc_s",iatsc_s," iatsc_e",iatsc_e + write (iout,*) "Psaxs" + do i=1,nsaxs + write (iout,'(i5,e15.5)') i, Psaxs(i) + enddo +#endif + Esaxs_constr = 0.0d0 + do k=1,nsaxs + Pcalc(k)=0.0d0 + do j=1,nres + do l=1,3 + PgradC(k,l,j)=0.0d0 + PgradX(k,l,j)=0.0d0 + enddo + enddo + enddo + do i=iatsc_s,iatsc_e + if (itype(i).eq.ntyp1) cycle + do iint=1,nint_gr(i) + do j=istart(i,iint),iend(i,iint) + if (itype(j).eq.ntyp1) cycle +#ifdef ALLSAXS + dijCACA=dist(i,j) + dijCASC=dist(i,j+nres) + dijSCCA=dist(i+nres,j) + dijSCSC=dist(i+nres,j+nres) + sigma2CACA=2.0d0/(pstok**2) + sigma2CASC=4.0d0/(pstok**2+restok(itype(j))**2) + sigma2SCCA=4.0d0/(pstok**2+restok(itype(i))**2) + sigma2SCSC=4.0d0/(restok(itype(j))**2+restok(itype(i))**2) + do k=1,nsaxs + dk = distsaxs(k) + expCACA = dexp(-0.5d0*sigma2CACA*(dijCACA-dk)**2) + if (itype(j).ne.10) then + expCASC = dexp(-0.5d0*sigma2CASC*(dijCASC-dk)**2) + else + endif + expCASC = 0.0d0 + if (itype(i).ne.10) then + expSCCA = dexp(-0.5d0*sigma2SCCA*(dijSCCA-dk)**2) + else + expSCCA = 0.0d0 + endif + if (itype(i).ne.10 .and. itype(j).ne.10) then + expSCSC = dexp(-0.5d0*sigma2SCSC*(dijSCSC-dk)**2) + else + expSCSC = 0.0d0 + endif + Pcalc(k) = Pcalc(k)+expCACA+expCASC+expSCCA+expSCSC +#ifdef DEBUG + write(iout,*) "i j k Pcalc",i,j,Pcalc(k) +#endif + CACAgrad = sigma2CACA*(dijCACA-dk)*expCACA + CASCgrad = sigma2CASC*(dijCASC-dk)*expCASC + SCCAgrad = sigma2SCCA*(dijSCCA-dk)*expSCCA + SCSCgrad = sigma2SCSC*(dijSCSC-dk)*expSCSC + do l=1,3 +c CA CA + aux = CACAgrad*(C(l,j)-C(l,i))/dijCACA + PgradC(k,l,i) = PgradC(k,l,i)-aux + PgradC(k,l,j) = PgradC(k,l,j)+aux +c CA SC + if (itype(j).ne.10) then + aux = CASCgrad*(C(l,j+nres)-C(l,i))/dijCASC + PgradC(k,l,i) = PgradC(k,l,i)-aux + PgradC(k,l,j) = PgradC(k,l,j)+aux + PgradX(k,l,j) = PgradX(k,l,j)+aux + endif +c SC CA + if (itype(i).ne.10) then + aux = SCCAgrad*(C(l,j)-C(l,i+nres))/dijSCCA + PgradX(k,l,i) = PgradX(k,l,i)-aux + PgradC(k,l,i) = PgradC(k,l,i)-aux + PgradC(k,l,j) = PgradC(k,l,j)+aux + endif +c SC SC + if (itype(i).ne.10 .and. itype(j).ne.10) then + aux = SCSCgrad*(C(l,j+nres)-C(l,i+nres))/dijSCSC + PgradC(k,l,i) = PgradC(k,l,i)-aux + PgradC(k,l,j) = PgradC(k,l,j)+aux + PgradX(k,l,i) = PgradX(k,l,i)-aux + PgradX(k,l,j) = PgradX(k,l,j)+aux + endif + enddo ! l + enddo ! k +#else + dijCACA=dist(i,j) + sigma2CACA=scal_rad**2*0.25d0/ + & (restok(itype(j))**2+restok(itype(i))**2) + rrr = 2.0d0/dsqrt(sigma2CACA) + do k=1,nsaxs + dk = distsaxs(k) +c write (2,*) "ijk",i,j,k + sss2 = sscale2(dijCACA,rrr,dk,0.3d0) + if (sss2.ne.0.0d0) then + ssgrad2 = sscalgrad2(dijCACA,rrr,dk,0.3d0) + if (energy_dec) write(iout,'(a4,3i5,4f10.4)') + & 'saxs',i,j,k,dijCACA,rrr,dk,sss2 + expCACA = dexp(-0.5d0*sigma2CACA*(dijCACA-dk)**2)*sss2 + Pcalc(k) = Pcalc(k)+expCACA +#ifdef DEBUG + write(iout,*) "i j k Pcalc",i,j,Pcalc(k) +#endif + CACAgrad = sigma2CACA*(dijCACA-dk)*expCACA+ + & ssgrad2*expCACA/sss2 + do l=1,3 +c CA CA + aux = CACAgrad*(C(l,j)-C(l,i))/dijCACA + PgradC(k,l,i) = PgradC(k,l,i)-aux + PgradC(k,l,j) = PgradC(k,l,j)+aux + enddo ! l + endif ! sss + enddo ! k +#endif + enddo ! j + enddo ! iint + enddo ! i +#ifdef MPI + if (nfgtasks.gt.1) then + call MPI_Reduce(Pcalc(1),Pcalc_(1),nsaxs,MPI_DOUBLE_PRECISION, + & MPI_SUM,king,FG_COMM,IERR) + if (fg_rank.eq.king) then + do k=1,nsaxs + Pcalc(k) = Pcalc_(k) + enddo + endif + call MPI_Reduce(PgradC(k,1,1),PgradC_(k,1,1),3*maxsaxs*nres, + & MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR) + if (fg_rank.eq.king) then + do i=1,nres + do l=1,3 + do k=1,nsaxs + PgradC(k,l,i) = PgradC_(k,l,i) + enddo + enddo + enddo + endif +#ifdef ALLSAXS + call MPI_Reduce(PgradX(k,1,1),PgradX_(k,1,1),3*maxsaxs*nres, + & MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR) + if (fg_rank.eq.king) then + do i=1,nres + do l=1,3 + do k=1,nsaxs + PgradX(k,l,i) = PgradX_(k,l,i) + enddo + enddo + enddo + endif +#endif + endif +#endif +#ifdef MPI + if (fg_rank.eq.king) then +#endif + Cnorm = 0.0d0 + do k=1,nsaxs + Cnorm = Cnorm + Pcalc(k) + enddo + Esaxs_constr = dlog(Cnorm) + do k=1,nsaxs + Esaxs_constr = Esaxs_constr - Psaxs(k)*dlog(Pcalc(k)) +#ifdef DEBUG + write (iout,*) "k",k," Esaxs_constr",Esaxs_constr +#endif + enddo +#ifdef DEBUG + write (iout,*) "Cnorm",Cnorm," Esaxs_constr",Esaxs_constr +#endif + do i=nnt,nct + do l=1,3 + auxC=0.0d0 + auxC1=0.0d0 + auxX=0.0d0 + auxX1=0.d0 + do k=1,nsaxs + auxC = auxC +Psaxs(k)*PgradC(k,l,i)/Pcalc(k) + auxC1 = auxC1+PgradC(k,l,i) +#ifdef ALLSAXS + auxX = auxX +Psaxs(k)*PgradX(k,l,i)/Pcalc(k) + auxX1 = auxX1+PgradX(k,l,i) +#endif + enddo + gsaxsC(l,i) = auxC - auxC1/Cnorm +#ifdef ALLSAXS + gsaxsX(l,i) = auxX - auxX1/Cnorm +#endif +c write (iout,*) "l i",l,i," gradC",wsaxs*(auxC - auxC1/Cnorm), +c * " gradX",wsaxs*(auxX - auxX1/Cnorm) + enddo + enddo +#ifdef MPI + endif +#endif + return + end +c---------------------------------------------------------------------------- + subroutine e_saxsC(Esaxs_constr) + implicit none + include 'DIMENSIONS' +#ifdef MPI + include "mpif.h" + include "COMMON.SETUP" + integer IERR +#endif + include 'COMMON.SBRIDGE' + include 'COMMON.CHAIN' + include 'COMMON.GEO' + include 'COMMON.DERIV' + include 'COMMON.LOCAL' + include 'COMMON.INTERACT' + include 'COMMON.VAR' + include 'COMMON.IOUNITS' + include 'COMMON.MD' + include 'COMMON.CONTROL' + include 'COMMON.NAMES' + include 'COMMON.TIME1' + include 'COMMON.FFIELD' +c + double precision Esaxs_constr + integer i,iint,j,k,l + double precision PgradC(3,maxres),PgradX(3,maxres),Pcalc,logPtot +#ifdef MPI + double precision gsaxsc_(3,maxres),gsaxsx_(3,maxres),logPtot_ +#endif + double precision dk,dijCASPH,dijSCSPH, + & sigma2CA,sigma2SC,expCASPH,expSCSPH, + & CASPHgrad,SCSPHgrad,aux,auxC,auxC1, + & auxX,auxX1,Cnorm +c SAXS restraint penalty function +#ifdef DEBUG + write(iout,*) "------- SAXS penalty function start -------" + write (iout,*) "nsaxs",nsaxs + + do i=nnt,nct + print *,MyRank,"C",i,(C(j,i),j=1,3) + enddo + do i=nnt,nct + print *,MyRank,"CSaxs",i,(Csaxs(j,i),j=1,3) + enddo +#endif + Esaxs_constr = 0.0d0 + logPtot=0.0d0 + do j=isaxs_start,isaxs_end + Pcalc=0.0d0 + do i=1,nres + do l=1,3 + PgradC(l,i)=0.0d0 + PgradX(l,i)=0.0d0 + enddo + enddo + do i=nnt,nct + if (itype(i).eq.ntyp1) cycle + dijCASPH=0.0d0 + dijSCSPH=0.0d0 + do l=1,3 + dijCASPH=dijCASPH+(C(l,i)-Csaxs(l,j))**2 + enddo + if (itype(i).ne.10) then + do l=1,3 + dijSCSPH=dijSCSPH+(C(l,i+nres)-Csaxs(l,j))**2 + enddo + endif + sigma2CA=2.0d0/pstok**2 + sigma2SC=4.0d0/restok(itype(i))**2 + expCASPH = dexp(-0.5d0*sigma2CA*dijCASPH) + expSCSPH = dexp(-0.5d0*sigma2SC*dijSCSPH) + Pcalc = Pcalc+expCASPH+expSCSPH +#ifdef DEBUG + write(*,*) "processor i j Pcalc", + & MyRank,i,j,dijCASPH,dijSCSPH, Pcalc +#endif + CASPHgrad = sigma2CA*expCASPH + SCSPHgrad = sigma2SC*expSCSPH + do l=1,3 + aux = (C(l,i+nres)-Csaxs(l,j))*SCSPHgrad + PgradX(l,i) = PgradX(l,i) + aux + PgradC(l,i) = PgradC(l,i)+(C(l,i)-Csaxs(l,j))*CASPHgrad+aux + enddo ! l + enddo ! i + do i=nnt,nct + do l=1,3 + gsaxsc(l,i)=gsaxsc(l,i)+PgradC(l,i)/Pcalc + gsaxsx(l,i)=gsaxsx(l,i)+PgradX(l,i)/Pcalc + enddo + enddo + logPtot = logPtot - dlog(Pcalc) +c print *,"me",me,MyRank," j",j," logPcalc",-dlog(Pcalc), +c & " logPtot",logPtot + enddo ! j +#ifdef MPI + if (nfgtasks.gt.1) then +c write (iout,*) "logPtot before reduction",logPtot + call MPI_Reduce(logPtot,logPtot_,1,MPI_DOUBLE_PRECISION, + & MPI_SUM,king,FG_COMM,IERR) + logPtot = logPtot_ +c write (iout,*) "logPtot after reduction",logPtot + call MPI_Reduce(gsaxsC(1,1),gsaxsC_(1,1),3*nres, + & MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR) + if (fg_rank.eq.king) then + do i=1,nres + do l=1,3 + gsaxsC(l,i) = gsaxsC_(l,i) + enddo + enddo + endif + call MPI_Reduce(gsaxsX(1,1),gsaxsX_(1,1),3*nres, + & MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR) + if (fg_rank.eq.king) then + do i=1,nres + do l=1,3 + gsaxsX(l,i) = gsaxsX_(l,i) + enddo + enddo + endif + endif +#endif + Esaxs_constr = logPtot + return + end