2 c------------------------------------------------
3 c The driver for molecular dynamics subroutines
4 c------------------------------------------------
5 implicit real*8 (a-h,o-z)
11 include 'COMMON.SETUP'
12 include 'COMMON.CONTROL'
16 include 'COMMON.LANGEVIN'
18 include 'COMMON.LANGEVIN.lang0'
20 include 'COMMON.CHAIN'
21 include 'COMMON.DERIV'
23 include 'COMMON.LOCAL'
24 include 'COMMON.INTERACT'
25 include 'COMMON.IOUNITS'
26 include 'COMMON.NAMES'
27 include 'COMMON.TIME1'
28 double precision cm(3),L(3),vcm(3)
30 double precision v_work(maxres6),v_transf(maxres6)
39 if (ilen(tmpdir).gt.0)
40 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_"
41 & //liczba(:ilen(liczba))//'.rst')
43 if (ilen(tmpdir).gt.0)
44 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_"//'.rst')
51 write (iout,'(20(1h=),a20,20(1h=))') "MD calculation started"
57 c Determine the inverse of the inertia matrix.
58 call setup_MD_matrices
62 t_MDsetup = MPI_Wtime()-tt0
64 t_MDsetup = tcpu()-tt0
67 c Entering the MD loop
73 if (lang.eq.2 .or. lang.eq.3) then
77 call sd_verlet_p_setup
79 call sd_verlet_ciccotti_setup
83 pfric0_mat(i,j,0)=pfric_mat(i,j)
84 afric0_mat(i,j,0)=afric_mat(i,j)
85 vfric0_mat(i,j,0)=vfric_mat(i,j)
86 prand0_mat(i,j,0)=prand_mat(i,j)
87 vrand0_mat1(i,j,0)=vrand_mat1(i,j)
88 vrand0_mat2(i,j,0)=vrand_mat2(i,j)
97 & "LANG=2 or 3 NOT SUPPORTED. Recompile without -DLANG0"
99 call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE)
103 else if (lang.eq.1 .or. lang.eq.4) then
107 t_langsetup=MPI_Wtime()-tt0
110 t_langsetup=tcpu()-tt0
113 do itime=1,n_timestep
115 if (lang.gt.0 .and. surfarea .and.
116 & mod(itime,reset_fricmat).eq.0) then
117 if (lang.eq.2 .or. lang.eq.3) then
121 call sd_verlet_p_setup
123 call sd_verlet_ciccotti_setup
127 pfric0_mat(i,j,0)=pfric_mat(i,j)
128 afric0_mat(i,j,0)=afric_mat(i,j)
129 vfric0_mat(i,j,0)=vfric_mat(i,j)
130 prand0_mat(i,j,0)=prand_mat(i,j)
131 vrand0_mat1(i,j,0)=vrand_mat1(i,j)
132 vrand0_mat2(i,j,0)=vrand_mat2(i,j)
137 flag_stoch(i)=.false.
140 else if (lang.eq.1 .or. lang.eq.4) then
143 write (iout,'(a,i10)')
144 & "Friction matrix reset based on surface area, itime",itime
146 if (reset_vel .and. tbf .and. lang.eq.0
147 & .and. mod(itime,count_reset_vel).eq.0) then
149 write(iout,'(a,f20.2)')
150 & "Velocities reset to random values, time",totT
153 d_t_old(j,i)=d_t(j,i)
157 if (reset_moment .and. mod(itime,count_reset_moment).eq.0) then
161 d_t(j,0)=d_t(j,0)-vcm(j)
164 kinetic_T=2.0d0/(dimen3*Rb)*EK
165 scalfac=dsqrt(T_bath/kinetic_T)
166 write(iout,'(a,f20.2)') "Momenta zeroed out, time",totT
169 d_t_old(j,i)=scalfac*d_t(j,i)
175 c Time-reversible RESPA algorithm
176 c (Tuckerman et al., J. Chem. Phys., 97, 1990, 1992)
177 call RESPA_step(itime)
179 c Variable time step algorithm.
180 call velverlet_step(itime)
184 call brown_step(itime)
186 print *,"Brown dynamics not here!"
188 call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE)
194 if (mod(itime,ntwe).eq.0) call statout(itime)
207 if (itype(i).ne.10) then
210 v_work(ind)=d_t(j,i+nres)
215 write (66,'(80f10.5)')
216 & ((d_t(j,i),j=1,3),i=0,nres-1),((d_t(j,i+nres),j=1,3),i=1,nres)
220 v_transf(i)=v_transf(i)+gvec(j,i)*v_work(j)
222 v_transf(i)= v_transf(i)*dsqrt(geigen(i))
224 write (67,'(80f10.5)') (v_transf(i),i=1,ind)
227 if (mod(itime,ntwx).eq.0) then
228 write (tytul,'("time",f8.2)') totT
230 call pdbout(potE,tytul,ipdb)
235 if (rstcount.eq.1000.or.itime.eq.n_timestep) then
236 open(irest2,file=rest2name,status='unknown')
237 write(irest2,*) totT,EK,potE,totE,t_bath
239 write (irest2,'(3e15.5)') (d_t(j,i),j=1,3)
242 write (irest2,'(3e15.5)') (dc(j,i),j=1,3)
253 write (iout,'(//35(1h=),a10,35(1h=)/10(/a40,1pe15.5))')
255 & 'MD calculations setup:',t_MDsetup,
256 & 'Energy & gradient evaluation:',t_enegrad,
257 & 'Stochastic MD setup:',t_langsetup,
258 & 'Stochastic MD step setup:',t_sdsetup,
260 write (iout,'(/28(1h=),a25,27(1h=))')
261 & ' End of MD calculation '
263 write (iout,*) "time for etotal",t_etotal," elong",t_elong,
265 write (iout,*) "time_fric",time_fric," time_stoch",time_stoch,
266 & " time_fricmatmult",time_fricmatmult," time_fsample ",
271 c-------------------------------------------------------------------------------
272 subroutine velverlet_step(itime)
273 c-------------------------------------------------------------------------------
274 c Perform a single velocity Verlet step; the time step can be rescaled if
275 c increments in accelerations exceed the threshold
276 c-------------------------------------------------------------------------------
277 implicit real*8 (a-h,o-z)
281 integer ierror,ierrcode
283 include 'COMMON.SETUP'
284 include 'COMMON.CONTROL'
288 include 'COMMON.LANGEVIN'
290 include 'COMMON.LANGEVIN.lang0'
292 include 'COMMON.CHAIN'
293 include 'COMMON.DERIV'
295 include 'COMMON.LOCAL'
296 include 'COMMON.INTERACT'
297 include 'COMMON.IOUNITS'
298 include 'COMMON.NAMES'
299 include 'COMMON.TIME1'
300 include 'COMMON.MUCA'
301 double precision vcm(3),incr(3)
302 double precision cm(3),L(3)
303 integer ilen,count,rstcount
306 integer maxcount_scale /20/
308 double precision stochforcvec(MAXRES6)
309 common /stochcalc/ stochforcvec
312 double precision HNose1,HNose,HNose_nh,H,vtnp(maxres6)
313 double precision vtnp_(maxres6),vtnp_a(maxres6)
319 else if (lang.eq.2 .or. lang.eq.3) then
321 call stochastic_force(stochforcvec)
324 & "LANG=2 or 3 NOT SUPPORTED. Recompile without -DLANG0"
326 call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE)
333 icount_scale=icount_scale+1
334 if (icount_scale.gt.maxcount_scale) then
336 & "ERROR: too many attempts at scaling down the time step. ",
337 & "amax=",amax,"epdrift=",epdrift,
338 & "damax=",damax,"edriftmax=",edriftmax,
342 call MPI_Abort(MPI_COMM_WORLD,IERROR,IERRCODE)
346 c First step of the velocity Verlet algorithm
351 else if (lang.eq.3) then
353 call sd_verlet1_ciccotti
355 else if (lang.eq.1) then
363 call nhcint(EK,scale_nh,wdti,wdti2,wdti4,wdti8)
366 d_t_old(j,i)=d_t_old(j,i)*scale_nh
372 c Build the chain from the newly calculated coordinates
374 if (rattle) call rattle1
376 if (large.and. mod(itime,ntwe).eq.0) then
377 write (iout,*) "Cartesian and internal coordinates: step 1"
382 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(dc(j,i),j=1,3),
383 & (dc(j,i+nres),j=1,3)
385 write (iout,*) "Accelerations"
387 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
388 & (d_a(j,i+nres),j=1,3)
390 write (iout,*) "Velocities, step 1"
392 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
393 & (d_t(j,i+nres),j=1,3)
402 c Calculate energy and forces
404 call etotal(potEcomp)
407 t_etotal=t_etotal+MPI_Wtime()-tt0
409 t_etotal=t_etotal+tcpu()-tt0
413 potE=potEcomp(0)-potEcomp(20)
415 c Get the new accelerations
418 t_enegrad=t_enegrad+MPI_Wtime()-tt0
420 t_enegrad=t_enegrad+tcpu()-tt0
422 c Determine maximum acceleration and scale down the timestep if needed
424 amax=amax/(itime_scal+1)**2
425 call predict_edrift(epdrift)
426 if (amax/(itime_scal+1).gt.damax .or. epdrift.gt.edriftmax) then
427 c Maximum acceleration or maximum predicted energy drift exceeded, rescale the time step
429 ifac_time=dmax1(dlog(amax/damax),dlog(epdrift/edriftmax))
431 itime_scal=itime_scal+ifac_time
432 c fac_time=dmin1(damax/amax,0.5d0)
433 fac_time=0.5d0**ifac_time
434 d_time=d_time*fac_time
435 if (lang.eq.2 .or. lang.eq.3) then
437 c write (iout,*) "Calling sd_verlet_setup: 1"
438 c Rescale the stochastic forces and recalculate or restore
439 c the matrices of tinker integrator
440 if (itime_scal.gt.maxflag_stoch) then
441 if (large) write (iout,'(a,i5,a)')
442 & "Calculate matrices for stochastic step;",
443 & " itime_scal ",itime_scal
445 call sd_verlet_p_setup
447 call sd_verlet_ciccotti_setup
449 write (iout,'(2a,i3,a,i3,1h.)')
450 & "Warning: cannot store matrices for stochastic",
451 & " integration because the index",itime_scal,
452 & " is greater than",maxflag_stoch
453 write (iout,'(2a)')"Increase MAXFLAG_STOCH or use direct",
454 & " integration Langevin algorithm for better efficiency."
455 else if (flag_stoch(itime_scal)) then
456 if (large) write (iout,'(a,i5,a,l1)')
457 & "Restore matrices for stochastic step; itime_scal ",
458 & itime_scal," flag ",flag_stoch(itime_scal)
461 pfric_mat(i,j)=pfric0_mat(i,j,itime_scal)
462 afric_mat(i,j)=afric0_mat(i,j,itime_scal)
463 vfric_mat(i,j)=vfric0_mat(i,j,itime_scal)
464 prand_mat(i,j)=prand0_mat(i,j,itime_scal)
465 vrand_mat1(i,j)=vrand0_mat1(i,j,itime_scal)
466 vrand_mat2(i,j)=vrand0_mat2(i,j,itime_scal)
470 if (large) write (iout,'(2a,i5,a,l1)')
471 & "Calculate & store matrices for stochastic step;",
472 & " itime_scal ",itime_scal," flag ",flag_stoch(itime_scal)
474 call sd_verlet_p_setup
476 call sd_verlet_ciccotti_setup
478 flag_stoch(ifac_time)=.true.
481 pfric0_mat(i,j,itime_scal)=pfric_mat(i,j)
482 afric0_mat(i,j,itime_scal)=afric_mat(i,j)
483 vfric0_mat(i,j,itime_scal)=vfric_mat(i,j)
484 prand0_mat(i,j,itime_scal)=prand_mat(i,j)
485 vrand0_mat1(i,j,itime_scal)=vrand_mat1(i,j)
486 vrand0_mat2(i,j,itime_scal)=vrand_mat2(i,j)
490 fac_time=1.0d0/dsqrt(fac_time)
492 stochforcvec(i)=fac_time*stochforcvec(i)
495 else if (lang.eq.1) then
496 c Rescale the accelerations due to stochastic forces
497 fac_time=1.0d0/dsqrt(fac_time)
499 d_as_work(i)=d_as_work(i)*fac_time
502 if (large) write (iout,'(a,i10,a,f8.6,a,i3,a,i3)')
503 & "itime",itime," Timestep scaled down to ",
504 & d_time," ifac_time",ifac_time," itime_scal",itime_scal
506 c Second step of the velocity Verlet algorithm
511 else if (lang.eq.3) then
513 call sd_verlet2_ciccotti
515 else if (lang.eq.1) then
525 call nhcint(EK,scale_nh,wdti,wdti2,wdti4,wdti8)
528 d_t(j,i)=d_t(j,i)*scale_nh
533 if (rattle) call rattle2
535 if (d_time.ne.d_time0) then
538 if (lang.eq.2 .or. lang.eq.3) then
539 if (large) write (iout,'(a)')
540 & "Restore original matrices for stochastic step"
541 c write (iout,*) "Calling sd_verlet_setup: 2"
542 c Restore the matrices of tinker integrator if the time step has been restored
545 pfric_mat(i,j)=pfric0_mat(i,j,0)
546 afric_mat(i,j)=afric0_mat(i,j,0)
547 vfric_mat(i,j)=vfric0_mat(i,j,0)
548 prand_mat(i,j)=prand0_mat(i,j,0)
549 vrand_mat1(i,j)=vrand0_mat1(i,j,0)
550 vrand_mat2(i,j)=vrand0_mat2(i,j,0)
559 c Calculate the kinetic and the total energy and the kinetic temperature
560 if (tnp .or. tnp1) then
563 d_t_old(j,i)=d_t(j,i)
564 d_t(j,i)=d_t(j,i)/s_np
572 c write (iout,*) "step",itime," EK",EK," EK1",EK1
574 c Couple the system to Berendsen bath if needed
575 if (tbf .and. lang.eq.0) then
578 kinetic_T=2.0d0/(dimen3*Rb)*EK
579 c Backup the coordinates, velocities, and accelerations
583 if(.not.(tnp .or. tnp1)) d_t_old(j,i)=d_t(j,i)
584 d_a_old(j,i)=d_a(j,i)
588 if (mod(itime,ntwe).eq.0) then
590 if(tnp .or. tnp1) then
591 HNose1=Hnose(EK,s_np,potE,pi_np,Q_np,t_bath,dimen3)
593 cd write (iout,'(a,10f)') "hhh",EK,s_np,potE,pi_np,H0
594 cd & ,EK+potE+pi_np**2/(2*Q_np)+dimen3*0.001986d0*t_bath*log(s_np)
595 cd write (iout,*) "HHH H=",H,abs(HNose1-H0)/H0
600 HNose1=Hnose_nh(EK,potE)
603 cd write (iout,*) "HHH H=",H,abs(HNose1-H0)/H0
619 if (itype(i).ne.10) then
623 vtnp(itnp)=d_t(j,inres)
628 c Transform velocities from UNRES coordinate space to cartesian and Gvec
635 vtnp_(i)=vtnp_(i)+Gvec(j,i)*vtnp(j)
636 vtnp_a(i)=vtnp_a(i)+A(i,j)*vtnp(j)
638 vtnp_(i)=vtnp_(i)*dsqrt(geigen(i))
642 write (iout,'("WWW",i3,3f10.5)') i,vtnp(i),vtnp_(i),vtnp_a(i)
645 write (iout,*) "Velocities, step 2"
647 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
648 & (d_t(j,i+nres),j=1,3)
655 c-------------------------------------------------------------------------------
656 subroutine RESPA_step(itime)
657 c-------------------------------------------------------------------------------
658 c Perform a single RESPA step.
659 c-------------------------------------------------------------------------------
660 implicit real*8 (a-h,o-z)
664 integer IERROR,ERRCODE
666 include 'COMMON.SETUP'
667 include 'COMMON.CONTROL'
671 include 'COMMON.LANGEVIN'
673 include 'COMMON.LANGEVIN.lang0'
675 include 'COMMON.CHAIN'
676 include 'COMMON.DERIV'
678 include 'COMMON.LOCAL'
679 include 'COMMON.INTERACT'
680 include 'COMMON.IOUNITS'
681 include 'COMMON.NAMES'
682 include 'COMMON.TIME1'
683 double precision energia_short(0:n_ene),
684 & energia_long(0:n_ene)
685 double precision cm(3),L(3),vcm(3),incr(3)
686 double precision dc_old0(3,0:maxres2),d_t_old0(3,0:maxres2),
687 & d_a_old0(3,0:maxres2)
688 integer ilen,count,rstcount
691 integer maxcount_scale /10/
692 common /gucio/ cm,energia_short
693 double precision stochforcvec(MAXRES6)
694 common /stochcalc/ stochforcvec
697 double precision vtnp(maxres6), vtnp_(maxres6), vtnp_a(maxres6)
698 common /cipiszcze/ itt
701 if (large.and. mod(itime,ntwe).eq.0) then
702 write (iout,*) "***************** RESPA itime",itime
703 write (iout,*) "Cartesian and internal coordinates: step 0"
705 call pdbout(0.0d0,"cipiszcze",iout)
707 write (iout,*) "Accelerations from long-range forces"
709 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
710 & (d_a(j,i+nres),j=1,3)
712 write (iout,*) "Velocities, step 0"
714 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
715 & (d_t(j,i+nres),j=1,3)
720 c Perform the initial RESPA step (increment velocities)
721 c write (iout,*) "*********************** RESPA ini"
727 if (tnh.and..not.xiresp) then
728 call nhcint(EK,scale_nh,wdti,wdti2,wdti4,wdti8)
731 d_t(j,i)=d_t(j,i)*scale_nh
738 cd if(tnp .or. tnp1) then
739 cd write (iout,'(a,3f)') "EE1 NP S, pi",totT, s_np, pi_np
743 if (mod(itime,ntwe).eq.0 .and. large) then
744 write (iout,*) "Velocities, end"
746 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
747 & (d_t(j,i+nres),j=1,3)
751 c Compute the short-range forces
757 C 7/2/2009 commented out
759 c call etotal_short(energia_short)
760 if (tnp.or.tnp1) potE=energia_short(0)
763 C 7/2/2009 Copy accelerations due to short-lange forces from previous MD step
766 d_a(j,i)=d_a_short(j,i)
770 if (large.and. mod(itime,ntwe).eq.0) then
771 write (iout,*) "energia_short",energia_short(0)
772 write (iout,*) "Accelerations from short-range forces"
774 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
775 & (d_a(j,i+nres),j=1,3)
780 t_enegrad=t_enegrad+MPI_Wtime()-tt0
782 t_enegrad=t_enegrad+tcpu()-tt0
787 if(.not.(tnp .or. tnp1)) d_t_old(j,i)=d_t(j,i)
788 d_a_old(j,i)=d_a(j,i)
791 c 6/30/08 A-MTS: attempt at increasing the split number
794 dc_old0(j,i)=dc_old(j,i)
795 d_t_old0(j,i)=d_t_old(j,i)
796 d_a_old0(j,i)=d_a_old(j,i)
799 if (ntime_split.gt.ntime_split0) ntime_split=ntime_split/2
800 if (ntime_split.lt.ntime_split0) ntime_split=ntime_split0
807 c write (iout,*) "itime",itime," ntime_split",ntime_split
808 c Split the time step
809 d_time=d_time0/ntime_split
810 c Perform the short-range RESPA steps (velocity Verlet increments of
811 c positions and velocities using short-range forces)
812 c write (iout,*) "*********************** RESPA split"
813 do itsplit=1,ntime_split
816 else if (lang.eq.2 .or. lang.eq.3) then
818 call stochastic_force(stochforcvec)
821 & "LANG=2 or 3 NOT SUPPORTED. Recompile without -DLANG0"
823 call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE)
828 c First step of the velocity Verlet algorithm
833 else if (lang.eq.3) then
835 call sd_verlet1_ciccotti
837 else if (lang.eq.1) then
840 call tnp1_respa_i_step1
842 call tnp_respa_i_step1
844 if (tnh.and.xiresp) then
846 call nhcint(EK,scale_nh,wdtii,wdtii2,wdtii4,wdtii8)
849 d_t_old(j,i)=d_t_old(j,i)*scale_nh
852 cd write(iout,*) "SSS1",itsplit,EK,scale_nh
856 c Build the chain from the newly calculated coordinates
858 if (rattle) call rattle1
860 if (large.and. mod(itime,ntwe).eq.0) then
861 write (iout,*) "***** ITSPLIT",itsplit
862 write (iout,*) "Cartesian and internal coordinates: step 1"
863 call pdbout(0.0d0,"cipiszcze",iout)
866 write (iout,*) "Velocities, step 1"
868 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
869 & (d_t(j,i+nres),j=1,3)
878 c Calculate energy and forces
880 call etotal_short(energia_short)
882 potE=energia_short(0)
885 t_eshort=t_eshort+MPI_Wtime()-tt0
887 t_eshort=t_eshort+tcpu()-tt0
891 c Get the new accelerations
893 C 7/2/2009 Copy accelerations due to short-lange forces to an auxiliary array
896 d_a_short(j,i)=d_a(j,i)
900 if (large.and. mod(itime,ntwe).eq.0) then
901 write (iout,*)"energia_short",energia_short(0)
902 write (iout,*) "Accelerations from short-range forces"
904 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
905 & (d_a(j,i+nres),j=1,3)
910 c Determine maximum acceleration and scale down the timestep if needed
912 amax=amax/ntime_split**2
913 call predict_edrift(epdrift)
914 if (ntwe.gt.0 .and. large .and. mod(itime,ntwe).eq.0)
915 & write (iout,*) "amax",amax," damax",damax,
916 & " epdrift",epdrift," epdriftmax",epdriftmax
917 c Exit loop and try with increased split number if the change of
918 c acceleration is too big
919 if (amax.gt.damax .or. epdrift.gt.edriftmax) then
920 if (ntime_split.lt.maxtime_split) then
922 ntime_split=ntime_split*2
925 dc_old(j,i)=dc_old0(j,i)
926 d_t_old(j,i)=d_t_old0(j,i)
927 d_a_old(j,i)=d_a_old0(j,i)
930 write (iout,*) "acceleration/energy drift too large",amax,
931 & epdrift," split increased to ",ntime_split," itime",itime,
936 & "Uh-hu. Bumpy landscape. Maximum splitting number",
938 & " already reached!!! Trying to carry on!"
942 t_enegrad=t_enegrad+MPI_Wtime()-tt0
944 t_enegrad=t_enegrad+tcpu()-tt0
946 c Second step of the velocity Verlet algorithm
951 else if (lang.eq.3) then
953 call sd_verlet2_ciccotti
955 else if (lang.eq.1) then
958 call tnp1_respa_i_step2
960 call tnp_respa_i_step2
963 if (tnh.and.xiresp) then
965 call nhcint(EK,scale_nh,wdtii,wdtii2,wdtii4,wdtii8)
968 d_t(j,i)=d_t(j,i)*scale_nh
971 cd write(iout,*) "SSS2",itsplit,EK,scale_nh
974 if (rattle) call rattle2
975 c Backup the coordinates, velocities, and accelerations
976 if (tnp .or. tnp1) then
979 d_t_old(j,i)=d_t(j,i)
980 if (tnp) d_t(j,i)=d_t(j,i)/s_np
981 if (tnp1) d_t(j,i)=d_t(j,i)/s_np
989 if(.not.(tnp .or. tnp1)) d_t_old(j,i)=d_t(j,i)
990 d_a_old(j,i)=d_a(j,i)
997 c Restore the time step
999 c Compute long-range forces
1006 call etotal_long(energia_long)
1007 E_long=energia_long(0)
1008 potE=energia_short(0)+energia_long(0)
1011 t_elong=t_elong+MPI_Wtime()-tt0
1013 t_elong=t_elong+tcpu()-tt0
1019 t_enegrad=t_enegrad+MPI_Wtime()-tt0
1021 t_enegrad=t_enegrad+tcpu()-tt0
1023 c Compute accelerations from long-range forces
1025 if (large.and. mod(itime,ntwe).eq.0) then
1026 write (iout,*) "energia_long",energia_long(0)
1027 write (iout,*) "Cartesian and internal coordinates: step 2"
1029 call pdbout(0.0d0,"cipiszcze",iout)
1031 write (iout,*) "Accelerations from long-range forces"
1033 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
1034 & (d_a(j,i+nres),j=1,3)
1036 write (iout,*) "Velocities, step 2"
1038 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
1039 & (d_t(j,i+nres),j=1,3)
1043 c Compute the final RESPA step (increment velocities)
1044 c write (iout,*) "*********************** RESPA fin"
1046 call tnp_respa_step2
1048 call tnp_respa_step2
1051 if (tnh.and..not.xiresp) then
1053 call nhcint(EK,scale_nh,wdti,wdti2,wdti4,wdti8)
1056 d_t(j,i)=d_t(j,i)*scale_nh
1062 if (tnp .or. tnp1) then
1065 d_t(j,i)=d_t_old(j,i)/s_np
1070 c Compute the complete potential energy
1072 potEcomp(i)=energia_short(i)+energia_long(i)
1074 potE=potEcomp(0)-potEcomp(20)
1075 c potE=energia_short(0)+energia_long(0)
1077 c Calculate the kinetic and the total energy and the kinetic temperature
1080 c Couple the system to Berendsen bath if needed
1081 if (tbf .and. lang.eq.0) then
1084 kinetic_T=2.0d0/(dimen3*Rb)*EK
1085 c Backup the coordinates, velocities, and accelerations
1087 if (mod(itime,ntwe).eq.0 .and. large) then
1088 write (iout,*) "Velocities, end"
1090 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
1091 & (d_t(j,i+nres),j=1,3)
1095 if (mod(itime,ntwe).eq.0) then
1097 if(tnp .or. tnp1) then
1099 write (iout,'(a3,7f)') "TTT",EK,s_np,potE,pi_np,Csplit,
1100 & E_long,energia_short(0)
1102 write (iout,'(a3,7f20.10)') "TTT",EK,s_np,potE,pi_np,Csplit,
1103 & E_long,energia_short(0)
1105 HNose1=Hnose(EK,s_np,potE,pi_np,Q_np,t_bath,dimen3)
1107 cd write (iout,'(a,10f)') "hhh",EK,s_np,potE,pi_np,H0
1108 cd & ,EK+potE+pi_np**2/(2*Q_np)+dimen3*0.001986d0*t_bath*log(s_np)
1109 cd write (iout,*) "HHH H=",H,abs(HNose1-H0)/H0
1111 cd write (iout,'(a,3f)') "EE2 NP S, pi",totT, s_np, pi_np
1115 HNose1=Hnose_nh(EK,potE)
1117 cd write (iout,*) "HHH H=",H,abs(HNose1-H0)/H0
1135 if (itype(i).ne.10) then
1139 vtnp(itnp)=d_t(j,inres)
1144 c Transform velocities from UNRES coordinate space to cartesian and Gvec
1151 vtnp_(i)=vtnp_(i)+Gvec(j,i)*vtnp(j)
1152 vtnp_a(i)=vtnp_a(i)+A(i,j)*vtnp(j)
1154 vtnp_(i)=vtnp_(i)*dsqrt(geigen(i))
1158 write (iout,'("WWW",i3,3f10.5)') i,vtnp(i),vtnp_(i),vtnp_a(i)
1168 c---------------------------------------------------------------------
1169 subroutine RESPA_vel
1170 c First and last RESPA step (incrementing velocities using long-range
1172 implicit real*8 (a-h,o-z)
1173 include 'DIMENSIONS'
1174 include 'COMMON.CONTROL'
1175 include 'COMMON.VAR'
1177 include 'COMMON.CHAIN'
1178 include 'COMMON.DERIV'
1179 include 'COMMON.GEO'
1180 include 'COMMON.LOCAL'
1181 include 'COMMON.INTERACT'
1182 include 'COMMON.IOUNITS'
1183 include 'COMMON.NAMES'
1185 d_t(j,0)=d_t(j,0)+0.5d0*d_a(j,0)*d_time
1189 d_t(j,i)=d_t(j,i)+0.5d0*d_a(j,i)*d_time
1193 if (itype(i).ne.10) then
1196 d_t(j,inres)=d_t(j,inres)+0.5d0*d_a(j,inres)*d_time
1202 c-----------------------------------------------------------------
1204 c Applying velocity Verlet algorithm - step 1 to coordinates
1205 implicit real*8 (a-h,o-z)
1206 include 'DIMENSIONS'
1207 include 'COMMON.CONTROL'
1208 include 'COMMON.VAR'
1210 include 'COMMON.CHAIN'
1211 include 'COMMON.DERIV'
1212 include 'COMMON.GEO'
1213 include 'COMMON.LOCAL'
1214 include 'COMMON.INTERACT'
1215 include 'COMMON.IOUNITS'
1216 include 'COMMON.NAMES'
1217 double precision adt,adt2
1220 write (iout,*) "VELVERLET1 START: DC"
1222 write (iout,'(i3,3f10.5,5x,3f10.5)') i,(dc(j,i),j=1,3),
1223 & (dc(j,i+nres),j=1,3)
1227 adt=d_a_old(j,0)*d_time
1229 dc(j,0)=dc_old(j,0)+(d_t_old(j,0)+adt2)*d_time
1230 d_t_new(j,0)=d_t_old(j,0)+adt2
1231 d_t(j,0)=d_t_old(j,0)+adt
1235 adt=d_a_old(j,i)*d_time
1237 dc(j,i)=dc_old(j,i)+(d_t_old(j,i)+adt2)*d_time
1238 d_t_new(j,i)=d_t_old(j,i)+adt2
1239 d_t(j,i)=d_t_old(j,i)+adt
1243 if (itype(i).ne.10) then
1246 adt=d_a_old(j,inres)*d_time
1248 dc(j,inres)=dc_old(j,inres)+(d_t_old(j,inres)+adt2)*d_time
1249 d_t_new(j,inres)=d_t_old(j,inres)+adt2
1250 d_t(j,inres)=d_t_old(j,inres)+adt
1255 write (iout,*) "VELVERLET1 END: DC"
1257 write (iout,'(i3,3f10.5,5x,3f10.5)') i,(dc(j,i),j=1,3),
1258 & (dc(j,i+nres),j=1,3)
1263 c---------------------------------------------------------------------
1265 c Step 2 of the velocity Verlet algorithm: update velocities
1266 implicit real*8 (a-h,o-z)
1267 include 'DIMENSIONS'
1268 include 'COMMON.CONTROL'
1269 include 'COMMON.VAR'
1271 include 'COMMON.CHAIN'
1272 include 'COMMON.DERIV'
1273 include 'COMMON.GEO'
1274 include 'COMMON.LOCAL'
1275 include 'COMMON.INTERACT'
1276 include 'COMMON.IOUNITS'
1277 include 'COMMON.NAMES'
1279 d_t(j,0)=d_t_new(j,0)+0.5d0*d_a(j,0)*d_time
1283 d_t(j,i)=d_t_new(j,i)+0.5d0*d_a(j,i)*d_time
1287 if (itype(i).ne.10) then
1290 d_t(j,inres)=d_t_new(j,inres)+0.5d0*d_a(j,inres)*d_time
1296 c-----------------------------------------------------------------
1297 subroutine sddir_precalc
1298 c Applying velocity Verlet algorithm - step 1 to coordinates
1299 implicit real*8 (a-h,o-z)
1300 include 'DIMENSIONS'
1304 include 'COMMON.CONTROL'
1305 include 'COMMON.VAR'
1308 include 'COMMON.LANGEVIN'
1310 include 'COMMON.LANGEVIN.lang0'
1312 include 'COMMON.CHAIN'
1313 include 'COMMON.DERIV'
1314 include 'COMMON.GEO'
1315 include 'COMMON.LOCAL'
1316 include 'COMMON.INTERACT'
1317 include 'COMMON.IOUNITS'
1318 include 'COMMON.NAMES'
1319 include 'COMMON.TIME1'
1320 double precision stochforcvec(MAXRES6)
1321 common /stochcalc/ stochforcvec
1323 c Compute friction and stochastic forces
1332 time_fric=time_fric+MPI_Wtime()-time00
1335 time_fric=time_fric+tcpu()-time00
1338 call stochastic_force(stochforcvec)
1340 time_stoch=time_stoch+MPI_Wtime()-time00
1342 time_stoch=time_stoch+tcpu()-time00
1345 c Compute the acceleration due to friction forces (d_af_work) and stochastic
1346 c forces (d_as_work)
1348 call ginv_mult(fric_work, d_af_work)
1349 call ginv_mult(stochforcvec, d_as_work)
1352 c---------------------------------------------------------------------
1353 subroutine sddir_verlet1
1354 c Applying velocity Verlet algorithm - step 1 to velocities
1355 implicit real*8 (a-h,o-z)
1356 include 'DIMENSIONS'
1357 include 'COMMON.CONTROL'
1358 include 'COMMON.VAR'
1361 include 'COMMON.LANGEVIN'
1363 include 'COMMON.LANGEVIN.lang0'
1365 include 'COMMON.CHAIN'
1366 include 'COMMON.DERIV'
1367 include 'COMMON.GEO'
1368 include 'COMMON.LOCAL'
1369 include 'COMMON.INTERACT'
1370 include 'COMMON.IOUNITS'
1371 include 'COMMON.NAMES'
1372 c Revised 3/31/05 AL: correlation between random contributions to
1373 c position and velocity increments included.
1374 double precision sqrt13 /0.57735026918962576451d0/ ! 1/sqrt(3)
1375 double precision adt,adt2
1377 c Add the contribution from BOTH friction and stochastic force to the
1378 c coordinates, but ONLY the contribution from the friction forces to velocities
1381 adt=(d_a_old(j,0)+d_af_work(j))*d_time
1382 adt2=0.5d0*adt+sqrt13*d_as_work(j)*d_time
1383 dc(j,0)=dc_old(j,0)+(d_t_old(j,0)+adt2)*d_time
1384 d_t_new(j,0)=d_t_old(j,0)+0.5d0*adt
1385 d_t(j,0)=d_t_old(j,0)+adt
1390 adt=(d_a_old(j,i)+d_af_work(ind+j))*d_time
1391 adt2=0.5d0*adt+sqrt13*d_as_work(ind+j)*d_time
1392 dc(j,i)=dc_old(j,i)+(d_t_old(j,i)+adt2)*d_time
1393 d_t_new(j,i)=d_t_old(j,i)+0.5d0*adt
1394 d_t(j,i)=d_t_old(j,i)+adt
1399 if (itype(i).ne.10) then
1402 adt=(d_a_old(j,inres)+d_af_work(ind+j))*d_time
1403 adt2=0.5d0*adt+sqrt13*d_as_work(ind+j)*d_time
1404 dc(j,inres)=dc_old(j,inres)+(d_t_old(j,inres)+adt2)*d_time
1405 d_t_new(j,inres)=d_t_old(j,inres)+0.5d0*adt
1406 d_t(j,inres)=d_t_old(j,inres)+adt
1413 c---------------------------------------------------------------------
1414 subroutine sddir_verlet2
1415 c Calculating the adjusted velocities for accelerations
1416 implicit real*8 (a-h,o-z)
1417 include 'DIMENSIONS'
1418 include 'COMMON.CONTROL'
1419 include 'COMMON.VAR'
1422 include 'COMMON.LANGEVIN'
1424 include 'COMMON.LANGEVIN.lang0'
1426 include 'COMMON.CHAIN'
1427 include 'COMMON.DERIV'
1428 include 'COMMON.GEO'
1429 include 'COMMON.LOCAL'
1430 include 'COMMON.INTERACT'
1431 include 'COMMON.IOUNITS'
1432 include 'COMMON.NAMES'
1433 double precision stochforcvec(MAXRES6),d_as_work1(MAXRES6)
1434 double precision cos60 /0.5d0/, sin60 /0.86602540378443864676d0/
1435 c Revised 3/31/05 AL: correlation between random contributions to
1436 c position and velocity increments included.
1437 c The correlation coefficients are calculated at low-friction limit.
1438 c Also, friction forces are now not calculated with new velocities.
1440 c call friction_force
1441 call stochastic_force(stochforcvec)
1443 c Compute the acceleration due to friction forces (d_af_work) and stochastic
1444 c forces (d_as_work)
1446 call ginv_mult(stochforcvec, d_as_work1)
1452 d_t(j,0)=d_t_new(j,0)+(0.5d0*(d_a(j,0)+d_af_work(j))
1453 & +sin60*d_as_work(j)+cos60*d_as_work1(j))*d_time
1458 d_t(j,i)=d_t_new(j,i)+(0.5d0*(d_a(j,i)+d_af_work(ind+j))
1459 & +sin60*d_as_work(ind+j)+cos60*d_as_work1(ind+j))*d_time
1464 if (itype(i).ne.10) then
1467 d_t(j,inres)=d_t_new(j,inres)+(0.5d0*(d_a(j,inres)
1468 & +d_af_work(ind+j))+sin60*d_as_work(ind+j)
1469 & +cos60*d_as_work1(ind+j))*d_time
1476 c---------------------------------------------------------------------
1477 subroutine max_accel
1479 c Find the maximum difference in the accelerations of the the sites
1480 c at the beginning and the end of the time step.
1482 implicit real*8 (a-h,o-z)
1483 include 'DIMENSIONS'
1484 include 'COMMON.CONTROL'
1485 include 'COMMON.VAR'
1487 include 'COMMON.CHAIN'
1488 include 'COMMON.DERIV'
1489 include 'COMMON.GEO'
1490 include 'COMMON.LOCAL'
1491 include 'COMMON.INTERACT'
1492 include 'COMMON.IOUNITS'
1493 double precision aux(3),accel(3),accel_old(3),dacc
1495 c aux(j)=d_a(j,0)-d_a_old(j,0)
1496 accel_old(j)=d_a_old(j,0)
1503 c 7/3/08 changed to asymmetric difference
1505 c accel(j)=aux(j)+0.5d0*(d_a(j,i)-d_a_old(j,i))
1506 accel_old(j)=accel_old(j)+0.5d0*d_a_old(j,i)
1507 accel(j)=accel(j)+0.5d0*d_a(j,i)
1508 c if (dabs(accel(j)).gt.amax) amax=dabs(accel(j))
1509 if (dabs(accel(j)).gt.dabs(accel_old(j))) then
1510 dacc=dabs(accel(j)-accel_old(j))
1511 if (dacc.gt.amax) amax=dacc
1519 accel_old(j)=d_a_old(j,0)
1524 accel_old(j)=accel_old(j)+d_a_old(j,1)
1525 accel(j)=accel(j)+d_a(j,1)
1529 if (itype(i).ne.10) then
1531 c accel(j)=accel(j)+d_a(j,i+nres)-d_a_old(j,i+nres)
1532 accel_old(j)=accel_old(j)+d_a_old(j,i+nres)
1533 accel(j)=accel(j)+d_a(j,i+nres)
1537 c if (dabs(accel(j)).gt.amax) amax=dabs(accel(j))
1538 if (dabs(accel(j)).gt.dabs(accel_old(j))) then
1539 dacc=dabs(accel(j)-accel_old(j))
1540 if (dacc.gt.amax) amax=dacc
1544 accel_old(j)=accel_old(j)+d_a_old(j,i)
1545 accel(j)=accel(j)+d_a(j,i)
1546 c aux(j)=aux(j)+d_a(j,i)-d_a_old(j,i)
1551 c---------------------------------------------------------------------
1552 subroutine predict_edrift(epdrift)
1554 c Predict the drift of the potential energy
1556 implicit real*8 (a-h,o-z)
1557 include 'DIMENSIONS'
1558 include 'COMMON.CONTROL'
1559 include 'COMMON.VAR'
1561 include 'COMMON.CHAIN'
1562 include 'COMMON.DERIV'
1563 include 'COMMON.GEO'
1564 include 'COMMON.LOCAL'
1565 include 'COMMON.INTERACT'
1566 include 'COMMON.IOUNITS'
1567 include 'COMMON.MUCA'
1568 double precision epdrift,epdriftij
1569 c Drift of the potential energy
1575 epdriftij=dabs((d_a(j,i)-d_a_old(j,i))*gcart(j,i))
1576 if (lmuca) epdriftij=epdriftij*factor
1577 c write (iout,*) "back",i,j,epdriftij
1578 if (epdriftij.gt.epdrift) epdrift=epdriftij
1582 if (itype(i).ne.10) then
1585 & dabs((d_a(j,i+nres)-d_a_old(j,i+nres))*gxcart(j,i))
1586 if (lmuca) epdriftij=epdriftij*factor
1587 c write (iout,*) "side",i,j,epdriftij
1588 if (epdriftij.gt.epdrift) epdrift=epdriftij
1592 epdrift=0.5d0*epdrift*d_time*d_time
1593 c write (iout,*) "epdrift",epdrift
1596 c-----------------------------------------------------------------------
1597 subroutine verlet_bath
1599 c Coupling to the thermostat by using the Berendsen algorithm
1601 implicit real*8 (a-h,o-z)
1602 include 'DIMENSIONS'
1603 include 'COMMON.CONTROL'
1604 include 'COMMON.VAR'
1606 include 'COMMON.CHAIN'
1607 include 'COMMON.DERIV'
1608 include 'COMMON.GEO'
1609 include 'COMMON.LOCAL'
1610 include 'COMMON.INTERACT'
1611 include 'COMMON.IOUNITS'
1612 include 'COMMON.NAMES'
1613 double precision T_half,fact
1615 T_half=2.0d0/(dimen3*Rb)*EK
1616 fact=dsqrt(1.0d0+(d_time/tau_bath)*(t_bath/T_half-1.0d0))
1617 c write(iout,*) "T_half", T_half
1618 c write(iout,*) "EK", EK
1619 c write(iout,*) "fact", fact
1621 d_t(j,0)=fact*d_t(j,0)
1625 d_t(j,i)=fact*d_t(j,i)
1629 if (itype(i).ne.10) then
1632 d_t(j,inres)=fact*d_t(j,inres)
1638 c---------------------------------------------------------
1640 c Set up the initial conditions of a MD simulation
1641 implicit real*8 (a-h,o-z)
1642 include 'DIMENSIONS'
1646 integer IERROR,ERRCODE
1648 include 'COMMON.SETUP'
1649 include 'COMMON.CONTROL'
1650 include 'COMMON.VAR'
1653 include 'COMMON.LANGEVIN'
1655 include 'COMMON.LANGEVIN.lang0'
1657 include 'COMMON.CHAIN'
1658 include 'COMMON.DERIV'
1659 include 'COMMON.GEO'
1660 include 'COMMON.LOCAL'
1661 include 'COMMON.INTERACT'
1662 include 'COMMON.IOUNITS'
1663 include 'COMMON.NAMES'
1664 include 'COMMON.REMD'
1665 real*8 energia_long(0:n_ene),
1666 & energia_short(0:n_ene),vcm(3),incr(3),E_short
1667 double precision cm(3),L(3),xv,sigv,lowb,highb
1668 double precision varia(maxvar)
1676 c write(iout,*) "d_time", d_time
1677 c Compute the standard deviations of stochastic forces for Langevin dynamics
1678 c if the friction coefficients do not depend on surface area
1679 if (lang.gt.0 .and. .not.surfarea) then
1681 stdforcp(i)=stdfp*dsqrt(gamp)
1684 stdforcsc(i)=stdfsc(itype(i))*dsqrt(gamsc(itype(i)))
1687 c Open the pdb file for snapshotshots
1690 if (ilen(tmpdir).gt.0)
1691 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD"//
1692 & liczba(:ilen(liczba))//".pdb")
1694 & file=prefix(:ilen(prefix))//"_MD"//liczba(:ilen(liczba))
1698 if (ilen(tmpdir).gt.0 .and. (me.eq.king .or. .not.traj1file))
1699 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD"//
1700 & liczba(:ilen(liczba))//".x")
1701 cartname=prefix(:ilen(prefix))//"_MD"//liczba(:ilen(liczba))
1704 if (ilen(tmpdir).gt.0 .and. (me.eq.king .or. .not.traj1file))
1705 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD"//
1706 & liczba(:ilen(liczba))//".cx")
1707 cartname=prefix(:ilen(prefix))//"_MD"//liczba(:ilen(liczba))
1713 if (ilen(tmpdir).gt.0)
1714 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD.pdb")
1715 open(ipdb,file=prefix(:ilen(prefix))//"_MD.pdb")
1717 if (ilen(tmpdir).gt.0)
1718 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD.cx")
1719 cartname=prefix(:ilen(prefix))//"_MD.cx"
1723 write (qstr,'(256(1h ))')
1726 iq = qinfrag(i,iset)*10
1727 iw = wfrag(i,iset)/100
1729 if(me.eq.king.or..not.out1file)
1730 & write (iout,*) "Frag",qinfrag(i,iset),wfrag(i,iset),iq,iw
1731 write (qstr(ipos:ipos+6),'(2h_f,i1,1h_,i1,1h_,i1)') i,iq,iw
1736 iq = qinpair(i,iset)*10
1737 iw = wpair(i,iset)/100
1739 if(me.eq.king.or..not.out1file)
1740 & write (iout,*) "Pair",i,qinpair(i,iset),wpair(i,iset),iq,iw
1741 write (qstr(ipos:ipos+6),'(2h_p,i1,1h_,i1,1h_,i1)') i,iq,iw
1745 c pdbname=pdbname(:ilen(pdbname)-4)//qstr(:ipos-1)//'.pdb'
1747 c cartname=cartname(:ilen(cartname)-2)//qstr(:ipos-1)//'.x'
1749 c cartname=cartname(:ilen(cartname)-3)//qstr(:ipos-1)//'.cx'
1751 c statname=statname(:ilen(statname)-5)//qstr(:ipos-1)//'.stat'
1755 if (restart1file) then
1757 & inquire(file=mremd_rst_name,exist=file_exist)
1758 write (*,*) me," Before broadcast: file_exist",file_exist
1760 call MPI_Bcast(file_exist,1,MPI_LOGICAL,king,CG_COMM,
1762 write (*,*) me," After broadcast: file_exist",file_exist
1764 c inquire(file=mremd_rst_name,exist=file_exist)
1765 if(me.eq.king.or..not.out1file)
1766 & write(iout,*) "Initial state read by master and distributed"
1768 if (ilen(tmpdir).gt.0)
1769 & call copy_to_tmp(pref_orig(:ilen(pref_orig))//'_'
1770 & //liczba(:ilen(liczba))//'.rst')
1771 inquire(file=rest2name,exist=file_exist)
1774 if(.not.restart1file) then
1775 if(me.eq.king.or..not.out1file)
1776 & write(iout,*) "Initial state will be read from file ",
1777 & rest2name(:ilen(rest2name))
1780 call rescale_weights(t_bath)
1782 if(me.eq.king.or..not.out1file)then
1783 if (restart1file) then
1784 write(iout,*) "File ",mremd_rst_name(:ilen(mremd_rst_name)),
1787 write(iout,*) "File ",rest2name(:ilen(rest2name)),
1790 write(iout,*) "Initial velocities randomly generated"
1796 c Generate initial velocities
1797 if(me.eq.king.or..not.out1file)
1798 & write(iout,*) "Initial velocities randomly generated"
1802 c rest2name = prefix(:ilen(prefix))//'.rst'
1803 if(me.eq.king.or..not.out1file)then
1804 write (iout,*) "Initial velocities"
1806 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
1807 & (d_t(j,i+nres),j=1,3)
1810 c Zeroing the total angular momentum of the system
1811 write(iout,*) "Calling the zero-angular
1812 & momentum subroutine"
1815 c Getting the potential energy and forces and velocities and accelerations
1817 c write (iout,*) "velocity of the center of the mass:"
1818 c write (iout,*) (vcm(j),j=1,3)
1820 d_t(j,0)=d_t(j,0)-vcm(j)
1822 c Removing the velocity of the center of mass
1824 if(me.eq.king.or..not.out1file)then
1825 write (iout,*) "vcm right after adjustment:"
1826 write (iout,*) (vcm(j),j=1,3)
1831 if(iranconf.ne.0) then
1833 print *, 'Calling OVERLAP_SC'
1834 call overlap_sc(fail)
1838 call sc_move(2,nres-1,10,1d10,nft_sc,etot)
1839 print *,'SC_move',nft_sc,etot
1840 if(me.eq.king.or..not.out1file)
1841 & write(iout,*) 'SC_move',nft_sc,etot
1845 print *, 'Calling MINIM_DC'
1846 call minim_dc(etot,iretcode,nfun)
1848 call geom_to_var(nvar,varia)
1849 print *,'Calling MINIMIZE.'
1850 call minimize(etot,varia,iretcode,nfun)
1851 call var_to_geom(nvar,varia)
1853 if(me.eq.king.or..not.out1file)
1854 & write(iout,*) 'SUMSL return code is',iretcode,' eval ',nfun
1857 call chainbuild_cart
1862 kinetic_T=2.0d0/(dimen3*Rb)*EK
1863 if(me.eq.king.or..not.out1file)then
1873 call etotal(potEcomp)
1876 t_etotal=t_etotal+MPI_Wtime()-tt0
1878 t_etotal=t_etotal+tcpu()-tt0
1883 if(tnp .or. tnp1) then
1886 HNose1=Hnose(EK,s_np,potE,pi_np,Q_np,t_bath,dimen3)
1888 write(iout,*) 'H0= ',H0
1892 HNose1=Hnose_nh(EK,potE)
1894 write (iout,*) 'H0= ',H0
1900 if(me.eq.king.or..not.out1file)
1901 & write(iout,*) 'HMC',hmc_etot,potE,EK
1912 if (amax*d_time .gt. dvmax) then
1913 d_time=d_time*dvmax/amax
1914 if(me.eq.king.or..not.out1file) write (iout,*)
1915 & "Time step reduced to",d_time,
1916 & " because of too large initial acceleration."
1918 if(me.eq.king.or..not.out1file)then
1919 write(iout,*) "Potential energy and its components"
1920 call enerprint(potEcomp)
1921 c write(iout,*) (potEcomp(i),i=0,n_ene)
1923 potE=potEcomp(0)-potEcomp(20)
1926 if (ntwe.ne.0) call statout(itime)
1927 if(me.eq.king.or..not.out1file)
1928 & write (iout,'(/a/3(a25,1pe14.5/))') "Initial:",
1929 & " Kinetic energy",EK," potential energy",potE,
1930 & " total energy",totE," maximum acceleration ",
1933 write (iout,*) "Initial coordinates"
1935 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(c(j,i),j=1,3),
1936 & (c(j,i+nres),j=1,3)
1938 write (iout,*) "Initial dC"
1940 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(dc(j,i),j=1,3),
1941 & (dc(j,i+nres),j=1,3)
1943 write (iout,*) "Initial velocities"
1944 write (iout,"(13x,' backbone ',23x,' side chain')")
1946 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3),
1947 & (d_t(j,i+nres),j=1,3)
1949 write (iout,*) "Initial accelerations"
1951 c write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
1952 write (iout,'(i3,3f15.10,3x,3f15.10)') i,(d_a(j,i),j=1,3),
1953 & (d_a(j,i+nres),j=1,3)
1959 d_t_old(j,i)=d_t(j,i)
1960 d_a_old(j,i)=d_a(j,i)
1962 c write (iout,*) "dc_old",i,(dc_old(j,i),j=1,3)
1971 call etotal_short(energia_short)
1974 t_eshort=t_eshort+MPI_Wtime()-tt0
1976 t_eshort=t_eshort+tcpu()-tt0
1980 if(tnp .or. tnp1) then
1981 E_short=energia_short(0)
1982 HNose1=Hnose(EK,s_np,E_short,pi_np,Q_np,t_bath,dimen3)
1985 c_new_var_csplit Csplit=H0-E_long
1986 c Csplit = H0-energia_short(0)
1987 write(iout,*) 'Csplit= ',Csplit
1993 if(.not.out1file .and. large) then
1994 write (iout,*) "energia_long",energia_long(0),
1995 & " energia_short",energia_short(0),
1996 & " total",energia_long(0)+energia_short(0)
1997 write (iout,*) "Initial fast-force accelerations"
1999 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
2000 & (d_a(j,i+nres),j=1,3)
2003 C 7/2/2009 Copy accelerations due to short-lange forces to an auxiliary array
2006 d_a_short(j,i)=d_a(j,i)
2015 call etotal_long(energia_long)
2018 t_elong=t_elong+MPI_Wtime()-tt0
2020 t_elong=t_elong+tcpu()-tt0
2025 if(.not.out1file .and. large) then
2026 write (iout,*) "energia_long",energia_long(0)
2027 write (iout,*) "Initial slow-force accelerations"
2029 write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3),
2030 & (d_a(j,i+nres),j=1,3)
2034 t_enegrad=t_enegrad+MPI_Wtime()-tt0
2036 t_enegrad=t_enegrad+tcpu()-tt0
2044 c-----------------------------------------------------------
2045 subroutine random_vel
2046 implicit real*8 (a-h,o-z)
2047 include 'DIMENSIONS'
2048 include 'COMMON.CONTROL'
2049 include 'COMMON.VAR'
2052 include 'COMMON.LANGEVIN'
2054 include 'COMMON.LANGEVIN.lang0'
2056 include 'COMMON.CHAIN'
2057 include 'COMMON.DERIV'
2058 include 'COMMON.GEO'
2059 include 'COMMON.LOCAL'
2060 include 'COMMON.INTERACT'
2061 include 'COMMON.IOUNITS'
2062 include 'COMMON.NAMES'
2063 include 'COMMON.TIME1'
2064 double precision xv,sigv,lowb,highb
2065 c Generate random velocities from Gaussian distribution of mean 0 and std of KT/m
2066 c First generate velocities in the eigenspace of the G matrix
2067 c write (iout,*) "Calling random_vel dimen dimen3",dimen,dimen3
2069 c write (iout,*) "RANDOM_VEL dimen",dimen
2075 sigv=dsqrt((Rb*t_bath)/geigen(i))
2078 d_t_work_new(ii)=anorm_distr(xv,sigv,lowb,highb)
2079 c write (iout,*) "i",i," ii",ii," geigen",geigen(i),
2080 c & " d_t_work_new",d_t_work_new(ii)
2090 c Ek1=Ek1+0.5d0*geigen(i)*d_t_work_new(ii)**2
2093 c write (iout,*) "Ek from eigenvectors",Ek1
2095 c Transform velocities to UNRES coordinate space
2101 d_t_work(ind)=d_t_work(ind)
2102 & +Gvec(i,j)*d_t_work_new((j-1)*3+k+1)
2104 c write (iout,*) "i",i," ind",ind," d_t_work",d_t_work(ind)
2108 c Transfer to the d_t vector
2110 d_t(j,0)=d_t_work(j)
2116 d_t(j,i)=d_t_work(ind)
2120 c write (iout,*) "d_t",i,(d_t(j,i),j=1,3)
2124 if (itype(i).ne.10) then
2127 d_t(j,i+nres)=d_t_work(ind)
2132 c write (iout,*) "Kinetic energy",Ek,EK1," kinetic temperature",
2133 c & 2.0d0/(dimen3*Rb)*EK,2.0d0/(dimen3*Rb)*EK1
2138 c-----------------------------------------------------------
2139 subroutine sd_verlet_p_setup
2140 c Sets up the parameters of stochastic Verlet algorithm
2141 implicit real*8 (a-h,o-z)
2142 include 'DIMENSIONS'
2146 include 'COMMON.CONTROL'
2147 include 'COMMON.VAR'
2150 include 'COMMON.LANGEVIN'
2152 include 'COMMON.LANGEVIN.lang0'
2154 include 'COMMON.CHAIN'
2155 include 'COMMON.DERIV'
2156 include 'COMMON.GEO'
2157 include 'COMMON.LOCAL'
2158 include 'COMMON.INTERACT'
2159 include 'COMMON.IOUNITS'
2160 include 'COMMON.NAMES'
2161 include 'COMMON.TIME1'
2162 double precision emgdt(MAXRES6),
2163 & pterm,vterm,rho,rhoc,vsig,
2164 & pfric_vec(MAXRES6),vfric_vec(MAXRES6),
2165 & afric_vec(MAXRES6),prand_vec(MAXRES6),
2166 & vrand_vec1(MAXRES6),vrand_vec2(MAXRES6)
2167 logical lprn /.false./
2168 double precision zero /1.0d-8/, gdt_radius /0.05d0/
2169 double precision ktm
2176 c AL 8/17/04 Code adapted from tinker
2178 c Get the frictional and random terms for stochastic dynamics in the
2179 c eigenspace of mass-scaled UNRES friction matrix
2182 gdt = fricgam(i) * d_time
2184 c Stochastic dynamics reduces to simple MD for zero friction
2186 if (gdt .le. zero) then
2187 pfric_vec(i) = 1.0d0
2188 vfric_vec(i) = d_time
2189 afric_vec(i) = 0.5d0 * d_time * d_time
2190 prand_vec(i) = 0.0d0
2191 vrand_vec1(i) = 0.0d0
2192 vrand_vec2(i) = 0.0d0
2194 c Analytical expressions when friction coefficient is large
2197 if (gdt .ge. gdt_radius) then
2200 vfric_vec(i) = (1.0d0-egdt) / fricgam(i)
2201 afric_vec(i) = (d_time-vfric_vec(i)) / fricgam(i)
2202 pterm = 2.0d0*gdt - 3.0d0 + (4.0d0-egdt)*egdt
2203 vterm = 1.0d0 - egdt**2
2204 rho = (1.0d0-egdt)**2 / sqrt(pterm*vterm)
2206 c Use series expansions when friction coefficient is small
2217 afric_vec(i) = (gdt2/2.0d0 - gdt3/6.0d0 + gdt4/24.0d0
2218 & - gdt5/120.0d0 + gdt6/720.0d0
2219 & - gdt7/5040.0d0 + gdt8/40320.0d0
2220 & - gdt9/362880.0d0) / fricgam(i)**2
2221 vfric_vec(i) = d_time - fricgam(i)*afric_vec(i)
2222 pfric_vec(i) = 1.0d0 - fricgam(i)*vfric_vec(i)
2223 pterm = 2.0d0*gdt3/3.0d0 - gdt4/2.0d0
2224 & + 7.0d0*gdt5/30.0d0 - gdt6/12.0d0
2225 & + 31.0d0*gdt7/1260.0d0 - gdt8/160.0d0
2226 & + 127.0d0*gdt9/90720.0d0
2227 vterm = 2.0d0*gdt - 2.0d0*gdt2 + 4.0d0*gdt3/3.0d0
2228 & - 2.0d0*gdt4/3.0d0 + 4.0d0*gdt5/15.0d0
2229 & - 4.0d0*gdt6/45.0d0 + 8.0d0*gdt7/315.0d0
2230 & - 2.0d0*gdt8/315.0d0 + 4.0d0*gdt9/2835.0d0
2231 rho = sqrt(3.0d0) * (0.5d0 - 3.0d0*gdt/16.0d0
2232 & - 17.0d0*gdt2/1280.0d0
2233 & + 17.0d0*gdt3/6144.0d0
2234 & + 40967.0d0*gdt4/34406400.0d0
2235 & - 57203.0d0*gdt5/275251200.0d0
2236 & - 1429487.0d0*gdt6/13212057600.0d0)
2239 c Compute the scaling factors of random terms for the nonzero friction case
2241 ktm = 0.5d0*d_time/fricgam(i)
2242 psig = dsqrt(ktm*pterm) / fricgam(i)
2243 vsig = dsqrt(ktm*vterm)
2244 rhoc = dsqrt(1.0d0 - rho*rho)
2246 vrand_vec1(i) = vsig * rho
2247 vrand_vec2(i) = vsig * rhoc
2252 & "pfric_vec, vfric_vec, afric_vec, prand_vec, vrand_vec1,",
2255 write (iout,'(i5,6e15.5)') i,pfric_vec(i),vfric_vec(i),
2256 & afric_vec(i),prand_vec(i),vrand_vec1(i),vrand_vec2(i)
2260 c Transform from the eigenspace of mass-scaled friction matrix to UNRES variables
2263 call eigtransf(dimen,maxres2,mt3,mt2,pfric_vec,pfric_mat)
2264 call eigtransf(dimen,maxres2,mt3,mt2,vfric_vec,vfric_mat)
2265 call eigtransf(dimen,maxres2,mt3,mt2,afric_vec,afric_mat)
2266 call eigtransf(dimen,maxres2,mt3,mt1,prand_vec,prand_mat)
2267 call eigtransf(dimen,maxres2,mt3,mt1,vrand_vec1,vrand_mat1)
2268 call eigtransf(dimen,maxres2,mt3,mt1,vrand_vec2,vrand_mat2)
2271 t_sdsetup=t_sdsetup+MPI_Wtime()
2273 t_sdsetup=t_sdsetup+tcpu()-tt0
2277 c-------------------------------------------------------------
2278 subroutine eigtransf1(n,ndim,ab,d,c)
2281 double precision ab(ndim,ndim,n),c(ndim,n),d(ndim)
2287 c(i,j)=c(i,j)+ab(k,j,i)*d(k)
2293 c-------------------------------------------------------------
2294 subroutine eigtransf(n,ndim,a,b,d,c)
2297 double precision a(ndim,n),b(ndim,n),c(ndim,n),d(ndim)
2303 c(i,j)=c(i,j)+a(i,k)*b(k,j)*d(k)
2309 c-------------------------------------------------------------
2310 subroutine sd_verlet1
2311 c Applying stochastic velocity Verlet algorithm - step 1 to velocities
2312 implicit real*8 (a-h,o-z)
2313 include 'DIMENSIONS'
2314 include 'COMMON.CONTROL'
2315 include 'COMMON.VAR'
2318 include 'COMMON.LANGEVIN'
2320 include 'COMMON.LANGEVIN.lang0'
2322 include 'COMMON.CHAIN'
2323 include 'COMMON.DERIV'
2324 include 'COMMON.GEO'
2325 include 'COMMON.LOCAL'
2326 include 'COMMON.INTERACT'
2327 include 'COMMON.IOUNITS'
2328 include 'COMMON.NAMES'
2329 double precision stochforcvec(MAXRES6)
2330 common /stochcalc/ stochforcvec
2331 logical lprn /.false./
2333 c write (iout,*) "dc_old"
2335 c write (iout,'(i5,3f10.5,5x,3f10.5)')
2336 c & i,(dc_old(j,i),j=1,3),(dc_old(j,i+nres),j=1,3)
2339 dc_work(j)=dc_old(j,0)
2340 d_t_work(j)=d_t_old(j,0)
2341 d_a_work(j)=d_a_old(j,0)
2346 dc_work(ind+j)=dc_old(j,i)
2347 d_t_work(ind+j)=d_t_old(j,i)
2348 d_a_work(ind+j)=d_a_old(j,i)
2353 if (itype(i).ne.10) then
2355 dc_work(ind+j)=dc_old(j,i+nres)
2356 d_t_work(ind+j)=d_t_old(j,i+nres)
2357 d_a_work(ind+j)=d_a_old(j,i+nres)
2365 & "pfric_mat, vfric_mat, afric_mat, prand_mat, vrand_mat1,",
2369 write (iout,'(2i5,6e15.5)') i,j,pfric_mat(i,j),
2370 & vfric_mat(i,j),afric_mat(i,j),
2371 & prand_mat(i,j),vrand_mat1(i,j),vrand_mat2(i,j)
2379 dc_work(i)=dc_work(i)+vfric_mat(i,j)*d_t_work(j)
2380 & +afric_mat(i,j)*d_a_work(j)+prand_mat(i,j)*stochforcvec(j)
2381 ddt1=ddt1+pfric_mat(i,j)*d_t_work(j)
2382 ddt2=ddt2+vfric_mat(i,j)*d_a_work(j)
2384 d_t_work_new(i)=ddt1+0.5d0*ddt2
2385 d_t_work(i)=ddt1+ddt2
2390 d_t(j,0)=d_t_work(j)
2395 dc(j,i)=dc_work(ind+j)
2396 d_t(j,i)=d_t_work(ind+j)
2401 if (itype(i).ne.10) then
2404 dc(j,inres)=dc_work(ind+j)
2405 d_t(j,inres)=d_t_work(ind+j)
2412 c--------------------------------------------------------------------------
2413 subroutine sd_verlet2
2414 c Calculating the adjusted velocities for accelerations
2415 implicit real*8 (a-h,o-z)
2416 include 'DIMENSIONS'
2417 include 'COMMON.CONTROL'
2418 include 'COMMON.VAR'
2421 include 'COMMON.LANGEVIN'
2423 include 'COMMON.LANGEVIN.lang0'
2425 include 'COMMON.CHAIN'
2426 include 'COMMON.DERIV'
2427 include 'COMMON.GEO'
2428 include 'COMMON.LOCAL'
2429 include 'COMMON.INTERACT'
2430 include 'COMMON.IOUNITS'
2431 include 'COMMON.NAMES'
2432 double precision stochforcvec(MAXRES6),stochforcvecV(MAXRES6)
2433 common /stochcalc/ stochforcvec
2435 c Compute the stochastic forces which contribute to velocity change
2437 call stochastic_force(stochforcvecV)
2444 ddt1=ddt1+vfric_mat(i,j)*d_a_work(j)
2445 ddt2=ddt2+vrand_mat1(i,j)*stochforcvec(j)+
2446 & vrand_mat2(i,j)*stochforcvecV(j)
2448 d_t_work(i)=d_t_work_new(i)+0.5d0*ddt1+ddt2
2452 d_t(j,0)=d_t_work(j)
2457 d_t(j,i)=d_t_work(ind+j)
2462 if (itype(i).ne.10) then
2465 d_t(j,inres)=d_t_work(ind+j)
2472 c-----------------------------------------------------------
2473 subroutine sd_verlet_ciccotti_setup
2474 c Sets up the parameters of stochastic velocity Verlet algorithmi; Ciccotti's
2476 implicit real*8 (a-h,o-z)
2477 include 'DIMENSIONS'
2481 include 'COMMON.CONTROL'
2482 include 'COMMON.VAR'
2485 include 'COMMON.LANGEVIN'
2487 include 'COMMON.LANGEVIN.lang0'
2489 include 'COMMON.CHAIN'
2490 include 'COMMON.DERIV'
2491 include 'COMMON.GEO'
2492 include 'COMMON.LOCAL'
2493 include 'COMMON.INTERACT'
2494 include 'COMMON.IOUNITS'
2495 include 'COMMON.NAMES'
2496 include 'COMMON.TIME1'
2497 double precision emgdt(MAXRES6),
2498 & pterm,vterm,rho,rhoc,vsig,
2499 & pfric_vec(MAXRES6),vfric_vec(MAXRES6),
2500 & afric_vec(MAXRES6),prand_vec(MAXRES6),
2501 & vrand_vec1(MAXRES6),vrand_vec2(MAXRES6)
2502 logical lprn /.false./
2503 double precision zero /1.0d-8/, gdt_radius /0.05d0/
2504 double precision ktm
2511 c AL 8/17/04 Code adapted from tinker
2513 c Get the frictional and random terms for stochastic dynamics in the
2514 c eigenspace of mass-scaled UNRES friction matrix
2517 write (iout,*) "i",i," fricgam",fricgam(i)
2518 gdt = fricgam(i) * d_time
2520 c Stochastic dynamics reduces to simple MD for zero friction
2522 if (gdt .le. zero) then
2523 pfric_vec(i) = 1.0d0
2524 vfric_vec(i) = d_time
2525 afric_vec(i) = 0.5d0*d_time*d_time
2526 prand_vec(i) = afric_vec(i)
2527 vrand_vec2(i) = vfric_vec(i)
2529 c Analytical expressions when friction coefficient is large
2534 vfric_vec(i) = dexp(-0.5d0*gdt)*d_time
2535 afric_vec(i) = 0.5d0*dexp(-0.25d0*gdt)*d_time*d_time
2536 prand_vec(i) = afric_vec(i)
2537 vrand_vec2(i) = vfric_vec(i)
2539 c Compute the scaling factors of random terms for the nonzero friction case
2541 c ktm = 0.5d0*d_time/fricgam(i)
2542 c psig = dsqrt(ktm*pterm) / fricgam(i)
2543 c vsig = dsqrt(ktm*vterm)
2544 c prand_vec(i) = psig*afric_vec(i)
2545 c vrand_vec2(i) = vsig*vfric_vec(i)
2550 & "pfric_vec, vfric_vec, afric_vec, prand_vec, vrand_vec1,",
2553 write (iout,'(i5,6e15.5)') i,pfric_vec(i),vfric_vec(i),
2554 & afric_vec(i),prand_vec(i),vrand_vec1(i),vrand_vec2(i)
2558 c Transform from the eigenspace of mass-scaled friction matrix to UNRES variables
2560 call eigtransf(dimen,maxres2,mt3,mt2,pfric_vec,pfric_mat)
2561 call eigtransf(dimen,maxres2,mt3,mt2,vfric_vec,vfric_mat)
2562 call eigtransf(dimen,maxres2,mt3,mt2,afric_vec,afric_mat)
2563 call eigtransf(dimen,maxres2,mt3,mt1,prand_vec,prand_mat)
2564 call eigtransf(dimen,maxres2,mt3,mt1,vrand_vec2,vrand_mat2)
2566 t_sdsetup=t_sdsetup+MPI_Wtime()
2568 t_sdsetup=t_sdsetup+tcpu()-tt0
2572 c-------------------------------------------------------------
2573 subroutine sd_verlet1_ciccotti
2574 c Applying stochastic velocity Verlet algorithm - step 1 to velocities
2575 implicit real*8 (a-h,o-z)
2576 include 'DIMENSIONS'
2580 include 'COMMON.CONTROL'
2581 include 'COMMON.VAR'
2584 include 'COMMON.LANGEVIN'
2586 include 'COMMON.LANGEVIN.lang0'
2588 include 'COMMON.CHAIN'
2589 include 'COMMON.DERIV'
2590 include 'COMMON.GEO'
2591 include 'COMMON.LOCAL'
2592 include 'COMMON.INTERACT'
2593 include 'COMMON.IOUNITS'
2594 include 'COMMON.NAMES'
2595 double precision stochforcvec(MAXRES6)
2596 common /stochcalc/ stochforcvec
2597 logical lprn /.false./
2599 c write (iout,*) "dc_old"
2601 c write (iout,'(i5,3f10.5,5x,3f10.5)')
2602 c & i,(dc_old(j,i),j=1,3),(dc_old(j,i+nres),j=1,3)
2605 dc_work(j)=dc_old(j,0)
2606 d_t_work(j)=d_t_old(j,0)
2607 d_a_work(j)=d_a_old(j,0)
2612 dc_work(ind+j)=dc_old(j,i)
2613 d_t_work(ind+j)=d_t_old(j,i)
2614 d_a_work(ind+j)=d_a_old(j,i)
2619 if (itype(i).ne.10) then
2621 dc_work(ind+j)=dc_old(j,i+nres)
2622 d_t_work(ind+j)=d_t_old(j,i+nres)
2623 d_a_work(ind+j)=d_a_old(j,i+nres)
2632 & "pfric_mat, vfric_mat, afric_mat, prand_mat, vrand_mat1,",
2636 write (iout,'(2i5,6e15.5)') i,j,pfric_mat(i,j),
2637 & vfric_mat(i,j),afric_mat(i,j),
2638 & prand_mat(i,j),vrand_mat1(i,j),vrand_mat2(i,j)
2646 dc_work(i)=dc_work(i)+vfric_mat(i,j)*d_t_work(j)
2647 & +afric_mat(i,j)*d_a_work(j)+prand_mat(i,j)*stochforcvec(j)
2648 ddt1=ddt1+pfric_mat(i,j)*d_t_work(j)
2649 ddt2=ddt2+vfric_mat(i,j)*d_a_work(j)
2651 d_t_work_new(i)=ddt1+0.5d0*ddt2
2652 d_t_work(i)=ddt1+ddt2
2657 d_t(j,0)=d_t_work(j)
2662 dc(j,i)=dc_work(ind+j)
2663 d_t(j,i)=d_t_work(ind+j)
2668 if (itype(i).ne.10) then
2671 dc(j,inres)=dc_work(ind+j)
2672 d_t(j,inres)=d_t_work(ind+j)
2679 c--------------------------------------------------------------------------
2680 subroutine sd_verlet2_ciccotti
2681 c Calculating the adjusted velocities for accelerations
2682 implicit real*8 (a-h,o-z)
2683 include 'DIMENSIONS'
2684 include 'COMMON.CONTROL'
2685 include 'COMMON.VAR'
2688 include 'COMMON.LANGEVIN'
2690 include 'COMMON.LANGEVIN.lang0'
2692 include 'COMMON.CHAIN'
2693 include 'COMMON.DERIV'
2694 include 'COMMON.GEO'
2695 include 'COMMON.LOCAL'
2696 include 'COMMON.INTERACT'
2697 include 'COMMON.IOUNITS'
2698 include 'COMMON.NAMES'
2699 double precision stochforcvec(MAXRES6),stochforcvecV(MAXRES6)
2700 common /stochcalc/ stochforcvec
2702 c Compute the stochastic forces which contribute to velocity change
2704 call stochastic_force(stochforcvecV)
2711 ddt1=ddt1+vfric_mat(i,j)*d_a_work(j)
2712 c ddt2=ddt2+vrand_mat2(i,j)*stochforcvecV(j)
2713 ddt2=ddt2+vrand_mat2(i,j)*stochforcvec(j)
2715 d_t_work(i)=d_t_work_new(i)+0.5d0*ddt1+ddt2
2719 d_t(j,0)=d_t_work(j)
2724 d_t(j,i)=d_t_work(ind+j)
2729 if (itype(i).ne.10) then
2732 d_t(j,inres)=d_t_work(ind+j)
2740 c------------------------------------------------------
2741 double precision function HNose(ek,s,e,pi,Q,t_bath,dimenl)
2743 double precision ek,s,e,pi,Q,t_bath,Rb
2746 HNose=ek+e+pi**2/(2*Q)+dimenl*Rb*t_bath*log(s)
2747 c print '(6f15.5,i5,a2,2f15.5)',ek,s,e,pi,Q,t_bath,dimenl,"--",
2748 c & pi**2/(2*Q),dimenl*Rb*t_bath*log(s)
2751 c-----------------------------------------------------------------
2752 double precision function HNose_nh(eki,e)
2753 implicit real*8 (a-h,o-z)
2754 include 'DIMENSIONS'
2756 HNose_nh=eki+e+dimen3*Rb*t_bath*xlogs(1)+qmass(1)*vlogs(1)**2/2
2758 HNose_nh=HNose_nh+qmass(i)*vlogs(i)**2/2+Rb*t_bath*xlogs(i)
2760 c write(4,'(5e15.5)')
2761 c & vlogs(1),xlogs(1),HNose,eki,e
2764 c-----------------------------------------------------------------
2765 SUBROUTINE NHCINT(akin,scale,wdti,wdti2,wdti4,wdti8)
2766 implicit real*8 (a-h,o-z)
2767 include 'DIMENSIONS'
2769 double precision akin,gnkt,dt,aa,gkt,scale
2770 double precision wdti(maxyosh),wdti2(maxyosh),
2771 & wdti4(maxyosh),wdti8(maxyosh)
2772 integer i,iresn,iyosh,inos,nnos1
2781 C THIS ROUTINE DOES THE NOSE-HOOVER PART OF THE
2782 C INTEGRATION FROM t=0 TO t=DT/2
2783 C GET THE TOTAL KINETIC ENERGY
2785 c CALL GETKINP(MASS,VX,VY,VZ,AKIN)
2787 GLOGS(1) = (AKIN - GNKT)/QMASS(1)
2788 C START THE MULTIPLE TIME STEP PROCEDURE
2791 C UPDATE THE THERMOSTAT VELOCITIES
2792 VLOGS(NNOS) = VLOGS(NNOS) + GLOGS(NNOS)*WDTI4(IYOSH)
2794 AA = EXP(-WDTI8(IYOSH)*VLOGS(NNOS1-INOS) )
2795 VLOGS(NNOS-INOS) = VLOGS(NNOS-INOS)*AA*AA
2796 & + WDTI4(IYOSH)*GLOGS(NNOS-INOS)*AA
2798 C UPDATE THE PARTICLE VELOCITIES
2799 AA = EXP(-WDTI2(IYOSH)*VLOGS(1) )
2802 GLOGS(1) = (SCALE*SCALE*AKIN - GNKT)/QMASS(1)
2803 C UPDATE THE THERMOSTAT POSITIONS
2805 XLOGS(INOS) = XLOGS(INOS) + VLOGS(INOS)*WDTI2(IYOSH)
2807 C UPDATE THE THERMOSTAT VELOCITIES
2809 AA = EXP(-WDTI8(IYOSH)*VLOGS(INOS+1) )
2810 VLOGS(INOS) = VLOGS(INOS)*AA*AA
2811 & + WDTI4(IYOSH)*GLOGS(INOS)*AA
2812 GLOGS(INOS+1) = (QMASS(INOS)*VLOGS(INOS)*VLOGS(INOS)
2813 & -GKT)/QMASS(INOS+1)
2815 VLOGS(NNOS) = VLOGS(NNOS) + GLOGS(NNOS)*WDTI4(IYOSH)
2818 C UPDATE THE PARTICLE VELOCITIES
2819 c outside of this subroutine
2821 c VX(I) = VX(I)*SCALE
2822 c VY(I) = VY(I)*SCALE
2823 c VZ(I) = VZ(I)*SCALE
2827 c-----------------------------------------------------------------
2828 subroutine tnp1_respa_i_step1
2829 c Applying Nose-Poincare algorithm - step 1 to coordinates
2830 c JPSJ 70 75 (2001) S. Nose
2832 c d_t is not updated here
2834 implicit real*8 (a-h,o-z)
2835 include 'DIMENSIONS'
2836 include 'COMMON.CONTROL'
2837 include 'COMMON.VAR'
2839 include 'COMMON.CHAIN'
2840 include 'COMMON.DERIV'
2841 include 'COMMON.GEO'
2842 include 'COMMON.LOCAL'
2843 include 'COMMON.INTERACT'
2844 include 'COMMON.IOUNITS'
2845 include 'COMMON.NAMES'
2846 double precision adt,adt2,tmp
2848 tmp=1+pi_np/(2*Q_np)*0.5*d_time
2851 s12_dt=d_time/s12_np
2852 d_time_s12=d_time*0.5*s12_np
2855 d_t_new(j,0)=d_t_old(j,0)+d_a_old(j,0)*d_time_s12
2856 dc(j,0)=dc_old(j,0)+d_t_new(j,0)*s12_dt
2860 d_t_new(j,i)=d_t_old(j,i)+d_a_old(j,i)*d_time_s12
2861 dc(j,i)=dc_old(j,i)+d_t_new(j,i)*s12_dt
2865 if (itype(i).ne.10) then
2868 d_t_new(j,inres)=d_t_old(j,inres)+d_a_old(j,inres)*d_time_s12
2869 dc(j,inres)=dc_old(j,inres)+d_t_new(j,inres)*s12_dt
2875 c---------------------------------------------------------------------
2876 subroutine tnp1_respa_i_step2
2877 c Step 2 of the velocity Verlet algorithm: update velocities
2878 implicit real*8 (a-h,o-z)
2879 include 'DIMENSIONS'
2880 include 'COMMON.CONTROL'
2881 include 'COMMON.VAR'
2883 include 'COMMON.CHAIN'
2884 include 'COMMON.DERIV'
2885 include 'COMMON.GEO'
2886 include 'COMMON.LOCAL'
2887 include 'COMMON.INTERACT'
2888 include 'COMMON.IOUNITS'
2889 include 'COMMON.NAMES'
2891 double precision d_time_s12
2895 d_t(j,i)=d_t_new(j,i)
2902 d_time_s12=0.5d0*s12_np*d_time
2905 d_t(j,0)=d_t_new(j,0)+d_a(j,0)*d_time_s12
2909 d_t(j,i)=d_t_new(j,i)+d_a(j,i)*d_time_s12
2913 if (itype(i).ne.10) then
2916 d_t(j,inres)=d_t_new(j,inres)+d_a(j,inres)*d_time_s12
2921 pistar=pistar+(EK-0.5*(E_old+potE)
2922 & -dimen3*Rb*t_bath*log(s12_np)+Csplit-dimen3*Rb*t_bath)*d_time
2923 tmp=1+pistar/(2*Q_np)*0.5*d_time
2929 c-------------------------------------------------------
2931 subroutine tnp1_step1
2932 c Applying Nose-Poincare algorithm - step 1 to coordinates
2933 c JPSJ 70 75 (2001) S. Nose
2935 c d_t is not updated here
2937 implicit real*8 (a-h,o-z)
2938 include 'DIMENSIONS'
2939 include 'COMMON.CONTROL'
2940 include 'COMMON.VAR'
2942 include 'COMMON.CHAIN'
2943 include 'COMMON.DERIV'
2944 include 'COMMON.GEO'
2945 include 'COMMON.LOCAL'
2946 include 'COMMON.INTERACT'
2947 include 'COMMON.IOUNITS'
2948 include 'COMMON.NAMES'
2949 double precision adt,adt2,tmp
2951 tmp=1+pi_np/(2*Q_np)*0.5*d_time
2954 s12_dt=d_time/s12_np
2955 d_time_s12=d_time*0.5*s12_np
2958 d_t_new(j,0)=d_t_old(j,0)+d_a_old(j,0)*d_time_s12
2959 dc(j,0)=dc_old(j,0)+d_t_new(j,0)*s12_dt
2963 d_t_new(j,i)=d_t_old(j,i)+d_a_old(j,i)*d_time_s12
2964 dc(j,i)=dc_old(j,i)+d_t_new(j,i)*s12_dt
2968 if (itype(i).ne.10) then
2971 d_t_new(j,inres)=d_t_old(j,inres)+d_a_old(j,inres)*d_time_s12
2972 dc(j,inres)=dc_old(j,inres)+d_t_new(j,inres)*s12_dt
2978 c---------------------------------------------------------------------
2979 subroutine tnp1_step2
2980 c Step 2 of the velocity Verlet algorithm: update velocities
2981 implicit real*8 (a-h,o-z)
2982 include 'DIMENSIONS'
2983 include 'COMMON.CONTROL'
2984 include 'COMMON.VAR'
2986 include 'COMMON.CHAIN'
2987 include 'COMMON.DERIV'
2988 include 'COMMON.GEO'
2989 include 'COMMON.LOCAL'
2990 include 'COMMON.INTERACT'
2991 include 'COMMON.IOUNITS'
2992 include 'COMMON.NAMES'
2994 double precision d_time_s12
2998 d_t(j,i)=d_t_new(j,i)
3005 d_time_s12=0.5d0*s12_np*d_time
3008 d_t(j,0)=d_t_new(j,0)+d_a(j,0)*d_time_s12
3012 d_t(j,i)=d_t_new(j,i)+d_a(j,i)*d_time_s12
3016 if (itype(i).ne.10) then
3019 d_t(j,inres)=d_t_new(j,inres)+d_a(j,inres)*d_time_s12
3024 cd write(iout,*) 'pistar',pistar,EK,E_old,potE,s12_np
3025 pistar=pistar+(EK-0.5*(E_old+potE)
3026 & -dimen3*Rb*t_bath*log(s12_np)+H0-dimen3*Rb*t_bath)*d_time
3027 tmp=1+pistar/(2*Q_np)*0.5*d_time
3034 c-----------------------------------------------------------------
3035 subroutine tnp_respa_i_step1
3036 c Applying Nose-Poincare algorithm - step 1 to coordinates
3037 c J.Comput.Phys. 151 114 (1999) S.D.Bond B.J.Leimkuhler B.B.Laird
3039 c d_t is not updated here, it is destroyed
3041 implicit real*8 (a-h,o-z)
3042 include 'DIMENSIONS'
3043 include 'COMMON.CONTROL'
3044 include 'COMMON.VAR'
3046 include 'COMMON.CHAIN'
3047 include 'COMMON.DERIV'
3048 include 'COMMON.GEO'
3049 include 'COMMON.LOCAL'
3050 include 'COMMON.INTERACT'
3051 include 'COMMON.IOUNITS'
3052 include 'COMMON.NAMES'
3053 double precision C_np,d_time_s,tmp,d_time_ss
3055 d_time_s=d_time*0.5*s_np
3056 ct2 d_time_s=d_time*0.5*s12_np
3059 d_t_new(j,0)=d_t_old(j,0)+d_a_old(j,0)*d_time_s
3063 d_t_new(j,i)=d_t_old(j,i)+d_a_old(j,i)*d_time_s
3067 if (itype(i).ne.10) then
3070 d_t_new(j,inres)=d_t_old(j,inres)+d_a_old(j,inres)*d_time_s
3077 d_t(j,i)=d_t_new(j,i)
3084 C_np=0.5*d_time*(dimen3*Rb*t_bath*(1.0+log(s_np))-EK+potE-Csplit)
3087 pistar=-2.0*C_np/(1.0+sqrt(1.0-C_np*d_time/Q_np))
3088 tmp=0.5*d_time*pistar/Q_np
3089 s12_np=s_np*(1.0+tmp)/(1.0-tmp)
3091 d_time_ss=0.5*d_time*(1.0/s12_np+1.0/s_np)
3092 ct2 d_time_ss=d_time/s12_np
3093 c d_time_ss=0.5*d_time*(1.0/sold_np+1.0/s_np)
3096 dc(j,0)=dc_old(j,0)+d_t_new(j,0)*d_time_ss
3100 dc(j,i)=dc_old(j,i)+d_t_new(j,i)*d_time_ss
3104 if (itype(i).ne.10) then
3107 dc(j,inres)=dc_old(j,inres)+d_t_new(j,inres)*d_time_ss
3114 c---------------------------------------------------------------------
3116 subroutine tnp_respa_i_step2
3117 c Step 2 of the velocity Verlet algorithm: update velocities
3118 implicit real*8 (a-h,o-z)
3119 include 'DIMENSIONS'
3120 include 'COMMON.CONTROL'
3121 include 'COMMON.VAR'
3123 include 'COMMON.CHAIN'
3124 include 'COMMON.DERIV'
3125 include 'COMMON.GEO'
3126 include 'COMMON.LOCAL'
3127 include 'COMMON.INTERACT'
3128 include 'COMMON.IOUNITS'
3129 include 'COMMON.NAMES'
3131 double precision d_time_s
3133 EK=EK*(s_np/s12_np)**2
3134 HNose1=Hnose(EK,s12_np,potE,pistar,Q_np,t_bath,dimen3)
3135 pi_np=pistar+0.5*d_time*(2*EK-dimen3*Rb*t_bath
3138 cr print '(a,5f)','i_step2',EK,potE,HNose1,pi_np,E_long
3139 d_time_s=d_time*0.5*s12_np
3140 c d_time_s=d_time*0.5*s_np
3143 d_t(j,0)=d_t_new(j,0)+d_a(j,0)*d_time_s
3147 d_t(j,i)=d_t_new(j,i)+d_a(j,i)*d_time_s
3151 if (itype(i).ne.10) then
3154 d_t(j,inres)=d_t_new(j,inres)+d_a(j,inres)*d_time_s
3163 c-----------------------------------------------------------------
3164 subroutine tnp_respa_step1
3165 c Applying Nose-Poincare algorithm - step 1 to vel for RESPA
3166 c J.Comput.Phys. 151 114 (1999) S.D.Bond B.J.Leimkuhler B.B.Laird
3168 c d_t is not updated here, it is destroyed
3170 implicit real*8 (a-h,o-z)
3171 include 'DIMENSIONS'
3172 include 'COMMON.CONTROL'
3173 include 'COMMON.VAR'
3175 include 'COMMON.CHAIN'
3176 include 'COMMON.DERIV'
3177 include 'COMMON.GEO'
3178 include 'COMMON.LOCAL'
3179 include 'COMMON.INTERACT'
3180 include 'COMMON.IOUNITS'
3181 include 'COMMON.NAMES'
3182 double precision C_np,d_time_s,tmp,d_time_ss
3183 double precision energia(0:n_ene)
3185 d_time_s=d_time*0.5*s_np
3188 d_t_old(j,0)=d_t_old(j,0)+d_a(j,0)*d_time_s
3192 d_t_old(j,i)=d_t_old(j,i)+d_a(j,i)*d_time_s
3196 if (itype(i).ne.10) then
3199 d_t_old(j,inres)=d_t_old(j,inres)+d_a(j,inres)*d_time_s
3205 c C_np=0.5*d_time*(dimen3*Rb*t_bath*(1.0+log(s_np))-EK+potE-H0)
3208 c pistar=-2.0*C_np/(1.0+sqrt(1.0-C_np*d_time/Q_np))
3209 c tmp=0.5*d_time*pistar/Q_np
3210 c s12_np=s_np*(1.0+tmp)/(1.0-tmp)
3211 c write(iout,*) 'tnp_respa_step1',s_np,s12_np,EK,potE,C_np,pistar,tmp
3217 c-------------------------------------
3218 c test of reviewer's comment
3219 pi_np=pi_np-0.5*d_time*(E_long+Csplit-H0)
3220 cr print '(a,3f)','1 pi_np,s_np',pi_np,s_np,E_long
3221 c-------------------------------------
3225 c---------------------------------------------------------------------
3226 subroutine tnp_respa_step2
3227 c Step 2 of the velocity Verlet algorithm: update velocities for RESPA
3228 implicit real*8 (a-h,o-z)
3229 include 'DIMENSIONS'
3230 include 'COMMON.CONTROL'
3231 include 'COMMON.VAR'
3233 include 'COMMON.CHAIN'
3234 include 'COMMON.DERIV'
3235 include 'COMMON.GEO'
3236 include 'COMMON.LOCAL'
3237 include 'COMMON.INTERACT'
3238 include 'COMMON.IOUNITS'
3239 include 'COMMON.NAMES'
3241 double precision d_time_s
3247 ct HNose1=Hnose(EK,s12_np,potE,pistar,Q_np,t_bath,dimen3)
3248 ct pi_np=pistar+0.5*d_time*(2*EK-dimen3*Rb*t_bath)
3249 ct & -0.5*d_time*(HNose1-H0)
3251 c-------------------------------------
3252 c test of reviewer's comment
3253 pi_np=pi_np-0.5*d_time*(E_long+Csplit-H0)
3254 cr print '(a,3f)','2 pi_np,s_np',pi_np,s_np,E_long
3255 c-------------------------------------
3256 d_time_s=d_time*0.5*s_np
3259 d_t_old(j,0)=d_t_old(j,0)+d_a(j,0)*d_time_s
3263 d_t_old(j,i)=d_t_old(j,i)+d_a(j,i)*d_time_s
3267 if (itype(i).ne.10) then
3270 d_t_old(j,inres)=d_t_old(j,inres)+d_a(j,inres)*d_time_s
3279 c---------------------------------------------------------------------
3280 subroutine tnp_step1
3281 c Applying Nose-Poincare algorithm - step 1 to coordinates
3282 c J.Comput.Phys. 151 114 (1999) S.D.Bond B.J.Leimkuhler B.B.Laird
3284 c d_t is not updated here, it is destroyed
3286 implicit real*8 (a-h,o-z)
3287 include 'DIMENSIONS'
3288 include 'COMMON.CONTROL'
3289 include 'COMMON.VAR'
3291 include 'COMMON.CHAIN'
3292 include 'COMMON.DERIV'
3293 include 'COMMON.GEO'
3294 include 'COMMON.LOCAL'
3295 include 'COMMON.INTERACT'
3296 include 'COMMON.IOUNITS'
3297 include 'COMMON.NAMES'
3298 double precision C_np,d_time_s,tmp,d_time_ss
3300 d_time_s=d_time*0.5*s_np
3303 d_t_new(j,0)=d_t_old(j,0)+d_a_old(j,0)*d_time_s
3307 d_t_new(j,i)=d_t_old(j,i)+d_a_old(j,i)*d_time_s
3311 if (itype(i).ne.10) then
3314 d_t_new(j,inres)=d_t_old(j,inres)+d_a_old(j,inres)*d_time_s
3321 d_t(j,i)=d_t_new(j,i)
3328 C_np=0.5*d_time*(dimen3*Rb*t_bath*(1.0+log(s_np))-EK+potE-H0)
3331 pistar=-2.0*C_np/(1.0+sqrt(1.0-C_np*d_time/Q_np))
3332 tmp=0.5*d_time*pistar/Q_np
3333 s12_np=s_np*(1.0+tmp)/(1.0-tmp)
3334 c write(iout,*) 'tnp_step1',s_np,s12_np,EK,potE,C_np,pistar,tmp
3336 d_time_ss=0.5*d_time*(1.0/s12_np+1.0/s_np)
3339 dc(j,0)=dc_old(j,0)+d_t_new(j,0)*d_time_ss
3343 dc(j,i)=dc_old(j,i)+d_t_new(j,i)*d_time_ss
3347 if (itype(i).ne.10) then
3350 dc(j,inres)=dc_old(j,inres)+d_t_new(j,inres)*d_time_ss
3357 c-----------------------------------------------------------------
3358 subroutine tnp_step2
3359 c Step 2 of the velocity Verlet algorithm: update velocities
3360 implicit real*8 (a-h,o-z)
3361 include 'DIMENSIONS'
3362 include 'COMMON.CONTROL'
3363 include 'COMMON.VAR'
3365 include 'COMMON.CHAIN'
3366 include 'COMMON.DERIV'
3367 include 'COMMON.GEO'
3368 include 'COMMON.LOCAL'
3369 include 'COMMON.INTERACT'
3370 include 'COMMON.IOUNITS'
3371 include 'COMMON.NAMES'
3373 double precision d_time_s
3375 EK=EK*(s_np/s12_np)**2
3376 HNose1=Hnose(EK,s12_np,potE,pistar,Q_np,t_bath,dimen3)
3377 pi_np=pistar+0.5*d_time*(2*EK-dimen3*Rb*t_bath)
3378 & -0.5*d_time*(HNose1-H0)
3380 cd write(iout,'(a,4f)') 'mmm',EK,potE,HNose1,pi_np
3381 d_time_s=d_time*0.5*s12_np
3384 d_t(j,0)=d_t_new(j,0)+d_a(j,0)*d_time_s
3388 d_t(j,i)=d_t_new(j,i)+d_a(j,i)*d_time_s
3392 if (itype(i).ne.10) then
3395 d_t(j,inres)=d_t_new(j,inres)+d_a(j,inres)*d_time_s
3405 subroutine hmc_test(itime)
3406 implicit real*8 (a-h,o-z)
3407 include 'DIMENSIONS'
3408 include 'COMMON.CONTROL'
3410 include 'COMMON.CHAIN'
3413 delta=-(potE+EK-hmc_etot)/(Rb*t_bath)
3414 if (delta .lt. -50.0d0) then
3419 xxx=ran_number(0.0d0,1.0d0)
3421 if (me.eq.king .or. .not. out1file)
3422 & write(iout,'(a8,i5,6f10.4)')
3423 & 'HMC',itime,potE+EK,potE,EK,hmc_etot,delta,xxx
3425 if (delta .le. xxx) then
3434 if (me.eq.king .or. .not. out1file)
3435 & write(iout,*) 'HMC accepting new'
3444 call chainbuild_cart
3448 d_t_old(j,i)=d_t(j,i)
3452 kinetic_T=2.0d0/(dimen3*Rb)*EK
3453 call etotal(potEcomp)
3456 if (me.eq.king .or. .not. out1file)
3457 & write(iout,'(a8,i5,3f10.4)')'HMC new',itime,potE+EK,potE,EK