subroutine MD c------------------------------------------------ c The driver for molecular dynamics subroutines c------------------------------------------------ implicit none include 'DIMENSIONS' #ifdef MPI include "mpif.h" integer IERROR,ERRCODE #endif include 'COMMON.SETUP' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.TIME1' include 'COMMON.HAIRPIN' double precision cm(3),L(3),vcm(3) #ifdef VOUT double precision v_work(maxres6),v_transf(maxres6) #endif integer ilen,rstcount external ilen character*50 tytul common /gucio/ cm integer itime logical ovrtim integer i,j,icount_scale,itime_scal integer nharp,iharp(4,maxres/3) double precision scalfac double precision tt0 c #ifdef MPI if (ilen(tmpdir).gt.0) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_" & //liczba(:ilen(liczba))//'.rst') #else if (ilen(tmpdir).gt.0) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_"//'.rst') #endif write (iout,*) "MD lang",lang t_MDsetup=0.0d0 t_langsetup=0.0d0 t_MD=0.0d0 t_enegrad=0.0d0 t_sdsetup=0.0d0 write (iout,'(20(1h=),a20,20(1h=))') "MD calculation started" #ifdef MPI tt0=MPI_Wtime() #else tt0 = tcpu() #endif c Determine the inverse of the inertia matrix. call setup_MD_matrices c Initialize MD call init_MD #ifdef MPI t_MDsetup = MPI_Wtime()-tt0 #else t_MDsetup = tcpu()-tt0 #endif rstcount=0 c Entering the MD loop #ifdef MPI tt0 = MPI_Wtime() #else tt0 = tcpu() #endif if (lang.eq.2 .or. lang.eq.3) then #ifndef LANG0 call setup_fricmat if (lang.eq.2) then call sd_verlet_p_setup else call sd_verlet_ciccotti_setup endif do i=1,dimen do j=1,dimen pfric0_mat(i,j,0)=pfric_mat(i,j) afric0_mat(i,j,0)=afric_mat(i,j) vfric0_mat(i,j,0)=vfric_mat(i,j) prand0_mat(i,j,0)=prand_mat(i,j) vrand0_mat1(i,j,0)=vrand_mat1(i,j) vrand0_mat2(i,j,0)=vrand_mat2(i,j) enddo enddo flag_stoch(0)=.true. do i=1,maxflag_stoch flag_stoch(i)=.false. enddo #else write (iout,*) & "LANG=2 or 3 NOT SUPPORTED. Recompile without -DLANG0" #ifdef MPI call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE) #endif stop #endif else if (lang.eq.1 .or. lang.eq.4) then call setup_fricmat endif #ifdef MPI t_langsetup=MPI_Wtime()-tt0 tt0=MPI_Wtime() #else t_langsetup=tcpu()-tt0 tt0=tcpu() #endif do itime=1,n_timestep if (ovrtim()) exit if (large.and. mod(itime,ntwe).eq.0) & write (iout,*) "itime",itime rstcount=rstcount+1 if (lang.gt.0 .and. surfarea .and. & mod(itime,reset_fricmat).eq.0) then if (lang.eq.2 .or. lang.eq.3) then #ifndef LANG0 call setup_fricmat if (lang.eq.2) then call sd_verlet_p_setup else call sd_verlet_ciccotti_setup endif do i=1,dimen do j=1,dimen pfric0_mat(i,j,0)=pfric_mat(i,j) afric0_mat(i,j,0)=afric_mat(i,j) vfric0_mat(i,j,0)=vfric_mat(i,j) prand0_mat(i,j,0)=prand_mat(i,j) vrand0_mat1(i,j,0)=vrand_mat1(i,j) vrand0_mat2(i,j,0)=vrand_mat2(i,j) enddo enddo flag_stoch(0)=.true. do i=1,maxflag_stoch flag_stoch(i)=.false. enddo #endif else if (lang.eq.1 .or. lang.eq.4) then call setup_fricmat endif write (iout,'(a,i10)') & "Friction matrix reset based on surface area, itime",itime endif if (reset_vel .and. tbf .and. lang.eq.0 & .and. mod(itime,count_reset_vel).eq.0) then call random_vel write(iout,'(a,f20.2)') & "Velocities reset to random values, time",totT do i=0,2*nres do j=1,3 d_t_old(j,i)=d_t(j,i) enddo enddo endif if (reset_moment .and. mod(itime,count_reset_moment).eq.0) then call inertia_tensor call vcm_vel(vcm) do j=1,3 d_t(j,0)=d_t(j,0)-vcm(j) enddo call kinetic(EK) kinetic_T=2.0d0/(dimen3*Rb)*EK scalfac=dsqrt(T_bath/kinetic_T) write(iout,'(a,f20.2)') "Momenta zeroed out, time",totT do i=0,2*nres do j=1,3 d_t_old(j,i)=scalfac*d_t(j,i) enddo enddo endif if (lang.ne.4) then if (RESPA) then c Time-reversible RESPA algorithm c (Tuckerman et al., J. Chem. Phys., 97, 1990, 1992) call RESPA_step(itime) else c Variable time step algorithm. call velverlet_step(itime) endif else #ifdef BROWN call brown_step(itime) #else print *,"Brown dynamics not here!" #ifdef MPI call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE) #endif stop #endif endif itime_mat=itime if (ntwe.ne.0) then if (mod(itime,ntwe).eq.0) then call statout(itime) C call enerprint(potEcomp) C print *,itime,'AFM',Eafmforc,etot endif #ifdef VOUT do j=1,3 v_work(j)=d_t(j,0) enddo ind=3 do i=nnt,nct-1 do j=1,3 ind=ind+1 v_work(ind)=d_t(j,i) enddo enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then do j=1,3 ind=ind+1 v_work(ind)=d_t(j,i+nres) enddo endif enddo write (66,'(80f10.5)') & ((d_t(j,i),j=1,3),i=0,nres-1),((d_t(j,i+nres),j=1,3),i=1,nres) do i=1,ind v_transf(i)=0.0d0 do j=1,ind v_transf(i)=v_transf(i)+gvec(j,i)*v_work(j) enddo v_transf(i)= v_transf(i)*dsqrt(geigen(i)) enddo write (67,'(80f10.5)') (v_transf(i),i=1,ind) #endif endif if (mod(itime,ntwx).eq.0) then c write(iout,*) 'time=',itime C call check_ecartint call returnbox write (tytul,'("time",f8.2)') totT if(mdpdb) then call hairpin(.true.,nharp,iharp) call secondary2(.true.) call pdbout(potE,tytul,ipdb) else call cartout(totT) endif endif if (rstcount.eq.1000.or.itime.eq.n_timestep) then open(irest2,file=rest2name,status='unknown') write(irest2,*) totT,EK,potE,totE,t_bath do i=0,2*nres-1 write (irest2,'(3e15.5)') (d_t(j,i),j=1,3) enddo do i=0,2*nres-1 write (irest2,'(3e15.5)') (dc(j,i),j=1,3) enddo close(irest2) rstcount=0 endif enddo #ifdef MPI t_MD=MPI_Wtime()-tt0 #else t_MD=tcpu()-tt0 #endif write (iout,'(//35(1h=),a10,35(1h=)/10(/a40,1pe15.5))') & ' Timing ', & 'MD calculations setup:',t_MDsetup, & 'Energy & gradient evaluation:',t_enegrad, & 'Stochastic MD setup:',t_langsetup, & 'Stochastic MD step setup:',t_sdsetup, & 'MD steps:',t_MD write (iout,'(/28(1h=),a25,27(1h=))') & ' End of MD calculation ' #ifdef TIMING_ENE write (iout,*) "time for etotal",t_etotal," elong",t_elong, & " eshort",t_eshort write (iout,*) "time_fric",time_fric," time_stoch",time_stoch, & " time_fricmatmult",time_fricmatmult," time_fsample ", & time_fsample #endif return end c------------------------------------------------------------------------------- subroutine velverlet_step(itime) c------------------------------------------------------------------------------- c Perform a single velocity Verlet step; the time step can be rescaled if c increments in accelerations exceed the threshold c------------------------------------------------------------------------------- implicit none include 'DIMENSIONS' #ifdef MPI include 'mpif.h' integer ierror,ierrcode,errcode #endif include 'COMMON.SETUP' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.TIME1' include 'COMMON.MUCA' double precision vcm(3),incr(3) double precision cm(3),L(3) integer ilen,count,rstcount external ilen character*50 tytul integer maxcount_scale /20/ common /gucio/ cm double precision stochforcvec(MAXRES6) common /stochcalc/ stochforcvec integer itime,icount_scale,itime_scal,ifac_time,i,j,itt logical scale double precision epdrift,fac_time double precision tt0 c scale=.true. icount_scale=0 if (lang.eq.1) then call sddir_precalc else if (lang.eq.2 .or. lang.eq.3) then #ifndef LANG0 call stochastic_force(stochforcvec) #else write (iout,*) & "LANG=2 or 3 NOT SUPPORTED. Recompile without -DLANG0" #ifdef MPI call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE) #endif stop #endif endif itime_scal=0 do while (scale) icount_scale=icount_scale+1 if (icount_scale.gt.maxcount_scale) then write (iout,*) & "ERROR: too many attempts at scaling down the time step. ", & "amax=",amax,"epdrift=",epdrift, & "damax=",damax,"edriftmax=",edriftmax, & "d_time=",d_time call flush(iout) #ifdef MPI call MPI_Abort(MPI_COMM_WORLD,IERROR,IERRCODE) #endif stop endif c First step of the velocity Verlet algorithm if (lang.eq.2) then #ifndef LANG0 call sd_verlet1 #endif else if (lang.eq.3) then #ifndef LANG0 call sd_verlet1_ciccotti #endif else if (lang.eq.1) then call sddir_verlet1 else call verlet1 endif c Build the chain from the newly calculated coordinates call chainbuild_cart if (rattle) call rattle1 if (ntwe.ne.0) then if (large.and. mod(itime,ntwe).eq.0) then write (iout,*) "Cartesian and internal coordinates: step 1" call cartprint call intout write (iout,*) "dC" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(dc(j,i),j=1,3), & (dc(j,i+nres),j=1,3) enddo write (iout,*) "Accelerations" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo write (iout,*) "Velocities, step 1" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif endif #ifdef MPI tt0 = MPI_Wtime() #else tt0 = tcpu() #endif c Calculate energy and forces call zerograd call etotal(potEcomp) ! AL 4/17/17: Reduce the steps if NaNs occurred. if (potEcomp(0).gt.0.99e20 .or. isnan(potEcomp(0)).gt.0) then d_time=d_time/2 cycle endif ! end change if (large.and. mod(itime,ntwe).eq.0) & call enerprint(potEcomp) #ifdef TIMING_ENE #ifdef MPI t_etotal=t_etotal+MPI_Wtime()-tt0 #else t_etotal=t_etotal+tcpu()-tt0 #endif #endif potE=potEcomp(0)-potEcomp(27) call cartgrad c Get the new accelerations call lagrangian #ifdef MPI t_enegrad=t_enegrad+MPI_Wtime()-tt0 #else t_enegrad=t_enegrad+tcpu()-tt0 #endif c Determine maximum acceleration and scale down the timestep if needed call max_accel amax=amax/(itime_scal+1)**2 call predict_edrift(epdrift) if (amax/(itime_scal+1).gt.damax .or. epdrift.gt.edriftmax) then c Maximum acceleration or maximum predicted energy drift exceeded, rescale the time step scale=.true. ifac_time=dmax1(dlog(amax/damax),dlog(epdrift/edriftmax)) & /dlog(2.0d0)+1 itime_scal=itime_scal+ifac_time c fac_time=dmin1(damax/amax,0.5d0) fac_time=0.5d0**ifac_time d_time=d_time*fac_time if (lang.eq.2 .or. lang.eq.3) then #ifndef LANG0 c write (iout,*) "Calling sd_verlet_setup: 1" c Rescale the stochastic forces and recalculate or restore c the matrices of tinker integrator if (itime_scal.gt.maxflag_stoch) then if (large) write (iout,'(a,i5,a)') & "Calculate matrices for stochastic step;", & " itime_scal ",itime_scal if (lang.eq.2) then call sd_verlet_p_setup else call sd_verlet_ciccotti_setup endif write (iout,'(2a,i3,a,i3,1h.)') & "Warning: cannot store matrices for stochastic", & " integration because the index",itime_scal, & " is greater than",maxflag_stoch write (iout,'(2a)')"Increase MAXFLAG_STOCH or use direct", & " integration Langevin algorithm for better efficiency." else if (flag_stoch(itime_scal)) then if (large) write (iout,'(a,i5,a,l1)') & "Restore matrices for stochastic step; itime_scal ", & itime_scal," flag ",flag_stoch(itime_scal) do i=1,dimen do j=1,dimen pfric_mat(i,j)=pfric0_mat(i,j,itime_scal) afric_mat(i,j)=afric0_mat(i,j,itime_scal) vfric_mat(i,j)=vfric0_mat(i,j,itime_scal) prand_mat(i,j)=prand0_mat(i,j,itime_scal) vrand_mat1(i,j)=vrand0_mat1(i,j,itime_scal) vrand_mat2(i,j)=vrand0_mat2(i,j,itime_scal) enddo enddo else if (large) write (iout,'(2a,i5,a,l1)') & "Calculate & store matrices for stochastic step;", & " itime_scal ",itime_scal," flag ",flag_stoch(itime_scal) if (lang.eq.2) then call sd_verlet_p_setup else call sd_verlet_ciccotti_setup endif flag_stoch(ifac_time)=.true. do i=1,dimen do j=1,dimen pfric0_mat(i,j,itime_scal)=pfric_mat(i,j) afric0_mat(i,j,itime_scal)=afric_mat(i,j) vfric0_mat(i,j,itime_scal)=vfric_mat(i,j) prand0_mat(i,j,itime_scal)=prand_mat(i,j) vrand0_mat1(i,j,itime_scal)=vrand_mat1(i,j) vrand0_mat2(i,j,itime_scal)=vrand_mat2(i,j) enddo enddo endif fac_time=1.0d0/dsqrt(fac_time) do i=1,dimen stochforcvec(i)=fac_time*stochforcvec(i) enddo #endif else if (lang.eq.1) then c Rescale the accelerations due to stochastic forces fac_time=1.0d0/dsqrt(fac_time) do i=1,dimen d_as_work(i)=d_as_work(i)*fac_time enddo endif if (large) write (iout,'(a,i10,a,f8.6,a,i3,a,i3)') & "itime",itime," Timestep scaled down to ", & d_time," ifac_time",ifac_time," itime_scal",itime_scal else c Second step of the velocity Verlet algorithm if (lang.eq.2) then #ifndef LANG0 call sd_verlet2 #endif else if (lang.eq.3) then #ifndef LANG0 call sd_verlet2_ciccotti #endif else if (lang.eq.1) then call sddir_verlet2 else call verlet2 endif if (rattle) call rattle2 totT=totT+d_time totTafm=totT C print *,totTafm,"TU?" if (d_time.ne.d_time0) then d_time=d_time0 #ifndef LANG0 if (lang.eq.2 .or. lang.eq.3) then if (large) write (iout,'(a)') & "Restore original matrices for stochastic step" c write (iout,*) "Calling sd_verlet_setup: 2" c Restore the matrices of tinker integrator if the time step has been restored do i=1,dimen do j=1,dimen pfric_mat(i,j)=pfric0_mat(i,j,0) afric_mat(i,j)=afric0_mat(i,j,0) vfric_mat(i,j)=vfric0_mat(i,j,0) prand_mat(i,j)=prand0_mat(i,j,0) vrand_mat1(i,j)=vrand0_mat1(i,j,0) vrand_mat2(i,j)=vrand0_mat2(i,j,0) enddo enddo endif #endif endif scale=.false. endif enddo c Calculate the kinetic and the total energy and the kinetic temperature call kinetic(EK) totE=EK+potE c diagnostics c call kinetic1(EK1) c write (iout,*) "step",itime," EK",EK," EK1",EK1 c end diagnostics c Couple the system to Berendsen bath if needed if (tbf .and. lang.eq.0) then call verlet_bath endif kinetic_T=2.0d0/(dimen3*Rb)*EK c Backup the coordinates, velocities, and accelerations do i=0,2*nres do j=1,3 dc_old(j,i)=dc(j,i) d_t_old(j,i)=d_t(j,i) d_a_old(j,i)=d_a(j,i) enddo enddo if (ntwe.ne.0) then if (mod(itime,ntwe).eq.0 .and. large) then write (iout,*) "Velocities, step 2" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif endif return end c------------------------------------------------------------------------------- subroutine RESPA_step(itime) c------------------------------------------------------------------------------- c Perform a single RESPA step. c------------------------------------------------------------------------------- implicit none include 'DIMENSIONS' #ifdef MPI include 'mpif.h' integer IERROR,ERRCODE #endif include 'COMMON.SETUP' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.TIME1' logical lprint_short common /shortcheck/ lprint_short double precision energia_short(0:n_ene), & energia_long(0:n_ene) double precision cm(3),L(3),vcm(3),incr(3) double precision dc_old0(3,0:maxres2),d_t_old0(3,0:maxres2), & d_a_old0(3,0:maxres2) integer i,j double precision fac_time logical PRINT_AMTS_MSG /.false./ integer ilen,count,rstcount external ilen character*50 tytul integer maxcount_scale /10/ common /gucio/ cm double precision stochforcvec(MAXRES6) common /stochcalc/ stochforcvec integer itime logical scale integer itt common /cipiszcze/ itt integer itsplit double precision epdrift,epdriftmax double precision tt0 itt=itime if (ntwe.ne.0) then if (large.and. mod(itime,ntwe).eq.0) then write (iout,*) "***************** RESPA itime",itime write (iout,*) "Cartesian and internal coordinates: step 0" c call cartprint call pdbout(0.0d0, & "cipiszcze ",iout) call intout write (iout,*) "Accelerations from long-range forces" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo write (iout,*) "Velocities, step 0" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif endif c c Perform the initial RESPA step (increment velocities) c write (iout,*) "*********************** RESPA ini" call RESPA_vel if (ntwe.ne.0) then if (mod(itime,ntwe).eq.0 .and. large) then write (iout,*) "Velocities, end" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif endif c Compute the short-range forces #ifdef MPI tt0 =MPI_Wtime() #else tt0 = tcpu() #endif C 7/2/2009 commented out c call zerograd c call etotal_short(energia_short) c call cartgrad c call lagrangian C 7/2/2009 Copy accelerations due to short-lange forces from previous MD step do i=0,2*nres do j=1,3 d_a(j,i)=d_a_short(j,i) enddo enddo if (ntwe.ne.0) then if (large.and. mod(itime,ntwe).eq.0) then write (iout,*) "energia_short",energia_short(0) write (iout,*) "Accelerations from short-range forces" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo endif endif #ifdef MPI t_enegrad=t_enegrad+MPI_Wtime()-tt0 #else t_enegrad=t_enegrad+tcpu()-tt0 #endif do i=0,2*nres do j=1,3 dc_old(j,i)=dc(j,i) d_t_old(j,i)=d_t(j,i) d_a_old(j,i)=d_a(j,i) enddo enddo c 6/30/08 A-MTS: attempt at increasing the split number do i=0,2*nres do j=1,3 dc_old0(j,i)=dc_old(j,i) d_t_old0(j,i)=d_t_old(j,i) d_a_old0(j,i)=d_a_old(j,i) enddo enddo if (ntime_split.gt.ntime_split0) ntime_split=ntime_split/2 if (ntime_split.lt.ntime_split0) ntime_split=ntime_split0 c scale=.true. d_time0=d_time do while (scale) scale=.false. c write (iout,*) "itime",itime," ntime_split",ntime_split c Split the time step d_time=d_time0/ntime_split c Perform the short-range RESPA steps (velocity Verlet increments of c positions and velocities using short-range forces) c write (iout,*) "*********************** RESPA split" do itsplit=1,ntime_split if (lang.eq.1) then call sddir_precalc else if (lang.eq.2 .or. lang.eq.3) then #ifndef LANG0 call stochastic_force(stochforcvec) #else write (iout,*) & "LANG=2 or 3 NOT SUPPORTED. Recompile without -DLANG0" #ifdef MPI call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE) #endif stop #endif endif c First step of the velocity Verlet algorithm if (lang.eq.2) then #ifndef LANG0 call sd_verlet1 #endif else if (lang.eq.3) then #ifndef LANG0 call sd_verlet1_ciccotti #endif else if (lang.eq.1) then call sddir_verlet1 else call verlet1 endif c Build the chain from the newly calculated coordinates call chainbuild_cart if (rattle) call rattle1 if (ntwe.ne.0) then if (large.and. mod(itime,ntwe).eq.0) then write (iout,*) "***** ITSPLIT",itsplit write (iout,*) "Cartesian and internal coordinates: step 1" call pdbout(0.0d0, & "cipiszcze ",iout) c call cartprint call intout write (iout,*) "Velocities, step 1" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif endif #ifdef MPI tt0 = MPI_Wtime() #else tt0 = tcpu() #endif c Calculate energy and forces c if (large.and. mod(itime,ntwe).eq.0) lprint_short=.true. call zerograd call etotal_short(energia_short) ! AL 4/17/17: Exit itime_split loop when energy goes infinite if (energia_short(0).gt.0.99e20 .or. isnan(energia_short(0)) ) & then if (PRINT_AMTS_MSG) write (iout,*) & "Infinities/NaNs in energia_short", & energia_short(0),"; increasing ntime_split to",ntime_split ntime_split=ntime_split*2 if (ntime_split.gt.maxtime_split) then #ifdef MPI write (iout,*) "Cannot rescue the run; aborting job.", & " Retry with a smaller time step" call flush(iout) call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE) #else write (iout,*) "Cannot rescue the run; terminating.", & " Retry with a smaller time step" #endif endif exit endif ! End change if (large.and. mod(itime,ntwe).eq.0) & call enerprint(energia_short) #ifdef TIMING_ENE #ifdef MPI t_eshort=t_eshort+MPI_Wtime()-tt0 #else t_eshort=t_eshort+tcpu()-tt0 #endif #endif call cartgrad c Get the new accelerations call lagrangian C 7/2/2009 Copy accelerations due to short-lange forces to an auxiliary array do i=0,2*nres do j=1,3 d_a_short(j,i)=d_a(j,i) enddo enddo if (ntwe.ne.0) then if (large.and. mod(itime,ntwe).eq.0) then write (iout,*)"energia_short",energia_short(0) write (iout,*) "Accelerations from short-range forces" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo endif endif c 6/30/08 A-MTS c Determine maximum acceleration and scale down the timestep if needed call max_accel amax=amax/ntime_split**2 call predict_edrift(epdrift) if (ntwe.gt.0 .and. large .and. mod(itime,ntwe).eq.0) & write (iout,*) "amax",amax," damax",damax, & " epdrift",epdrift," epdriftmax",epdriftmax c Exit loop and try with increased split number if the change of c acceleration is too big if (amax.gt.damax .or. epdrift.gt.edriftmax) then if (ntime_split.lt.maxtime_split) then scale=.true. ntime_split=ntime_split*2 do i=0,2*nres do j=1,3 dc_old(j,i)=dc_old0(j,i) d_t_old(j,i)=d_t_old0(j,i) d_a_old(j,i)=d_a_old0(j,i) enddo enddo if (PRINT_AMTS_MSG) then write (iout,*) "acceleration/energy drift too large",amax, & epdrift," split increased to ",ntime_split," itime",itime, & " itsplit",itsplit endif exit else write (iout,*) & "Uh-hu. Bumpy landscape. Maximum splitting number", & maxtime_split, & " already reached!!! Trying to carry on!" endif endif #ifdef MPI t_enegrad=t_enegrad+MPI_Wtime()-tt0 #else t_enegrad=t_enegrad+tcpu()-tt0 #endif c Second step of the velocity Verlet algorithm if (lang.eq.2) then #ifndef LANG0 call sd_verlet2 #endif else if (lang.eq.3) then #ifndef LANG0 call sd_verlet2_ciccotti #endif else if (lang.eq.1) then call sddir_verlet2 else call verlet2 endif if (rattle) call rattle2 c Backup the coordinates, velocities, and accelerations do i=0,2*nres do j=1,3 dc_old(j,i)=dc(j,i) d_t_old(j,i)=d_t(j,i) d_a_old(j,i)=d_a(j,i) enddo enddo enddo enddo ! while scale c Restore the time step d_time=d_time0 c Compute long-range forces #ifdef MPI tt0 =MPI_Wtime() #else tt0 = tcpu() #endif call zerograd call etotal_long(energia_long) ! AL 4/17/2017 Handling NaNs if (energia_long(0).gt.0.99e20 .or. isnan(energia_long(0))) then #ifdef MPI write (iout,*) & "Infinitied/NaNs in energia_long, Aborting MPI job." call enerprint(energia_long) call flush(iout) call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE) #else write (iout,*) "Infinitied/NaNs in energia_long, terminating." call enerprint(energia_long) stop #endif endif ! end change c lprint_short=.false. if (large.and. mod(itime,ntwe).eq.0) & call enerprint(energia_long) #ifdef TIMING_ENE #ifdef MPI t_elong=t_elong+MPI_Wtime()-tt0 #else t_elong=t_elong+tcpu()-tt0 #endif #endif call cartgrad call lagrangian #ifdef MPI t_enegrad=t_enegrad+MPI_Wtime()-tt0 #else t_enegrad=t_enegrad+tcpu()-tt0 #endif c Compute accelerations from long-range forces if (ntwe.ne.0) then if (large.and. mod(itime,ntwe).eq.0) then write (iout,*) "energia_long",energia_long(0) write (iout,*) "Cartesian and internal coordinates: step 2" c call cartprint call pdbout(0.0d0, & 'cipiszcze ',iout) call intout write (iout,*) "Accelerations from long-range forces" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo write (iout,*) "Velocities, step 2" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif endif c Compute the final RESPA step (increment velocities) c write (iout,*) "*********************** RESPA fin" call RESPA_vel c Compute the complete potential energy do i=0,n_ene potEcomp(i)=energia_short(i)+energia_long(i) enddo potE=potEcomp(0)-potEcomp(27) if (ntwe.ne.0) then if (large.and. mod(itime,ntwe).eq.0) then call enerprint(potEcomp) write (iout,*) "potE",potE endif endif c potE=energia_short(0)+energia_long(0) totT=totT+d_time totTafm=totT c Calculate the kinetic and the total energy and the kinetic temperature call kinetic(EK) totE=EK+potE c Couple the system to Berendsen bath if needed if (tbf .and. lang.eq.0) then call verlet_bath endif kinetic_T=2.0d0/(dimen3*Rb)*EK c Backup the coordinates, velocities, and accelerations if (ntwe.ne.0) then if (mod(itime,ntwe).eq.0 .and. large) then write (iout,*) "Velocities, end" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif endif return end c--------------------------------------------------------------------- subroutine RESPA_vel c First and last RESPA step (incrementing velocities using long-range c forces). implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' integer i,j,inres do j=1,3 d_t(j,0)=d_t(j,0)+0.5d0*d_a(j,0)*d_time enddo do i=nnt,nct-1 do j=1,3 d_t(j,i)=d_t(j,i)+0.5d0*d_a(j,i)*d_time enddo enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 d_t(j,inres)=d_t(j,inres)+0.5d0*d_a(j,inres)*d_time enddo endif enddo return end c----------------------------------------------------------------- subroutine verlet1 c Applying velocity Verlet algorithm - step 1 to coordinates implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' double precision adt,adt2 integer i,j,inres #ifdef DEBUG write (iout,*) "VELVERLET1 START: DC" do i=0,nres write (iout,'(i3,3f10.5,5x,3f10.5)') i,(dc(j,i),j=1,3), & (dc(j,i+nres),j=1,3) enddo #endif do j=1,3 adt=d_a_old(j,0)*d_time adt2=0.5d0*adt dc(j,0)=dc_old(j,0)+(d_t_old(j,0)+adt2)*d_time d_t_new(j,0)=d_t_old(j,0)+adt2 d_t(j,0)=d_t_old(j,0)+adt enddo do i=nnt,nct-1 C SPYTAC ADAMA C do i=0,nres #ifdef DEBUG write (iout,*) "i",i," d_a_old",(d_a_old(j,i),j=1,3) #endif do j=1,3 adt=d_a_old(j,i)*d_time adt2=0.5d0*adt dc(j,i)=dc_old(j,i)+(d_t_old(j,i)+adt2)*d_time d_t_new(j,i)=d_t_old(j,i)+adt2 d_t(j,i)=d_t_old(j,i)+adt enddo enddo do i=nnt,nct C do i=0,nres if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 adt=d_a_old(j,inres)*d_time adt2=0.5d0*adt dc(j,inres)=dc_old(j,inres)+(d_t_old(j,inres)+adt2)*d_time d_t_new(j,inres)=d_t_old(j,inres)+adt2 d_t(j,inres)=d_t_old(j,inres)+adt enddo endif enddo #ifdef DEBUG write (iout,*) "VELVERLET1 END: DC" do i=0,nres write (iout,'(i3,3f10.5,5x,3f10.5)') i,(dc(j,i),j=1,3), & (dc(j,i+nres),j=1,3) enddo #endif return end c--------------------------------------------------------------------- subroutine verlet2 c Step 2 of the velocity Verlet algorithm: update velocities implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' integer i,j,inres do j=1,3 d_t(j,0)=d_t_new(j,0)+0.5d0*d_a(j,0)*d_time enddo do i=nnt,nct-1 do j=1,3 d_t(j,i)=d_t_new(j,i)+0.5d0*d_a(j,i)*d_time enddo enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 d_t(j,inres)=d_t_new(j,inres)+0.5d0*d_a(j,inres)*d_time enddo endif enddo return end c----------------------------------------------------------------- subroutine sddir_precalc c Applying velocity Verlet algorithm - step 1 to coordinates implicit none include 'DIMENSIONS' #ifdef MPI include 'mpif.h' #endif include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.TIME1' double precision time00 double precision stochforcvec(MAXRES6) common /stochcalc/ stochforcvec integer i c c Compute friction and stochastic forces c #ifdef MPI time00=MPI_Wtime() #else time00=tcpu() #endif call friction_force c write (iout,*) "After friction_force" #ifdef MPI time_fric=time_fric+MPI_Wtime()-time00 time00=MPI_Wtime() #else time_fric=time_fric+tcpu()-time00 time00=tcpu() #endif call stochastic_force(stochforcvec) c write (iout,*) "After stochastic_force" #ifdef MPI time_stoch=time_stoch+MPI_Wtime()-time00 #else time_stoch=time_stoch+tcpu()-time00 #endif c c Compute the acceleration due to friction forces (d_af_work) and stochastic c forces (d_as_work) c #ifdef FIVEDIAG c write (iout,*) "friction accelerations" call fivediaginv_mult(dimen,fric_work, d_af_work) c write (iout,*) "stochastic acceleratios" call fivediaginv_mult(dimen,stochforcvec, d_as_work) #else call ginv_mult(fric_work, d_af_work) call ginv_mult(stochforcvec, d_as_work) #endif #ifdef DEBUG write (iout,*) "d_af_work" write (iout,'(3f10.5)') (d_af_work(i),i=1,dimen3) write (iout,*) "d_as_work" write (iout,'(3f10.5)') (d_as_work(i),i=1,dimen3) write (iout,*) "Leaving sddir_precalc" #endif return end c--------------------------------------------------------------------- subroutine sddir_verlet1 c Applying velocity Verlet algorithm - step 1 to velocities implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' c Revised 3/31/05 AL: correlation between random contributions to c position and velocity increments included. double precision sqrt13 /0.57735026918962576451d0/ ! 1/sqrt(3) double precision adt,adt2 integer i,j,ind,inres c c Add the contribution from BOTH friction and stochastic force to the c coordinates, but ONLY the contribution from the friction forces to velocities c do j=1,3 adt=(d_a_old(j,0)+d_af_work(j))*d_time adt2=0.5d0*adt+sqrt13*d_as_work(j)*d_time dc(j,0)=dc_old(j,0)+(d_t_old(j,0)+adt2)*d_time d_t_new(j,0)=d_t_old(j,0)+0.5d0*adt d_t(j,0)=d_t_old(j,0)+adt enddo ind=3 do i=nnt,nct-1 do j=1,3 adt=(d_a_old(j,i)+d_af_work(ind+j))*d_time adt2=0.5d0*adt+sqrt13*d_as_work(ind+j)*d_time dc(j,i)=dc_old(j,i)+(d_t_old(j,i)+adt2)*d_time d_t_new(j,i)=d_t_old(j,i)+0.5d0*adt d_t(j,i)=d_t_old(j,i)+adt enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 adt=(d_a_old(j,inres)+d_af_work(ind+j))*d_time adt2=0.5d0*adt+sqrt13*d_as_work(ind+j)*d_time dc(j,inres)=dc_old(j,inres)+(d_t_old(j,inres)+adt2)*d_time d_t_new(j,inres)=d_t_old(j,inres)+0.5d0*adt d_t(j,inres)=d_t_old(j,inres)+adt enddo ind=ind+3 endif enddo return end c--------------------------------------------------------------------- subroutine sddir_verlet2 c Calculating the adjusted velocities for accelerations implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' double precision stochforcvec(MAXRES6),d_as_work1(MAXRES6) double precision cos60 /0.5d0/, sin60 /0.86602540378443864676d0/ integer i,j,inres,ind c Revised 3/31/05 AL: correlation between random contributions to c position and velocity increments included. c The correlation coefficients are calculated at low-friction limit. c Also, friction forces are now not calculated with new velocities. c call friction_force call stochastic_force(stochforcvec) c c Compute the acceleration due to friction forces (d_af_work) and stochastic c forces (d_as_work) c #ifdef FIVEDIAG call fivediaginv_mult(maxres6,stochforcvec, d_as_work1) #else call ginv_mult(stochforcvec, d_as_work1) #endif c c Update velocities c do j=1,3 d_t(j,0)=d_t_new(j,0)+(0.5d0*(d_a(j,0)+d_af_work(j)) & +sin60*d_as_work(j)+cos60*d_as_work1(j))*d_time enddo ind=3 do i=nnt,nct-1 do j=1,3 d_t(j,i)=d_t_new(j,i)+(0.5d0*(d_a(j,i)+d_af_work(ind+j)) & +sin60*d_as_work(ind+j)+cos60*d_as_work1(ind+j))*d_time enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 d_t(j,inres)=d_t_new(j,inres)+(0.5d0*(d_a(j,inres) & +d_af_work(ind+j))+sin60*d_as_work(ind+j) & +cos60*d_as_work1(ind+j))*d_time enddo ind=ind+3 endif enddo return end c--------------------------------------------------------------------- subroutine max_accel c c Find the maximum difference in the accelerations of the the sites c at the beginning and the end of the time step. c implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' integer i,j double precision aux(3),accel(3),accel_old(3),dacc do j=1,3 c aux(j)=d_a(j,0)-d_a_old(j,0) accel_old(j)=d_a_old(j,0) accel(j)=d_a(j,0) enddo amax=0.0d0 do i=nnt,nct c Backbone if (i.lt.nct) then c 7/3/08 changed to asymmetric difference do j=1,3 c accel(j)=aux(j)+0.5d0*(d_a(j,i)-d_a_old(j,i)) accel_old(j)=accel_old(j)+0.5d0*d_a_old(j,i) accel(j)=accel(j)+0.5d0*d_a(j,i) c if (dabs(accel(j)).gt.amax) amax=dabs(accel(j)) if (dabs(accel(j)).gt.dabs(accel_old(j))) then dacc=dabs(accel(j)-accel_old(j)) c write (iout,*) i,dacc if (dacc.gt.amax) amax=dacc endif enddo endif enddo c Side chains do j=1,3 c accel(j)=aux(j) accel_old(j)=d_a_old(j,0) accel(j)=d_a(j,0) enddo if (nnt.eq.2) then do j=1,3 accel_old(j)=accel_old(j)+d_a_old(j,1) accel(j)=accel(j)+d_a(j,1) enddo endif do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then do j=1,3 c accel(j)=accel(j)+d_a(j,i+nres)-d_a_old(j,i+nres) accel_old(j)=accel_old(j)+d_a_old(j,i+nres) accel(j)=accel(j)+d_a(j,i+nres) enddo endif do j=1,3 c if (dabs(accel(j)).gt.amax) amax=dabs(accel(j)) if (dabs(accel(j)).gt.dabs(accel_old(j))) then dacc=dabs(accel(j)-accel_old(j)) c write (iout,*) "side-chain",i,dacc if (dacc.gt.amax) amax=dacc endif enddo do j=1,3 accel_old(j)=accel_old(j)+d_a_old(j,i) accel(j)=accel(j)+d_a(j,i) c aux(j)=aux(j)+d_a(j,i)-d_a_old(j,i) enddo enddo return end c--------------------------------------------------------------------- subroutine predict_edrift(epdrift) c c Predict the drift of the potential energy c implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.MUCA' double precision epdrift,epdriftij integer i,j c Drift of the potential energy epdrift=0.0d0 do i=nnt,nct c Backbone if (i.lt.nct) then do j=1,3 epdriftij=dabs((d_a(j,i)-d_a_old(j,i))*gcart(j,i)) if (lmuca) epdriftij=epdriftij*factor c write (iout,*) "back",i,j,epdriftij if (epdriftij.gt.epdrift) epdrift=epdriftij enddo endif c Side chains if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then do j=1,3 epdriftij= & dabs((d_a(j,i+nres)-d_a_old(j,i+nres))*gxcart(j,i)) if (lmuca) epdriftij=epdriftij*factor c write (iout,*) "side",i,j,epdriftij if (epdriftij.gt.epdrift) epdrift=epdriftij enddo endif enddo epdrift=0.5d0*epdrift*d_time*d_time c write (iout,*) "epdrift",epdrift return end c----------------------------------------------------------------------- subroutine verlet_bath c c Coupling to the thermostat by using the Berendsen algorithm c implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifdef LANG0 #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #else include 'COMMON.LANGEVIN' #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' double precision T_half,fact integer i,j ,inres T_half=2.0d0/(dimen3*Rb)*EK fact=dsqrt(1.0d0+(d_time/tau_bath)*(t_bath/T_half-1.0d0)) c write(iout,*) "T_half", T_half c write(iout,*) "EK", EK c write(iout,*) "fact", fact do j=1,3 d_t(j,0)=fact*d_t(j,0) enddo do i=nnt,nct-1 do j=1,3 d_t(j,i)=fact*d_t(j,i) enddo enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 d_t(j,inres)=fact*d_t(j,inres) enddo endif enddo return end c--------------------------------------------------------- subroutine init_MD c Set up the initial conditions of a MD simulation implicit none include 'DIMENSIONS' #ifdef MP include 'mpif.h' character*16 form integer IERROR,ERRCODE,error_msg,ierr,ierrcode #endif include 'COMMON.SETUP' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.QRESTR' #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.REMD' include 'COMMON.TIME1' #ifdef LBFGS character*9 status integer niter common /lbfgstat/ status,niter,nfun #endif integer n_model_try,list_model_try(max_template),k double precision tt0 real*8 energia_long(0:n_ene), & energia_short(0:n_ene),vcm(3),incr(3) double precision cm(3),L(3),xv,sigv,lowb,highb double precision varia(maxvar),energia(0:n_ene) character*256 qstr integer ilen external ilen character*50 tytul logical file_exist common /gucio/ cm integer i,ipos,iq,iw,j,iranmin,nft_sc,iretcode,nfun,itrial,itmp, & i_model,itime integer iran_num double precision etot logical fail integer i_start_models(0:nodes-1) write (iout,*) "init_MD INDPDB",indpdb d_time0=d_time c write(iout,*) "d_time", d_time c Compute the standard deviations of stochastic forces for Langevin dynamics c if the friction coefficients do not depend on surface area if (lang.gt.0 .and. .not.surfarea) then do i=nnt,nct-1 stdforcp(i)=stdfp*dsqrt(gamp) enddo do i=nnt,nct if (itype(i).ne.ntyp1) stdforcsc(i)=stdfsc(iabs(itype(i))) & *dsqrt(gamsc(iabs(itype(i)))) enddo endif c Open the pdb file for snapshotshots #ifdef MPI if(mdpdb) then if (ilen(tmpdir).gt.0) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD"// & liczba(:ilen(liczba))//".pdb") open(ipdb, & file=prefix(:ilen(prefix))//"_MD"//liczba(:ilen(liczba)) & //".pdb") else #ifdef NOXDR if (ilen(tmpdir).gt.0 .and. (me.eq.king .or. .not.traj1file)) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD"// & liczba(:ilen(liczba))//".x") cartname=prefix(:ilen(prefix))//"_MD"//liczba(:ilen(liczba)) & //".x" #else if (ilen(tmpdir).gt.0 .and. (me.eq.king .or. .not.traj1file)) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD"// & liczba(:ilen(liczba))//".cx") cartname=prefix(:ilen(prefix))//"_MD"//liczba(:ilen(liczba)) & //".cx" #endif endif #else if(mdpdb) then if (ilen(tmpdir).gt.0) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD.pdb") open(ipdb,file=prefix(:ilen(prefix))//"_MD.pdb") else if (ilen(tmpdir).gt.0) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_MD.cx") cartname=prefix(:ilen(prefix))//"_MD.cx" endif #endif if (usampl) then write (qstr,'(256(1h ))') ipos=1 do i=1,nfrag iq = qinfrag(i,iset)*10 iw = wfrag(i,iset)/100 if (iw.gt.0) then if(me.eq.king.or..not.out1file) & write (iout,*) "Frag",qinfrag(i,iset),wfrag(i,iset),iq,iw write (qstr(ipos:ipos+6),'(2h_f,i1,1h_,i1,1h_,i1)') i,iq,iw ipos=ipos+7 endif enddo do i=1,npair iq = qinpair(i,iset)*10 iw = wpair(i,iset)/100 if (iw.gt.0) then if(me.eq.king.or..not.out1file) & write (iout,*) "Pair",i,qinpair(i,iset),wpair(i,iset),iq,iw write (qstr(ipos:ipos+6),'(2h_p,i1,1h_,i1,1h_,i1)') i,iq,iw ipos=ipos+7 endif enddo c pdbname=pdbname(:ilen(pdbname)-4)//qstr(:ipos-1)//'.pdb' #ifdef NOXDR c cartname=cartname(:ilen(cartname)-2)//qstr(:ipos-1)//'.x' #else c cartname=cartname(:ilen(cartname)-3)//qstr(:ipos-1)//'.cx' #endif c statname=statname(:ilen(statname)-5)//qstr(:ipos-1)//'.stat' endif icg=1 write (iout,*) "REST ",rest if (rest) then if (restart1file) then if (me.eq.king) & inquire(file=mremd_rst_name,exist=file_exist) #ifdef MPI c write (*,*) me," Before broadcast: file_exist",file_exist call MPI_Bcast(file_exist,1,MPI_LOGICAL,king,CG_COMM, & IERR) c write (*,*) me," After broadcast: file_exist",file_exist c inquire(file=mremd_rst_name,exist=file_exist) #endif if(me.eq.king.or..not.out1file) & write(iout,*) "Initial state read by master and distributed" else if (ilen(tmpdir).gt.0) & call copy_to_tmp(pref_orig(:ilen(pref_orig))//'_' & //liczba(:ilen(liczba))//'.rst') inquire(file=rest2name,exist=file_exist) endif if(file_exist) then if(.not.restart1file) then if(me.eq.king.or..not.out1file) & write(iout,*) "Initial state will be read from file ", & rest2name(:ilen(rest2name)) call readrst endif call rescale_weights(t_bath) else rest=.false. if(me.eq.king.or..not.out1file)then if (restart1file) then write(iout,*) "File ",mremd_rst_name(:ilen(mremd_rst_name)), & " does not exist" else write(iout,*) "File ",rest2name(:ilen(rest2name)), & " does not exist" endif write(iout,*) "Initial velocities randomly generated" endif call random_vel totT=0.0d0 totTafm=totT endif else c Generate initial velocities if(me.eq.king.or..not.out1file) & write(iout,*) "Initial velocities randomly generated" call random_vel totT=0.0d0 CtotTafm is the variable for AFM time which eclipsed during totTafm=totT endif c rest2name = prefix(:ilen(prefix))//'.rst' if(me.eq.king.or..not.out1file)then write (iout,*) "Initial velocities" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo endif if (lang.eq.0) then c Zeroing the total angular momentum of the system write(iout,*) "Calling the zero-angular & momentum subroutine" call inertia_tensor c Getting the potential energy and forces and velocities and accelerations call vcm_vel(vcm) write (iout,*) "velocity of the center of the mass:" write (iout,*) (vcm(j),j=1,3) call flush(iout) do j=1,3 d_t(j,0)=d_t(j,0)-vcm(j) enddo c Removing the velocity of the center of mass call vcm_vel(vcm) if(me.eq.king.or..not.out1file)then write (iout,*) "vcm right after adjustment:" write (iout,*) (vcm(j),j=1,3) write (iout,*) "Initial velocities after adjustment" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo call flush(iout) endif endif c write (iout,*) "init_MD before initial structure REST ",rest if(start_from_model .and. (me.eq.king .or. .not. out1file)) & write(iout,*) 'START_FROM_MODELS is ON' if(start_from_model .and. rest) then if(me.eq.king .or. .not. out1file .and. fg_rank.eq.0) then write(iout,*) & 'START_FROM_MODELS is OFF because the run is restarted' write(iout,*) 'Remove restart keyword from input' endif endif c write (iout,*) "rest ",rest," start_from_model",start_from_model, c & " nmodel_start",nmodel_start," preminim",preminim if (.not.rest) then 122 continue if (iranconf.ne.0) then c 8/22/17 AL Loop to produce a low-energy random conformation do iranmin=1,10 call chainbuild if (overlapsc) then print *, 'Calling OVERLAP_SC' call overlap_sc(fail) endif if (searchsc) then call sc_move(2,nres-1,10,1d10,nft_sc,etot) print *,'SC_move',nft_sc,etot if(me.eq.king.or..not.out1file) & write(iout,*) 'SC_move',nft_sc,etot endif if (dccart) then if (me.eq.king.or..not.out1file) write(iout,*) & 'Minimizing random structure: Calling MINIM_DC' call minim_dc(etot,iretcode,nfun) else call geom_to_var(nvar,varia) if (me.eq.king.or..not.out1file) write (iout,*) & 'Minimizing random structure: Calling MINIMIZE.' call minimize(etot,varia,iretcode,nfun) call var_to_geom(nvar,varia) endif if (me.eq.king.or..not.out1file) & write(iout,*) 'SUMSL return code is',iretcode,' eval ',nfun if (isnan(etot) .or. etot.gt.1.0d4) then write (iout,*) "Energy too large",etot, & " trying another random conformation" do itrial=1,100 itmp=1 call gen_rand_conf(itmp,*30) goto 40 30 write (iout,*) 'Failed to generate random conformation', & ', itrial=',itrial write (*,*) 'Processor:',me, & ' Failed to generate random conformation', & ' itrial=',itrial call intout #ifdef AIX call flush_(iout) #else call flush(iout) #endif enddo write (iout,'(a,i3,a)') 'Processor:',me, & ' error in generating random conformation.' write (*,'(a,i3,a)') 'Processor:',me, & ' error in generating random conformation.' call flush(iout) #ifdef MPI call MPI_Abort(mpi_comm_world,error_msg,ierrcode) #else stop #endif 40 continue else goto 44 endif enddo write (iout,'(a,i3,a)') 'Processor:',me, & ' failed to generate a low-energy random conformation.' write (*,'(a,i3,a)') 'Processor:',me, & ' failed to generate a low-energy random conformation.' call flush(iout) #ifdef MPI call MPI_Abort(mpi_comm_world,error_msg,ierrcode) #else stop #endif 44 continue else if (preminim) then if (start_from_model) then n_model_try=0 fail=.true. list_model_try=0 do while (fail .and. n_model_try.lt.nmodel_start) write (iout,*) "n_model_try",n_model_try do i_model=iran_num(1,nmodel_start) do k=1,n_model_try if (i_model.eq.list_model_try(k)) exit enddo if (k.gt.n_model_try) exit enddo n_model_try=n_model_try+1 list_model_try(n_model_try)=i_model if (me.eq.king .or. .not. out1file) & write (iout,*) 'Trying to start from model ', & pdbfiles_chomo(i_model)(:ilen(pdbfiles_chomo(i_model))) do i=1,2*nres do j=1,3 c(j,i)=chomo(j,i,i_model) enddo enddo call int_from_cart(.true.,.false.) call sc_loc_geom(.false.) do i=1,nres-1 do j=1,3 dc(j,i)=c(j,i+1)-c(j,i) dc_norm(j,i)=dc(j,i)*vbld_inv(i+1) enddo enddo do i=2,nres-1 do j=1,3 dc(j,i+nres)=c(j,i+nres)-c(j,i) dc_norm(j,i+nres)=dc(j,i+nres)*vbld_inv(i+nres) enddo enddo if (me.eq.king.or..not.out1file) then write (iout,*) "Energies before removing overlaps" call etotal(energia(0)) call enerprint(energia(0)) endif ! Remove SC overlaps if requested if (overlapsc) then write (iout,*) 'Calling OVERLAP_SC' call overlap_sc(fail) if (fail) then write (iout,*) & "Failed to remove overlap from model",i_model cycle endif endif if (me.eq.king.or..not.out1file) then write (iout,*) "Energies after removing overlaps" call etotal(energia(0)) call enerprint(energia(0)) endif #ifdef SEARCHSC ! Search for better SC rotamers if requested if (searchsc) then call sc_move(2,nres-1,10,1d10,nft_sc,etot) print *,'SC_move',nft_sc,etot if (me.eq.king.or..not.out1file) & write(iout,*) 'SC_move',nft_sc,etot endif call etotal(energia(0)) #endif enddo call MPI_Gather(i_model,1,MPI_INTEGER,i_start_models(0), & 1,MPI_INTEGER,king,CG_COMM,IERROR) if (n_model_try.gt.nmodel_start .and. & (me.eq.king .or. out1file)) then write (iout,*) & "All models have irreparable overlaps. Trying randoms starts." iranconf=1 i_model=nmodel_start+1 goto 122 endif else ! Remove SC overlaps if requested if (overlapsc) then write (iout,*) 'Calling OVERLAP_SC' call overlap_sc(fail) if (fail) then write (iout,*) & "Failed to remove overlap" endif endif if (me.eq.king.or..not.out1file) then write (iout,*) "Energies after removing overlaps" call etotal(energia(0)) call enerprint(energia(0)) endif endif C 8/22/17 AL Minimize initial structure if (dccart) then if (me.eq.king.or..not.out1file) write(iout,*) & 'Minimizing initial PDB structure: Calling MINIM_DC' call minim_dc(etot,iretcode,nfun) else call geom_to_var(nvar,varia) if(me.eq.king.or..not.out1file) write (iout,*) & 'Minimizing initial PDB structure: Calling MINIMIZE.' call minimize(etot,varia,iretcode,nfun) call var_to_geom(nvar,varia) #ifdef LBFGS if (me.eq.king.or..not.out1file) & write(iout,*) 'LBFGS return code is ',status,' eval ',nfun if(me.eq.king.or..not.out1file) & write(iout,*) 'LBFGS return code is ',status,' eval ',nfun #else if (me.eq.king.or..not.out1file) & write(iout,*) 'SUMSL return code is',iretcode,' eval ',nfun if(me.eq.king.or..not.out1file) & write(iout,*) 'SUMSL return code is',iretcode,' eval ',nfun #endif endif endif if (nmodel_start.gt.0 .and. me.eq.king) then write (iout,'(a)') "Task Starting model" do i=0,nodes-1 if (i_start_models(i).gt.nmodel_start) then write (iout,'(i4,2x,a)') i,"RANDOM STRUCTURE" else write(iout,'(i4,2x,a)')i,pdbfiles_chomo(i_start_models(i)) & (:ilen(pdbfiles_chomo(i_start_models(i)))) endif enddo endif endif ! .not. rest call chainbuild_cart call kinetic(EK) if (tbf) then call verlet_bath endif kinetic_T=2.0d0/(dimen3*Rb)*EK write (iout,*) "Initial kinetic energy",EK," kinetic T",kinetic_T if(me.eq.king.or..not.out1file)then call cartprint call intout endif #ifdef MPI tt0=MPI_Wtime() #else tt0=tcpu() #endif call zerograd call etotal(potEcomp) if (large) call enerprint(potEcomp) #ifdef TIMING_ENE #ifdef MPI t_etotal=t_etotal+MPI_Wtime()-tt0 #else t_etotal=t_etotal+tcpu()-tt0 #endif #endif potE=potEcomp(0) c write (iout,*) "PotE-homology",potE-potEcomp(27) call cartgrad call lagrangian call max_accel if (amax*d_time .gt. dvmax) then d_time=d_time*dvmax/amax if(me.eq.king.or..not.out1file) write (iout,*) & "Time step reduced to",d_time, & " because of too large initial acceleration." endif if(me.eq.king.or..not.out1file)then write(iout,*) "Potential energy and its components" call enerprint(potEcomp) c write(iout,*) (potEcomp(i),i=0,n_ene) endif potE=potEcomp(0)-potEcomp(27) c write (iout,*) "PotE-homology",potE totE=EK+potE itime=0 itime_mat=itime if (ntwe.ne.0) call statout(itime) if(me.eq.king.or..not.out1file) & write (iout,'(/a/3(a25,1pe14.5/))') "Initial:", & " Kinetic energy",EK," Potential energy",potE, & " Total energy",totE," Maximum acceleration ", & amax if (large) then write (iout,*) "Initial coordinates" do i=1,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(c(j,i),j=1,3), & (c(j,i+nres),j=1,3) enddo write (iout,*) "Initial dC" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(dc(j,i),j=1,3), & (dc(j,i+nres),j=1,3) enddo write (iout,*) "Initial velocities" write (iout,"(13x,' backbone ',23x,' side chain')") do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_t(j,i),j=1,3), & (d_t(j,i+nres),j=1,3) enddo write (iout,*) "Initial accelerations" do i=0,nres c write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), write (iout,'(i3,3f15.10,3x,3f15.10)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo endif do i=0,2*nres do j=1,3 dc_old(j,i)=dc(j,i) d_t_old(j,i)=d_t(j,i) d_a_old(j,i)=d_a(j,i) enddo c write (iout,*) "dc_old",i,(dc_old(j,i),j=1,3) enddo if (RESPA) then #ifdef MPI tt0 =MPI_Wtime() #else tt0 = tcpu() #endif call zerograd call etotal_short(energia_short) if (large) call enerprint(potEcomp) #ifdef TIMING_ENE #ifdef MPI t_eshort=t_eshort+MPI_Wtime()-tt0 #else t_eshort=t_eshort+tcpu()-tt0 #endif #endif call cartgrad call lagrangian if(.not.out1file .and. large) then write (iout,*) "energia_long",energia_long(0), & " energia_short",energia_short(0), & " total",energia_long(0)+energia_short(0) write (iout,*) "Initial fast-force accelerations" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo endif C 7/2/2009 Copy accelerations due to short-lange forces to an auxiliary array do i=0,2*nres do j=1,3 d_a_short(j,i)=d_a(j,i) enddo enddo #ifdef MPI tt0=MPI_Wtime() #else tt0=tcpu() #endif call zerograd call etotal_long(energia_long) if (large) call enerprint(potEcomp) #ifdef TIMING_ENE #ifdef MPI t_elong=t_elong+MPI_Wtime()-tt0 #else t_elong=t_elong+tcpu()-tt0 #endif #endif call cartgrad call lagrangian if(.not.out1file .and. large) then write (iout,*) "energia_long",energia_long(0) write (iout,*) "Initial slow-force accelerations" do i=0,nres write (iout,'(i3,3f10.5,3x,3f10.5)') i,(d_a(j,i),j=1,3), & (d_a(j,i+nres),j=1,3) enddo endif #ifdef MPI t_enegrad=t_enegrad+MPI_Wtime()-tt0 #else t_enegrad=t_enegrad+tcpu()-tt0 #endif endif return end c----------------------------------------------------------- subroutine random_vel implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.TIME1' double precision xv,sigv,lowb,highb,vec_afm(3),Ek1,Ek2,Ek3,aux integer i,ii,j,k,l,ind double precision anorm_distr logical lprn /.false./ #ifdef FIVEDIAG integer ichain,n,innt,inct,ibeg,ierr double precision work(8*maxres6) integer iwork(maxres6) double precision Ghalf(mmaxres2_chain),Geigen(maxres2_chain), & Gvec(maxres2_chain,maxres2_chain) common /przechowalnia/Ghalf,Geigen,Gvec #ifdef DEBUG double precision inertia(maxres2_chain,maxres2_chain) #endif c Generate random velocities from Gaussian distribution of mean 0 and std of KT/m c First generate velocities in the eigenspace of the G matrix c write (iout,*) "Calling random_vel dimen dimen3",dimen,dimen3 c call flush(iout) #ifdef DEBUG write (iout,*) "Random_vel, fivediag" #endif d_t=0.0d0 Ek2=0.0d0 EK=0.0d0 Ek3=0.0d0 do ichain=1,nchain ind=0 ghalf=0.0d0 n=dimen_chain(ichain) innt=iposd_chain(ichain) inct=innt+n-1 #ifdef DEBUG write (iout,*) "Chain",ichain," n",n," start",innt do i=innt,inct if (i.lt.inct-1) then write (iout,'(2i3,3f10.5)') i,i-innt+1,DMorig(i),DU1orig(i), & DU2orig(i) else if (i.eq.inct-1) then write (iout,'(2i3,3f10.5)') i,i-innt+1,DMorig(i),DU1orig(i) else write (iout,'(2i3,3f10.5)') i,i-innt+1,DMorig(i) endif enddo #endif ghalf(ind+1)=dmorig(innt) ghalf(ind+2)=du1orig(innt) ghalf(ind+3)=dmorig(innt+1) ind=ind+3 do i=3,n ind=ind+i-3 c write (iout,*) "i",i," ind",ind," indu2",innt+i-2, c & " indu1",innt+i-1," indm",innt+i ghalf(ind+1)=du2orig(innt-1+i-2) ghalf(ind+2)=du1orig(innt-1+i-1) ghalf(ind+3)=dmorig(innt-1+i) c write (iout,'(3(a,i2,1x))') "DU2",innt-1+i-2, c & "DU1",innt-1+i-1,"DM ",innt-1+i ind=ind+3 enddo #ifdef DEBUG ind=0 do i=1,n do j=1,i ind=ind+1 inertia(i,j)=ghalf(ind) inertia(j,i)=ghalf(ind) enddo enddo #endif #ifdef DEBUG write (iout,*) "Chain ",ichain," ind",ind," dim",n*(n+1)/2 write (iout,*) "Five-diagonal inertia matrix, lower triangle" call matoutr(n,ghalf) #endif call gldiag(maxres2_chain,n,n,Ghalf,work,Geigen,Gvec,ierr,iwork) if (large) then write (iout,'(//a,i3)') & "Eigenvectors and eigenvalues of the G matrix chain",ichain call eigout(n,n,maxres2_chain,maxres2_chain,Gvec,Geigen) endif #ifdef DIAGCHECK c check diagonalization do i=1,n do j=1,n aux=0.0d0 do k=1,n do l=1,n aux=aux+gvec(k,i)*gvec(l,j)*inertia(k,l) enddo enddo if (i.eq.j) then write (iout,*) i,j,aux,geigen(i) else write (iout,*) i,j,aux endif enddo enddo #endif xv=0.0d0 ii=0 do i=1,n do k=1,3 ii=ii+1 sigv=dsqrt((Rb*t_bath)/geigen(i)) lowb=-5*sigv highb=5*sigv d_t_work_new(ii)=anorm_distr(xv,sigv,lowb,highb) EK=EK+0.5d0*geigen(i)*d_t_work_new(ii)**2 c write (iout,*) "i",i," ii",ii," geigen",geigen(i), c & " d_t_work_new",d_t_work_new(ii) enddo enddo do k=1,3 do i=1,n ind=(i-1)*3+k d_t_work(ind)=0.0d0 do j=1,n d_t_work(ind)=d_t_work(ind) & +Gvec(i,j)*d_t_work_new((j-1)*3+k) enddo c write (iout,*) "i",i," ind",ind," d_t_work",d_t_work(ind) c call flush(iout) enddo enddo #ifdef DEBUG aux=0.0d0 do k=1,3 do i=1,n do j=1,n aux=aux+inertia(i,j)*d_t_work(3*(i-1)+k)*d_t_work(3*(j-1)+k) enddo enddo enddo Ek3=Ek3+aux/2 #endif c Transfer to the d_t vector innt=chain_border(1,ichain) inct=chain_border(2,ichain) ind=0 c write (iout,*) "ichain",ichain," innt",innt," inct",inct do i=innt,inct do j=1,3 ind=ind+1 d_t(j,i)=d_t_work(ind) enddo if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then do j=1,3 ind=ind+1 d_t(j,i+nres)=d_t_work(ind) enddo endif enddo enddo if (large) then write (iout,*) write (iout,*) "Random velocities in the Calpha,SC space" do i=1,nres write (iout,'(a3,1h(,i5,1h),3f10.5,3x,3f10.5)') & restyp(itype(i)),i,(d_t(j,i),j=1,3),(d_t(j,i+nres),j=1,3) enddo endif call kinetic_CASC(Ek1) ! ! Transform the velocities to virtual-bond space ! #define WLOS #ifdef WLOS if (nnt.eq.1) then d_t(:,0)=d_t(:,1) endif do i=1,nres if (itype(i).eq.10 .or. itype(i).eq.ntyp1) then do j=1,3 d_t(j,i)=d_t(j,i+1)-d_t(j,i) enddo else do j=1,3 d_t(j,i+nres)=d_t(j,i+nres)-d_t(j,i) d_t(j,i)=d_t(j,i+1)-d_t(j,i) enddo end if enddo d_t(:,nres)=0.0d0 d_t(:,nct)=0.0d0 d_t(:,2*nres)=0.0d0 if (nnt.gt.1) then d_t(:,0)=d_t(:,1) d_t(:,1)=0.0d0 endif c d_a(:,0)=d_a(:,1) c d_a(:,1)=0.0d0 c write (iout,*) "Shifting accelerations" do ichain=2,nchain c write (iout,*) "ichain",chain_border1(1,ichain)-1, c & chain_border1(1,ichain) d_t(:,chain_border1(1,ichain)-1)=d_t(:,chain_border1(1,ichain)) d_t(:,chain_border1(1,ichain))=0.0d0 enddo c write (iout,*) "Adding accelerations" do ichain=2,nchain c write (iout,*) "chain",ichain,chain_border1(1,ichain)-1, c & chain_border(2,ichain-1) d_t(:,chain_border1(1,ichain)-1)= & d_t(:,chain_border1(1,ichain)-1)+d_t(:,chain_border(2,ichain-1)) d_t(:,chain_border(2,ichain-1))=0.0d0 enddo do ichain=2,nchain write (iout,*) "chain",ichain,chain_border1(1,ichain)-1, & chain_border(2,ichain-1) d_t(:,chain_border1(1,ichain)-1)= & d_t(:,chain_border1(1,ichain)-1)+d_t(:,chain_border(2,ichain-1)) d_t(:,chain_border(2,ichain-1))=0.0d0 enddo #else ibeg=0 c do j=1,3 c d_t(j,0)=d_t(j,nnt) c enddo do ichain=1,nchain innt=chain_border(1,ichain) inct=chain_border(2,ichain) c write (iout,*) "ichain",ichain," innt",innt," inct",inct c write (iout,*) "ibeg",ibeg do j=1,3 d_t(j,ibeg)=d_t(j,innt) enddo ibeg=inct+1 do i=innt,inct if (iabs(itype(i).eq.10)) then c write (iout,*) "i",i,(d_t(j,i),j=1,3),(d_t(j,i+1),j=1,3) do j=1,3 d_t(j,i)=d_t(j,i+1)-d_t(j,i) enddo else do j=1,3 d_t(j,i+nres)=d_t(j,i+nres)-d_t(j,i) d_t(j,i)=d_t(j,i+1)-d_t(j,i) enddo end if enddo enddo #endif if (large) then write (iout,*) write (iout,*) & "Random velocities in the virtual-bond-vector space" write (iout,'(3hORG,1h(,i5,1h),3f10.5)') 0,(d_t(j,0),j=1,3) do i=1,nres write (iout,'(a3,1h(,i5,1h),3f10.5,3x,3f10.5)') & restyp(itype(i)),i,(d_t(j,i),j=1,3),(d_t(j,i+nres),j=1,3) enddo write (iout,*) write (iout,*) "Kinetic energy from inertia matrix eigenvalues", & Ek write (iout,*) & "Kinetic temperatures from inertia matrix eigenvalues", & 2*Ek/(3*dimen*Rb) #ifdef DEBUG write (iout,*) "Kinetic energy from inertia matrix",Ek3 write (iout,*) "Kinetic temperatures from inertia", & 2*Ek3/(3*dimen*Rb) #endif write (iout,*) "Kinetic energy from velocities in CA-SC space", & Ek1 write (iout,*) & "Kinetic temperatures from velovities in CA-SC space", & 2*Ek1/(3*dimen*Rb) call kinetic(Ek1) write (iout,*) & "Kinetic energy from virtual-bond-vector velocities",Ek1 write (iout,*) & "Kinetic temperature from virtual-bond-vector velocities ", & 2*Ek1/(dimen3*Rb) endif #else xv=0.0d0 ii=0 do i=1,dimen do k=1,3 ii=ii+1 sigv=dsqrt((Rb*t_bath)/geigen(i)) lowb=-5*sigv highb=5*sigv d_t_work_new(ii)=anorm_distr(xv,sigv,lowb,highb) c write (iout,*) "i",i," ii",ii," geigen",geigen(i), c & " d_t_work_new",d_t_work_new(ii) enddo enddo C if (SELFGUIDE.gt.0) then C distance=0.0 C do j=1,3 C vec_afm(j)=c(j,afmend)-c(j,afmbeg) C distance=distance+vec_afm(j)**2 C enddo C distance=dsqrt(distance) C do j=1,3 C d_t_work_new(j+(afmbeg-1)*3)=-velAFMconst*vec_afm(j)/distance C d_t_work_new(j+(afmend-1)*3)=velAFMconst*vec_afm(j)/distance C write(iout,*) "myvel",d_t_work_new(j+(afmbeg-1)*3), C & d_t_work_new(j+(afmend-1)*3) C enddo C endif c diagnostics c Ek1=0.0d0 c ii=0 c do i=1,dimen c do k=1,3 c ii=ii+1 c Ek1=Ek1+0.5d0*geigen(i)*d_t_work_new(ii)**2 c enddo c enddo c write (iout,*) "Ek from eigenvectors",Ek1 c end diagnostics c Transform velocities to UNRES coordinate space do k=0,2 do i=1,dimen ind=(i-1)*3+k+1 d_t_work(ind)=0.0d0 do j=1,dimen d_t_work(ind)=d_t_work(ind) & +Gvec(i,j)*d_t_work_new((j-1)*3+k+1) enddo c write (iout,*) "i",i," ind",ind," d_t_work",d_t_work(ind) c call flush(iout) enddo enddo c Transfer to the d_t vector do j=1,3 d_t(j,0)=d_t_work(j) enddo ind=3 do i=nnt,nct-1 do j=1,3 ind=ind+1 d_t(j,i)=d_t_work(ind) enddo enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then do j=1,3 ind=ind+1 d_t(j,i+nres)=d_t_work(ind) enddo endif enddo c call kinetic(EK) c write (iout,*) "Kinetic energy",Ek,EK1," kinetic temperature", c & 2.0d0/(dimen3*Rb)*EK,2.0d0/(dimen3*Rb)*EK1 c call flush(iout) #endif return end #ifndef LANG0 c----------------------------------------------------------- subroutine sd_verlet_p_setup c Sets up the parameters of stochastic Verlet algorithm implicit none include 'DIMENSIONS' #ifdef MPI include 'mpif.h' #endif include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.TIME1' double precision emgdt(MAXRES6), & pterm,vterm,rho,rhoc,vsig, & pfric_vec(MAXRES6),vfric_vec(MAXRES6), & afric_vec(MAXRES6),prand_vec(MAXRES6), & vrand_vec1(MAXRES6),vrand_vec2(MAXRES6) logical lprn /.false./ double precision zero /1.0d-8/, gdt_radius /0.05d0/ double precision ktm #ifdef MPI tt0 = MPI_Wtime() #else tt0 = tcpu() #endif c c AL 8/17/04 Code adapted from tinker c c Get the frictional and random terms for stochastic dynamics in the c eigenspace of mass-scaled UNRES friction matrix c do i = 1, dimen gdt = fricgam(i) * d_time c c Stochastic dynamics reduces to simple MD for zero friction c if (gdt .le. zero) then pfric_vec(i) = 1.0d0 vfric_vec(i) = d_time afric_vec(i) = 0.5d0 * d_time * d_time prand_vec(i) = 0.0d0 vrand_vec1(i) = 0.0d0 vrand_vec2(i) = 0.0d0 c c Analytical expressions when friction coefficient is large c else if (gdt .ge. gdt_radius) then egdt = dexp(-gdt) pfric_vec(i) = egdt vfric_vec(i) = (1.0d0-egdt) / fricgam(i) afric_vec(i) = (d_time-vfric_vec(i)) / fricgam(i) pterm = 2.0d0*gdt - 3.0d0 + (4.0d0-egdt)*egdt vterm = 1.0d0 - egdt**2 rho = (1.0d0-egdt)**2 / sqrt(pterm*vterm) c c Use series expansions when friction coefficient is small c else gdt2 = gdt * gdt gdt3 = gdt * gdt2 gdt4 = gdt2 * gdt2 gdt5 = gdt2 * gdt3 gdt6 = gdt3 * gdt3 gdt7 = gdt3 * gdt4 gdt8 = gdt4 * gdt4 gdt9 = gdt4 * gdt5 afric_vec(i) = (gdt2/2.0d0 - gdt3/6.0d0 + gdt4/24.0d0 & - gdt5/120.0d0 + gdt6/720.0d0 & - gdt7/5040.0d0 + gdt8/40320.0d0 & - gdt9/362880.0d0) / fricgam(i)**2 vfric_vec(i) = d_time - fricgam(i)*afric_vec(i) pfric_vec(i) = 1.0d0 - fricgam(i)*vfric_vec(i) pterm = 2.0d0*gdt3/3.0d0 - gdt4/2.0d0 & + 7.0d0*gdt5/30.0d0 - gdt6/12.0d0 & + 31.0d0*gdt7/1260.0d0 - gdt8/160.0d0 & + 127.0d0*gdt9/90720.0d0 vterm = 2.0d0*gdt - 2.0d0*gdt2 + 4.0d0*gdt3/3.0d0 & - 2.0d0*gdt4/3.0d0 + 4.0d0*gdt5/15.0d0 & - 4.0d0*gdt6/45.0d0 + 8.0d0*gdt7/315.0d0 & - 2.0d0*gdt8/315.0d0 + 4.0d0*gdt9/2835.0d0 rho = sqrt(3.0d0) * (0.5d0 - 3.0d0*gdt/16.0d0 & - 17.0d0*gdt2/1280.0d0 & + 17.0d0*gdt3/6144.0d0 & + 40967.0d0*gdt4/34406400.0d0 & - 57203.0d0*gdt5/275251200.0d0 & - 1429487.0d0*gdt6/13212057600.0d0) end if c c Compute the scaling factors of random terms for the nonzero friction case c ktm = 0.5d0*d_time/fricgam(i) psig = dsqrt(ktm*pterm) / fricgam(i) vsig = dsqrt(ktm*vterm) rhoc = dsqrt(1.0d0 - rho*rho) prand_vec(i) = psig vrand_vec1(i) = vsig * rho vrand_vec2(i) = vsig * rhoc end if end do if (lprn) then write (iout,*) & "pfric_vec, vfric_vec, afric_vec, prand_vec, vrand_vec1,", & " vrand_vec2" do i=1,dimen write (iout,'(i5,6e15.5)') i,pfric_vec(i),vfric_vec(i), & afric_vec(i),prand_vec(i),vrand_vec1(i),vrand_vec2(i) enddo endif c c Transform from the eigenspace of mass-scaled friction matrix to UNRES variables c call eigtransf(dimen,maxres2,mt3,mt2,pfric_vec,pfric_mat) call eigtransf(dimen,maxres2,mt3,mt2,vfric_vec,vfric_mat) call eigtransf(dimen,maxres2,mt3,mt2,afric_vec,afric_mat) call eigtransf(dimen,maxres2,mt3,mt1,prand_vec,prand_mat) call eigtransf(dimen,maxres2,mt3,mt1,vrand_vec1,vrand_mat1) call eigtransf(dimen,maxres2,mt3,mt1,vrand_vec2,vrand_mat2) #ifdef MPI t_sdsetup=t_sdsetup+MPI_Wtime() #else t_sdsetup=t_sdsetup+tcpu()-tt0 #endif return end c------------------------------------------------------------- subroutine eigtransf1(n,ndim,ab,d,c) implicit none integer n,ndim double precision ab(ndim,ndim,n),c(ndim,n),d(ndim) integer i,j,k do i=1,n do j=1,n c(i,j)=0.0d0 do k=1,n c(i,j)=c(i,j)+ab(k,j,i)*d(k) enddo enddo enddo return end c------------------------------------------------------------- subroutine eigtransf(n,ndim,a,b,d,c) implicit none integer n,ndim double precision a(ndim,n),b(ndim,n),c(ndim,n),d(ndim) integer i,j,k do i=1,n do j=1,n c(i,j)=0.0d0 do k=1,n c(i,j)=c(i,j)+a(i,k)*b(k,j)*d(k) enddo enddo enddo return end c------------------------------------------------------------- subroutine sd_verlet1 c Applying stochastic velocity Verlet algorithm - step 1 to velocities implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' double precision stochforcvec(MAXRES6) common /stochcalc/ stochforcvec logical lprn /.false./ integer i,j,ind,inres c write (iout,*) "dc_old" c do i=0,nres c write (iout,'(i5,3f10.5,5x,3f10.5)') c & i,(dc_old(j,i),j=1,3),(dc_old(j,i+nres),j=1,3) c enddo do j=1,3 dc_work(j)=dc_old(j,0) d_t_work(j)=d_t_old(j,0) d_a_work(j)=d_a_old(j,0) enddo ind=3 do i=nnt,nct-1 do j=1,3 dc_work(ind+j)=dc_old(j,i) d_t_work(ind+j)=d_t_old(j,i) d_a_work(ind+j)=d_a_old(j,i) enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then do j=1,3 dc_work(ind+j)=dc_old(j,i+nres) d_t_work(ind+j)=d_t_old(j,i+nres) d_a_work(ind+j)=d_a_old(j,i+nres) enddo ind=ind+3 endif enddo #ifndef LANG0 if (lprn) then write (iout,*) & "pfric_mat, vfric_mat, afric_mat, prand_mat, vrand_mat1,", & " vrand_mat2" do i=1,dimen do j=1,dimen write (iout,'(2i5,6e15.5)') i,j,pfric_mat(i,j), & vfric_mat(i,j),afric_mat(i,j), & prand_mat(i,j),vrand_mat1(i,j),vrand_mat2(i,j) enddo enddo endif do i=1,dimen ddt1=0.0d0 ddt2=0.0d0 do j=1,dimen dc_work(i)=dc_work(i)+vfric_mat(i,j)*d_t_work(j) & +afric_mat(i,j)*d_a_work(j)+prand_mat(i,j)*stochforcvec(j) ddt1=ddt1+pfric_mat(i,j)*d_t_work(j) ddt2=ddt2+vfric_mat(i,j)*d_a_work(j) enddo d_t_work_new(i)=ddt1+0.5d0*ddt2 d_t_work(i)=ddt1+ddt2 enddo #endif do j=1,3 dc(j,0)=dc_work(j) d_t(j,0)=d_t_work(j) enddo ind=3 do i=nnt,nct-1 do j=1,3 dc(j,i)=dc_work(ind+j) d_t(j,i)=d_t_work(ind+j) enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10) then inres=i+nres do j=1,3 dc(j,inres)=dc_work(ind+j) d_t(j,inres)=d_t_work(ind+j) enddo ind=ind+3 endif enddo return end c-------------------------------------------------------------------------- subroutine sd_verlet2 c Calculating the adjusted velocities for accelerations implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif #ifndef LANG0 include 'COMMON.LANGEVIN' #else #ifdef FIVEDIAG include 'COMMON.LANGEVIN.lang0.5diag' #else include 'COMMON.LANGEVIN.lang0' #endif #endif include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' double precision stochforcvec(MAXRES6),stochforcvecV(MAXRES6) common /stochcalc/ stochforcvec integer i,j,ind,inres c c Compute the stochastic forces which contribute to velocity change c call stochastic_force(stochforcvecV) #ifndef LANG0 do i=1,dimen ddt1=0.0d0 ddt2=0.0d0 do j=1,dimen ddt1=ddt1+vfric_mat(i,j)*d_a_work(j) ddt2=ddt2+vrand_mat1(i,j)*stochforcvec(j)+ & vrand_mat2(i,j)*stochforcvecV(j) enddo d_t_work(i)=d_t_work_new(i)+0.5d0*ddt1+ddt2 enddo #endif do j=1,3 d_t(j,0)=d_t_work(j) enddo ind=3 do i=nnt,nct-1 do j=1,3 d_t(j,i)=d_t_work(ind+j) enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 d_t(j,inres)=d_t_work(ind+j) enddo ind=ind+3 endif enddo return end c----------------------------------------------------------- subroutine sd_verlet_ciccotti_setup c Sets up the parameters of stochastic velocity Verlet algorithmi; Ciccotti's c version implicit none include 'DIMENSIONS' #ifdef MPI include 'mpif.h' #endif include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.LANGEVIN' include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' include 'COMMON.TIME1' double precision emgdt(MAXRES6), & pterm,vterm,rho,rhoc,vsig, & pfric_vec(MAXRES6),vfric_vec(MAXRES6), & afric_vec(MAXRES6),prand_vec(MAXRES6), & vrand_vec1(MAXRES6),vrand_vec2(MAXRES6) logical lprn /.false./ double precision zero /1.0d-8/, gdt_radius /0.05d0/ double precision ktm integer i #ifdef MPI tt0 = MPI_Wtime() #else tt0 = tcpu() #endif c c AL 8/17/04 Code adapted from tinker c c Get the frictional and random terms for stochastic dynamics in the c eigenspace of mass-scaled UNRES friction matrix c do i = 1, dimen write (iout,*) "i",i," fricgam",fricgam(i) gdt = fricgam(i) * d_time c c Stochastic dynamics reduces to simple MD for zero friction c if (gdt .le. zero) then pfric_vec(i) = 1.0d0 vfric_vec(i) = d_time afric_vec(i) = 0.5d0*d_time*d_time prand_vec(i) = afric_vec(i) vrand_vec2(i) = vfric_vec(i) c c Analytical expressions when friction coefficient is large c else egdt = dexp(-gdt) pfric_vec(i) = egdt vfric_vec(i) = dexp(-0.5d0*gdt)*d_time afric_vec(i) = 0.5d0*dexp(-0.25d0*gdt)*d_time*d_time prand_vec(i) = afric_vec(i) vrand_vec2(i) = vfric_vec(i) c c Compute the scaling factors of random terms for the nonzero friction case c c ktm = 0.5d0*d_time/fricgam(i) c psig = dsqrt(ktm*pterm) / fricgam(i) c vsig = dsqrt(ktm*vterm) c prand_vec(i) = psig*afric_vec(i) c vrand_vec2(i) = vsig*vfric_vec(i) end if end do if (lprn) then write (iout,*) & "pfric_vec, vfric_vec, afric_vec, prand_vec, vrand_vec1,", & " vrand_vec2" do i=1,dimen write (iout,'(i5,6e15.5)') i,pfric_vec(i),vfric_vec(i), & afric_vec(i),prand_vec(i),vrand_vec1(i),vrand_vec2(i) enddo endif c c Transform from the eigenspace of mass-scaled friction matrix to UNRES variables c call eigtransf(dimen,maxres2,mt3,mt2,pfric_vec,pfric_mat) call eigtransf(dimen,maxres2,mt3,mt2,vfric_vec,vfric_mat) call eigtransf(dimen,maxres2,mt3,mt2,afric_vec,afric_mat) call eigtransf(dimen,maxres2,mt3,mt1,prand_vec,prand_mat) call eigtransf(dimen,maxres2,mt3,mt1,vrand_vec2,vrand_mat2) #ifdef MPI t_sdsetup=t_sdsetup+MPI_Wtime() #else t_sdsetup=t_sdsetup+tcpu()-tt0 #endif return end c------------------------------------------------------------- subroutine sd_verlet1_ciccotti c Applying stochastic velocity Verlet algorithm - step 1 to velocities implicit none include 'DIMENSIONS' #ifdef MPI include 'mpif.h' #endif include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.LANGEVIN' include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' double precision stochforcvec(MAXRES6) common /stochcalc/ stochforcvec logical lprn /.false./ integer i,j c write (iout,*) "dc_old" c do i=0,nres c write (iout,'(i5,3f10.5,5x,3f10.5)') c & i,(dc_old(j,i),j=1,3),(dc_old(j,i+nres),j=1,3) c enddo do j=1,3 dc_work(j)=dc_old(j,0) d_t_work(j)=d_t_old(j,0) d_a_work(j)=d_a_old(j,0) enddo ind=3 do i=nnt,nct-1 do j=1,3 dc_work(ind+j)=dc_old(j,i) d_t_work(ind+j)=d_t_old(j,i) d_a_work(ind+j)=d_a_old(j,i) enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then do j=1,3 dc_work(ind+j)=dc_old(j,i+nres) d_t_work(ind+j)=d_t_old(j,i+nres) d_a_work(ind+j)=d_a_old(j,i+nres) enddo ind=ind+3 endif enddo #ifndef LANG0 if (lprn) then write (iout,*) & "pfric_mat, vfric_mat, afric_mat, prand_mat, vrand_mat1,", & " vrand_mat2" do i=1,dimen do j=1,dimen write (iout,'(2i5,6e15.5)') i,j,pfric_mat(i,j), & vfric_mat(i,j),afric_mat(i,j), & prand_mat(i,j),vrand_mat1(i,j),vrand_mat2(i,j) enddo enddo endif do i=1,dimen ddt1=0.0d0 ddt2=0.0d0 do j=1,dimen dc_work(i)=dc_work(i)+vfric_mat(i,j)*d_t_work(j) & +afric_mat(i,j)*d_a_work(j)+prand_mat(i,j)*stochforcvec(j) ddt1=ddt1+pfric_mat(i,j)*d_t_work(j) ddt2=ddt2+vfric_mat(i,j)*d_a_work(j) enddo d_t_work_new(i)=ddt1+0.5d0*ddt2 d_t_work(i)=ddt1+ddt2 enddo #endif do j=1,3 dc(j,0)=dc_work(j) d_t(j,0)=d_t_work(j) enddo ind=3 do i=nnt,nct-1 do j=1,3 dc(j,i)=dc_work(ind+j) d_t(j,i)=d_t_work(ind+j) enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 dc(j,inres)=dc_work(ind+j) d_t(j,inres)=d_t_work(ind+j) enddo ind=ind+3 endif enddo return end c-------------------------------------------------------------------------- subroutine sd_verlet2_ciccotti c Calculating the adjusted velocities for accelerations implicit none include 'DIMENSIONS' include 'COMMON.CONTROL' include 'COMMON.VAR' include 'COMMON.MD' #ifdef FIVEDIAG include 'COMMON.LAGRANGE.5diag' #else include 'COMMON.LAGRANGE' #endif include 'COMMON.LANGEVIN' include 'COMMON.CHAIN' include 'COMMON.DERIV' include 'COMMON.GEO' include 'COMMON.LOCAL' include 'COMMON.INTERACT' include 'COMMON.IOUNITS' include 'COMMON.NAMES' double precision stochforcvec(MAXRES6),stochforcvecV(MAXRES6) common /stochcalc/ stochforcvec integer i,j c c Compute the stochastic forces which contribute to velocity change c call stochastic_force(stochforcvecV) #ifndef LANG0 do i=1,dimen ddt1=0.0d0 ddt2=0.0d0 do j=1,dimen ddt1=ddt1+vfric_mat(i,j)*d_a_work(j) c ddt2=ddt2+vrand_mat2(i,j)*stochforcvecV(j) ddt2=ddt2+vrand_mat2(i,j)*stochforcvec(j) enddo d_t_work(i)=d_t_work_new(i)+0.5d0*ddt1+ddt2 enddo #endif do j=1,3 d_t(j,0)=d_t_work(j) enddo ind=3 do i=nnt,nct-1 do j=1,3 d_t(j,i)=d_t_work(ind+j) enddo ind=ind+3 enddo do i=nnt,nct if (itype(i).ne.10 .and. itype(i).ne.ntyp1) then inres=i+nres do j=1,3 d_t(j,inres)=d_t_work(ind+j) enddo ind=ind+3 endif enddo return end #endif