1 subroutine gen_rand_conf(nstart,*)
2 C Generate random conformation or chain cut and regrowth.
3 implicit real*8 (a-h,o-z)
8 include 'COMMON.INTERACT'
9 include 'COMMON.IOUNITS'
12 include 'COMMON.CONTROL'
13 logical overlap,back,fail
14 cd print *,' CG Processor',me,' maxgen=',maxgen
16 cd write (iout,*) 'Gen_Rand_conf: nstart=',nstart
19 phi(4)=gen_phi(4,iabs(itype(2)),abs(itype(3)))
20 c write(iout,*)'phi(4)=',rad2deg*phi(4)
21 if (nstart.lt.3) theta(3)=gen_theta(itype(2),pi,phi(4))
22 c write(iout,*)'theta(3)=',rad2deg*theta(3)
26 do while (fail.and.nsi.le.maxsi)
27 call gen_side(it1,theta(3),alph(2),omeg(2),fail)
30 if (nsi.gt.maxsi) return1
45 do while (i.le.nres .and. niter.lt.maxgen)
48 write (iout,'(/80(1h*)/2a/80(1h*))')
49 & 'Generation procedure went down to ',
50 & 'chain beginning. Cannot continue...'
51 write (*,'(/80(1h*)/2a/80(1h*))')
52 & 'Generation procedure went down to ',
53 & 'chain beginning. Cannot continue...'
60 c print *,'Gen_Rand_Conf: i=',i,' it=',it,' it1=',it1,' it2=',it2,
61 c & ' nit=',nit,' niter=',niter,' maxgen=',maxgen
62 phi(i+1)=gen_phi(i+1,it1,it)
64 phi(i)=gen_phi(i+1,it2,it1)
65 c print *,'phi(',i,')=',phi(i)
66 theta(i-1)=gen_theta(it2,phi(i-1),phi(i))
70 do while (fail.and.nsi.le.maxsi)
71 call gen_side(it2,theta(i-1),alph(i-2),omeg(i-2),fail)
74 if (nsi.gt.maxsi) return1
76 call locate_next_res(i-1)
78 theta(i)=gen_theta(it1,phi(i),phi(i+1))
82 do while (fail.and.nsi.le.maxsi)
83 call gen_side(it1,theta(i),alph(i-1),omeg(i-1),fail)
86 if (nsi.gt.maxsi) return1
88 call locate_next_res(i)
89 if (overlap(i-1)) then
90 if (nit.lt.maxnit) then
100 & 'Cannot generate non-overlaping conformation. Increase MAXNIT.'
102 & 'Cannot generate non-overlaping conformation. Increase MAXNIT.'
113 if (niter.ge.maxgen) then
114 write (iout,'(a,2i5)')
115 & 'Too many trials in conformation generation',niter,maxgen
117 & 'Too many trials in conformation generation',niter,maxgen
122 c(j,nres+nres)=c(j,nres)
126 c-------------------------------------------------------------------------
127 logical function overlap(i)
128 implicit real*8 (a-h,o-z)
130 include 'COMMON.CHAIN'
131 include 'COMMON.INTERACT'
132 include 'COMMON.FFIELD'
136 if (iti.gt.ntyp) return
137 C Check for SC-SC overlaps.
138 cd print *,'nnt=',nnt,' nct=',nct
141 if (j.lt.i-1 .or. ipot.ne.4) then
142 rcomp=sigmaii(iti,itj)
147 if (dist(nres+i,nres+j).lt.redfac*rcomp) then
149 c print *,'overlap, SC-SC: i=',i,' j=',j,
150 c & ' dist=',dist(nres+i,nres+j),' rcomp=',
155 C Check for overlaps between the added peptide group and the preceding
159 c(j,maxres2+1)=0.5D0*(c(j,i)+c(j,i+1))
163 cd print *,'overlap, p-Sc: i=',i,' j=',j,
164 cd & ' dist=',dist(nres+j,maxres2+1)
165 if (dist(nres+j,maxres2+1).lt.4.0D0*redfac) then
170 C Check for overlaps between the added side chain and the preceding peptide
174 c(k,maxres2+1)=0.5D0*(c(k,j)+c(k,j+1))
176 cd print *,'overlap, SC-p: i=',i,' j=',j,
177 cd & ' dist=',dist(nres+i,maxres2+1)
178 if (dist(nres+i,maxres2+1).lt.4.0D0*redfac) then
183 C Check for p-p overlaps
185 c(j,maxres2+2)=0.5D0*(c(j,i)+c(j,i+1))
190 c(k,maxres2+2)=0.5D0*(c(k,j)+c(k,j+1))
192 cd print *,'overlap, p-p: i=',i,' j=',j,
193 cd & ' dist=',dist(maxres2+1,maxres2+2)
194 if(iteli.ne.0.and.itelj.ne.0)then
195 if (dist(maxres2+1,maxres2+2).lt.rpp(iteli,itelj)*redfac) then
203 c--------------------------------------------------------------------------
204 double precision function gen_phi(i,it1,it2)
205 implicit real*8 (a-h,o-z)
208 include 'COMMON.BOUNDS'
209 c gen_phi=ran_number(-pi,pi)
210 C 8/13/98 Generate phi using pre-defined boundaries
211 gen_phi=ran_number(phibound(1,i),phibound(2,i))
214 c---------------------------------------------------------------------------
215 double precision function gen_theta(it,gama,gama1)
216 implicit real*8 (a-h,o-z)
218 include 'COMMON.LOCAL'
220 double precision y(2),z(2)
221 double precision theta_max,theta_min
222 c print *,'gen_theta: it=',it
225 if (dabs(gama).gt.dwapi) then
232 if (dabs(gama1).gt.dwapi) then
239 thet_pred_mean=a0thet(it)
241 thet_pred_mean=thet_pred_mean+athet(k,it,1,1)*y(k)
242 & +bthet(k,it,1,1)*z(k)
246 sig=sig*thet_pred_mean+polthet(j,it)
248 sig=0.5D0/(sig*sig+sigc0(it))
250 &0.5D0*((gthet(2,it)-thet_pred_mean)/gthet(3,it))**2)
251 c print '(i5,5(1pe14.4))',it,(gthet(j,it),j=1,3)
252 c print '(5(1pe14.4))',thet_pred_mean,theta0(it),sig,sig0(it),ak
253 theta_temp=binorm(thet_pred_mean,theta0(it),sig,sig0(it),ak)
254 if (theta_temp.lt.theta_min) theta_temp=theta_min
255 if (theta_temp.gt.theta_max) theta_temp=theta_max
257 c print '(a)','Exiting GENTHETA.'
260 c-------------------------------------------------------------------------
261 subroutine gen_side(it,the,al,om,fail)
262 implicit real*8 (a-h,o-z)
265 include 'COMMON.LOCAL'
266 include 'COMMON.SETUP'
267 include 'COMMON.IOUNITS'
268 double precision MaxBoxLen /10.0D0/
269 double precision Ap_inv(3,3),a(3,3),z(3,maxlob),W1(maxlob),
270 & sumW(0:maxlob),y(2),cm(2),eig(2),box(2,2),work(100),detAp(maxlob)
271 double precision eig_limit /1.0D-8/
272 double precision Big /10.0D0/
273 double precision vec(3,3)
274 logical lprint,fail,lcheck
278 if (the.eq.0.0D0 .or. the.eq.pi) then
280 write (*,'(a,i4,a,i3,a,1pe14.5)')
281 & 'CG Processor:',me,' Error in GenSide: it=',it,' theta=',the
283 cd write (iout,'(a,i3,a,1pe14.5)')
284 cd & 'Error in GenSide: it=',it,' theta=',the
293 print '(a,i4,a)','CG Processor:',me,' Enter Gen_Side.'
294 write (iout,'(a,i4,a)') 'Processor:',me,' Enter Gen_Side.'
296 print *,'it=',it,' nlobit=',nlobit,' the=',the,' tant=',tant
297 write (iout,*) 'it=',it,' nlobit=',nlobit,' the=',the,
301 zz1=tant-censc(1,i,it)
304 a(k,l)=gaussc(k,l,i,it)
307 detApi=a(2,2)*a(3,3)-a(2,3)**2
308 Ap_inv(2,2)=a(3,3)/detApi
309 Ap_inv(2,3)=-a(2,3)/detApi
310 Ap_inv(3,2)=Ap_inv(2,3)
311 Ap_inv(3,3)=a(2,2)/detApi
313 write (*,'(/a,i2/)') 'Cluster #',i
314 write (*,'(3(1pe14.5),5x,1pe14.5)')
315 & ((a(l,k),l=1,3),censc(k,i,it),k=1,3)
316 write (iout,'(/a,i2/)') 'Cluster #',i
317 write (iout,'(3(1pe14.5),5x,1pe14.5)')
318 & ((a(l,k),l=1,3),censc(k,i,it),k=1,3)
323 W1i=W1i+a(k,1)*a(l,1)*Ap_inv(k,l)
327 W1(i)=dexp(bsc(i,it)-0.5D0*W1i*zz1*zz1)
328 c if (lprint) write(*,'(a,3(1pe15.5)/)')
329 c & 'detAp, W1, anormi',detApi,W1i,anormi
333 zk=zk+zz1*Ap_inv(k,l)*a(l,1)
337 detAp(i)=dsqrt(detApi)
341 print *,'W1:',(w1(i),i=1,nlobit)
342 print *,'detAp:',(detAp(i),i=1,nlobit)
345 print '(i2,3f10.5)',i,(rad2deg*z(j,i),j=2,3)
347 write (iout,*) 'W1:',(w1(i),i=1,nlobit)
348 write (iout,*) 'detAp:',(detAp(i),i=1,nlobit)
351 write (iout,'(i2,3f10.5)') i,(rad2deg*z(j,i),j=2,3)
355 C Writing the distribution just to check the procedure
361 fac=fac+W1(i)/detAp(i)
363 fac=1.0D0/(2.0D0*fac*pi)
364 cd print *,it,'fac=',fac
373 a(j-1,k-1)=gaussc(j,k,i,it)
385 wykl=wykl+a(j,k)*(y(j)-z(j+1,i))*(y(k)-z(k+1,i))
388 wart=wart+W1(i)*dexp(-0.5D0*wykl)
395 c print *,'y',y(1),y(2),' fac=',fac
397 write (20,'(2f10.3,1pd15.5)') y(1)*rad2deg,y(2)*rad2deg,wart
402 c print *,'it=',it,' sum=',sum*dV,' sum1=',sum1*dV
406 C Calculate the CM of the system
413 sumW(i)=sumW(i-1)+W1(i)
418 cm(1)=cm(1)+z(2,j)*W1(j)
419 cm(2)=cm(2)+W1(j)*(z(3,1)+pinorm(z(3,j)-z(3,1)))
421 cm(1)=cm(1)/sumW(nlobit)
422 cm(2)=cm(2)/sumW(nlobit)
423 if (cm(1).gt.Big .or. cm(1).lt.-Big .or.
424 & cm(2).gt.Big .or. cm(2).lt.-Big) then
425 cd write (iout,'(a)')
426 cd & 'Unexpected error in GenSide - CM coordinates too large.'
427 cd write (iout,'(i5,2(1pe14.5))') it,cm(1),cm(2)
429 cd & 'Unexpected error in GenSide - CM coordinates too large.'
430 cd write (*,'(i5,2(1pe14.5))') it,cm(1),cm(2)
434 cd print *,'CM:',cm(1),cm(2)
436 C Find the largest search distance from CM
442 a(j-1,k-1)=gaussc(j,k,i,it)
446 call f02faf('N','U',2,a,3,eig,work,100,ifail)
448 call djacob(2,3,10000,1.0d-10,a,vec,eig)
452 print *,'*************** CG Processor',me
453 print *,'CM:',cm(1),cm(2)
454 write (iout,*) '*************** CG Processor',me
455 write (iout,*) 'CM:',cm(1),cm(2)
456 print '(A,8f10.5)','Eigenvalues: ',(1.0/dsqrt(eig(k)),k=1,2)
457 write (iout,'(A,8f10.5)')
458 & 'Eigenvalues: ',(1.0/dsqrt(eig(k)),k=1,2)
461 if (eig(1).lt.eig_limit) then
463 & 'From Mult_Norm: Eigenvalues of A are too small.'
465 & 'From Mult_Norm: Eigenvalues of A are too small.'
472 radius=radius+pinorm(z(j+1,i)-cm(j))**2
474 radius=dsqrt(radius)+3.0D0/dsqrt(eig(1))
475 if (radius.gt.radmax) radmax=radius
477 if (radmax.gt.pi) radmax=pi
479 C Determine the boundaries of the search rectangle.
482 print '(a,4(1pe14.4))','W1: ',(W1(i),i=1,nlob(it) )
483 print '(a,4(1pe14.4))','radmax: ',radmax
485 box(1,1)=dmax1(cm(1)-radmax,0.0D0)
486 box(2,1)=dmin1(cm(1)+radmax,pi)
487 box(1,2)=cm(2)-radmax
488 box(2,2)=cm(2)+radmax
491 print *,'CG Processor',me,' Array BOX:'
495 print '(4(1pe14.4))',((box(k,j),k=1,2),j=1,2)
496 print '(a,4(1pe14.4))','sumW: ',(sumW(i),i=0,nlob(it) )
498 write (iout,*)'CG Processor',me,' Array BOX:'
500 write (iout,*)'Array BOX:'
502 write(iout,'(4(1pe14.4))') ((box(k,j),k=1,2),j=1,2)
503 write(iout,'(a,4(1pe14.4))')'sumW: ',(sumW(i),i=0,nlob(it) )
505 if (box(1,2).lt.-MaxBoxLen .or. box(2,2).gt.MaxBoxLen) then
507 write (iout,'(a,i4,a)') 'CG Processor:',me,': bad sampling box.'
508 write (*,'(a,i4,a)') 'CG Processor:',me,': bad sampling box.'
510 c write (iout,'(a)') 'Bad sampling box.'
515 which_lobe=ran_number(0.0D0,sumW(nlobit))
516 c print '(a,1pe14.4)','which_lobe=',which_lobe
518 if (sumW(i-1).le.which_lobe .and. sumW(i).ge.which_lobe) goto 1
521 c print *,'ilob=',ilob,' nlob=',nlob(it)
525 a(i-1,j-1)=gaussc(i,j,ilob,it)
528 cd print '(a,i4,a)','CG Processor',me,' Calling MultNorm1.'
529 call mult_norm1(3,2,a,cm,box,y,fail)
533 cd print *,'al=',al,' om=',om
537 c---------------------------------------------------------------------------
538 double precision function ran_number(x1,x2)
539 C Calculate a random real number from the range (x1,x2).
540 implicit real*8 (a-h,o-z)
542 double precision x1,x2,fctor
543 data fctor /2147483647.0D0/
546 include 'COMMON.SETUP'
547 ran_number=x1+(x2-x1)*prng_next(me)
550 ran_number=x1+(x2-x1)*ix/fctor
554 c--------------------------------------------------------------------------
555 integer function iran_num(n1,n2)
556 C Calculate a random integer number from the range (n1,n2).
557 implicit real*8 (a-h,o-z)
560 real fctor /2147483647.0/
563 include 'COMMON.SETUP'
564 ix=n1+(n2-n1+1)*prng_next(me)
570 ix=n1+(n2-n1+1)*(ix/fctor)
576 c--------------------------------------------------------------------------
577 double precision function binorm(x1,x2,sigma1,sigma2,ak)
578 implicit real*8 (a-h,o-z)
579 c print '(a)','Enter BINORM.'
580 alowb=dmin1(x1-3.0D0*sigma1,x2-3.0D0*sigma2)
581 aupb=dmax1(x1+3.0D0*sigma1,x2+3.0D0*sigma2)
582 seg=sigma1/(sigma1+ak*sigma2)
583 alen=ran_number(0.0D0,1.0D0)
584 if (alen.lt.seg) then
585 binorm=anorm_distr(x1,sigma1,alowb,aupb)
587 binorm=anorm_distr(x2,sigma2,alowb,aupb)
589 c print '(a)','Exiting BINORM.'
592 c-----------------------------------------------------------------------
593 c double precision function anorm_distr(x,sigma,alowb,aupb)
594 c implicit real*8 (a-h,o-z)
595 c print '(a)','Enter ANORM_DISTR.'
596 c 10 y=ran_number(alowb,aupb)
597 c expon=dexp(-0.5D0*((y-x)/sigma)**2)
598 c ran=ran_number(0.0D0,1.0D0)
599 c if (expon.lt.ran) goto 10
601 c print '(a)','Exiting ANORM_DISTR.'
604 c-----------------------------------------------------------------------
605 double precision function anorm_distr(x,sigma,alowb,aupb)
606 implicit real*8 (a-h,o-z)
607 c to make a normally distributed deviate with zero mean and unit variance
610 real fac,gset,rsq,v1,v2,ran1
614 1 v1=2.0d0*ran_number(0.0d0,1.0d0)-1.0d0
615 v2=2.0d0*ran_number(0.0d0,1.0d0)-1.0d0
617 if(rsq.ge.1.d0.or.rsq.eq.0.0d0) goto 1
618 fac=sqrt(-2.0d0*log(rsq)/rsq)
626 anorm_distr=x+gaussdev*sigma
629 c------------------------------------------------------------------------
630 subroutine mult_norm(lda,n,a,x,fail)
632 C Generate the vector X whose elements obey the multiple-normal distribution
633 C from exp(-0.5*X'AX). LDA is the leading dimension of the moment matrix A,
634 C n is the dimension of the problem. FAIL is set at .TRUE., if the smallest
635 C eigenvalue of the matrix A is close to 0.
637 implicit double precision (a-h,o-z)
638 double precision a(lda,n),x(n),eig(100),vec(3,3),work(100)
639 double precision eig_limit /1.0D-8/
642 c print '(a)','Enter MULT_NORM.'
644 C Find the smallest eigenvalue of the matrix A.
647 c print '(8f10.5)',(a(i,j),j=1,n)
650 call f02faf('V','U',2,a,lda,eig,work,100,ifail)
652 call djacob(2,lda,10000,1.0d-10,a,vec,eig)
654 c print '(8f10.5)',(eig(i),i=1,n)
657 c print '(8f10.5)',(a(i,j),j=1,n)
659 if (eig(1).lt.eig_limit) then
660 print *,'From Mult_Norm: Eigenvalues of A are too small.'
665 C Generate points following the normal distributions along the principal
666 C axes of the moment matrix. Store in WORK.
669 sigma=1.0D0/dsqrt(eig(i))
671 work(i)=anorm_distr(0.0D0,sigma,-alim,alim)
674 C Transform the vector of normal variables back to the original basis.
685 c------------------------------------------------------------------------
686 subroutine mult_norm1(lda,n,a,z,box,x,fail)
688 C Generate the vector X whose elements obey the multi-gaussian multi-dimensional
689 C distribution from sum_{i=1}^m W(i)exp[-0.5*X'(i)A(i)X(i)]. LDA is the
690 C leading dimension of the moment matrix A, n is the dimension of the
691 C distribution, nlob is the number of lobes. FAIL is set at .TRUE., if the
692 C smallest eigenvalue of the matrix A is close to 0.
694 implicit real*8 (a-h,o-z)
699 double precision a(lda,n),z(n),x(n),box(n,n)
700 double precision etmp
701 include 'COMMON.IOUNITS'
703 include 'COMMON.SETUP'
707 C Generate points following the normal distributions along the principal
708 C axes of the moment matrix. Store in WORK.
710 cd print *,'CG Processor',me,' entered MultNorm1.'
711 cd print '(2(1pe14.4),3x,1pe14.4)',((a(i,j),j=1,2),z(i),i=1,2)
713 cd print *,i,box(1,i),box(2,i)
717 if (istep.gt.10000) then
718 c write (iout,'(a,i4,2a)') 'CG Processor: ',me,': too many steps',
720 c write (*,'(a,i4,2a)') 'CG Processor: ',me,': too many steps',
722 c write (iout,*) 'box',box
723 c write (iout,*) 'a',a
724 c write (iout,*) 'z',z
729 x(i)=ran_number(box(1,i),box(2,i))
734 ww=ww+0.5D0*a(i,i)*xi*xi
736 ww=ww+a(i,j)*xi*pinorm(x(j)-z(j))
739 dec=ran_number(0.0D0,1.0D0)
740 c print *,(x(i),i=1,n),ww,dexp(-ww),dec
741 crc if (dec.gt.dexp(-ww)) goto 10
747 if (dec.gt.etmp) goto 10
748 cd print *,'CG Processor',me,' exitting MultNorm1.'
752 crc--------------------------------------
753 subroutine overlap_sc(scfail)
754 c Internal and cartesian coordinates must be consistent as input,
755 c and will be up-to-date on return.
756 c At the end of this procedure, scfail is true if there are
757 c overlapping residues left, or false otherwise (success)
758 implicit real*8 (a-h,o-z)
760 include 'COMMON.CHAIN'
761 include 'COMMON.INTERACT'
762 include 'COMMON.FFIELD'
764 include 'COMMON.SBRIDGE'
765 include 'COMMON.IOUNITS'
766 logical had_overlaps,fail,scfail
767 integer ioverlap(maxres),ioverlap_last
770 call overlap_sc_list(ioverlap,ioverlap_last)
771 if (ioverlap_last.gt.0) then
772 write (iout,*) '#OVERLAPing residues ',ioverlap_last
773 write (iout,'(20i4)') (ioverlap(k),k=1,ioverlap_last)
779 if (ioverlap_last.eq.0) exit
781 do ires=1,ioverlap_last
787 do while (fail.and.nsi.le.maxsi)
788 call gen_side(iti,theta(i+1),alph(i),omeg(i),fail)
796 call overlap_sc_list(ioverlap,ioverlap_last)
797 c write (iout,*) 'Overlaping residues ',ioverlap_last,
798 c & (ioverlap(j),j=1,ioverlap_last)
801 if (k.le.1000.and.ioverlap_last.eq.0) then
803 if (had_overlaps) then
804 write (iout,*) '#OVERLAPing all corrected after ',k,
805 & ' random generation'
809 write (iout,*) '#OVERLAPing NOT all corrected ',ioverlap_last
810 write (iout,'(20i4)') (ioverlap(j),j=1,ioverlap_last)
816 write (iout,'(a30,i5,a12,i4)')
817 & '#OVERLAP FAIL in gen_side after',maxsi,
823 subroutine overlap_sc_list(ioverlap,ioverlap_last)
824 implicit real*8 (a-h,o-z)
827 include 'COMMON.LOCAL'
828 include 'COMMON.IOUNITS'
829 include 'COMMON.CHAIN'
830 include 'COMMON.INTERACT'
831 include 'COMMON.FFIELD'
833 include 'COMMON.CALC'
835 integer ioverlap(maxres),ioverlap_last
839 C Check for SC-SC overlaps and mark residues
840 c print *,'>>overlap_sc nnt=',nnt,' nct=',nct
844 itypi1=abs(itype(i+1))
848 dxi=dc_norm(1,nres+i)
849 dyi=dc_norm(2,nres+i)
850 dzi=dc_norm(3,nres+i)
851 dsci_inv=dsc_inv(itypi)
854 do j=istart(i,iint),iend(i,iint)
857 dscj_inv=dsc_inv(itypj)
858 sig0ij=sigma(itypi,itypj)
859 chi1=chi(itypi,itypj)
860 chi2=chi(itypj,itypi)
867 alf12=0.5D0*(alf1+alf2)
869 rcomp=sigmaii(itypi,itypj)
871 rcomp=sigma(itypi,itypj)
873 c print '(2(a3,2i3),a3,2f10.5)',
874 c & ' i=',i,iti,' j=',j,itj,' d=',dist(nres+i,nres+j)
879 dxj=dc_norm(1,nres+j)
880 dyj=dc_norm(2,nres+j)
881 dzj=dc_norm(3,nres+j)
882 rrij=1.0D0/(xj*xj+yj*yj+zj*zj)
886 sig=sig0ij*dsqrt(sigsq)
887 rij_shift=1.0D0/rij-sig+sig0ij
889 ct if ( 1.0/rij .lt. redfac*rcomp .or.
890 ct & rij_shift.le.0.0D0 ) then
891 if ( rij_shift.le.0.0D0 ) then
892 cd write (iout,'(a,i3,a,i3,a,f10.5,a,3f10.5)')
893 cd & 'overlap SC-SC: i=',i,' j=',j,
894 cd & ' dist=',dist(nres+i,nres+j),' rcomp=',
895 cd & rcomp,1.0/rij,rij_shift
896 ioverlap_last=ioverlap_last+1
897 ioverlap(ioverlap_last)=i
898 do k=1,ioverlap_last-1
899 if (ioverlap(k).eq.i) ioverlap_last=ioverlap_last-1
901 ioverlap_last=ioverlap_last+1
902 ioverlap(ioverlap_last)=j
903 do k=1,ioverlap_last-1
904 if (ioverlap(k).eq.j) ioverlap_last=ioverlap_last-1