9914d9babd0afa1532ac2190f403c99e23ed82d2
[unres4.git] / source / unres / MREMD.F90
1       module MREMDyn
2 !-----------------------------------------------------------------------------
3       use io_units
4       use names
5       use MPI_data
6       use md_data
7       use remd_data
8       use geometry_data
9       use energy_data
10 !      use control_data, only:maxprocs
11       use MDyn
12
13       implicit none
14 !-----------------------------------------------------------------------------
15 ! commom.remd
16 !      common /remdrestart/
17       integer(kind=2),dimension(:),allocatable :: i2set !(0:maxprocs)
18       integer(kind=2),dimension(:),allocatable :: ifirst !(maxprocs)
19       integer(kind=2),dimension(:,:),allocatable :: nupa,&
20        ndowna !(0:maxprocs/4,0:maxprocs)
21       real(kind=4),dimension(:,:),allocatable :: t_restart1 !(5,maxprocs)
22       integer,dimension(:),allocatable :: iset_restart1 !(maxprocs)
23 !      common /traj1cache/
24       real(kind=4),dimension(:),allocatable :: totT_cache,EK_cache,&
25        potE_cache,t_bath_cache,Uconst_cache !(max_cache_traj)
26       real(kind=4),dimension(:,:),allocatable :: qfrag_cache !(50,max_cache_traj)
27       real(kind=4),dimension(:,:),allocatable :: qpair_cache !(100,max_cache_traj)
28       real(kind=4),dimension(:,:),allocatable :: ugamma_cache,&
29        utheta_cache,uscdiff_cache !(maxfrag_back,max_cache_traj)
30       real(kind=4),dimension(:,:,:),allocatable :: c_cache !(3,maxres2+2,max_cache_traj)
31       integer :: ntwx_cache,ii_write    !,max_cache_traj_use
32       integer,dimension(:),allocatable :: iset_cache !(max_cache_traj)
33 !-----------------------------------------------------------------------------
34 !       common /przechowalnia/
35       real(kind=4),dimension(:,:),allocatable :: d_restart1 !(3,2*nres*maxprocs)
36       real(kind=4),dimension(:,:),allocatable :: d_restart2 !(3,2*nres*maxprocs)
37       real(kind=4),dimension(:,:),allocatable :: p_c !(3,(nres2+2)*maxprocs)
38 !-----------------------------------------------------------------------------
39 !
40 !
41 !-----------------------------------------------------------------------------
42       contains
43 !-----------------------------------------------------------------------------
44
45 !-----------------------------------------------------------------------------
46 ! MREMD.F
47 !-----------------------------------------------------------------------------
48
49       subroutine MREMD
50
51       use comm_gucio
52       use control, only:tcpu,ovrtim
53       use io_base, only:ilen
54       use control_data
55       use geometry_data
56       use random, only: iran_num,ran_number
57       use compare, only:hairpin,secondary2
58       use io, only:cartout,statout
59 !      implicit real*8 (a-h,o-z)
60 !      include 'DIMENSIONS'
61       include 'mpif.h'
62 !      include 'COMMON.CONTROL'
63 !      include 'COMMON.VAR'
64 !      include 'COMMON.MD'
65 !#ifndef LANG0
66 !      include 'COMMON.LANGEVIN'
67 !#else
68 !      include 'COMMON.LANGEVIN.lang0'
69 !#endif
70 !      include 'COMMON.CHAIN'
71 !      include 'COMMON.DERIV'
72 !      include 'COMMON.GEO'
73 !      include 'COMMON.LOCAL'
74 !      include 'COMMON.INTERACT'
75 !      include 'COMMON.IOUNITS'
76 !      include 'COMMON.NAMES'
77 !      include 'COMMON.TIME1'
78 !      include 'COMMON.REMD'
79 !      include 'COMMON.SETUP'
80 !      include 'COMMON.MUCA'
81 !      include 'COMMON.HAIRPIN'
82       integer :: ERRCODE
83       real(kind=8),dimension(3) :: L,vcm
84       real(kind=8) :: energia(0:n_ene)
85
86       real(kind=8) :: remd_t_bath(Nprocs) !(maxprocs)
87       integer :: iremd_iset(Nprocs) !(maxprocs)
88       integer(kind=2) :: i_index(Nprocs,Nprocs,Nprocs,Nprocs)
89 ! (maxprocs/4,maxprocs/20,maxprocs/200,maxprocs/200)
90       real(kind=8) :: remd_ene(0:n_ene+4,Nprocs) !(0:n_ene+4,maxprocs)
91       integer :: iremd_acc(Nprocs),iremd_tot(Nprocs) !(maxprocs)
92       integer :: iremd_acc_usa(Nprocs),iremd_tot_usa(Nprocs) !(maxprocs)
93       integer :: rstcount       !el ilen,
94 !el      external ilen
95       character(len=50) :: tytul
96 !el      common /gucio/ cm
97       integer :: itime
98 !old      integer nup(0:maxprocs),ndown(0:maxprocs)
99       integer :: rep2i(0:Nprocs),ireqi(Nprocs) !(maxprocs)
100       integer :: icache_all(Nprocs) !(maxprocs)
101       integer :: status(MPI_STATUS_SIZE),statusi(MPI_STATUS_SIZE,Nprocs)! (MPI_STATUS_SIZE,maxprocs)
102       logical :: synflag, end_of_run, file_exist = .false.!, ovrtim
103
104       real(kind=8) :: delta,time00,time01,time001,time02,time03,time04,&
105                  time05,time06,time07,time08,tt0,scalfac,ene_iex_iex,&
106                  ene_i_i,ene_iex_i,ene_i_iex,xxx,tmp,econstr_temp_i,&
107                  econstr_temp_iex
108       integer :: k,il,il1,i,j,nharp,ii,ierr,itime_master,irr,iex,&
109             i_set_temp,itmp,i_temp,i_mult,i_iset,i_mset,i_dir,i_temp1,&
110             i_mult1,i_iset1,i_mset1,ierror
111       integer,dimension(4,nres/3) :: iharp      !(4,nres/3)(4,maxres/3)
112 !deb      imin_itime_old=0
113       integer :: nres2  !el
114       WRITE(iout,*) "JUST AFTER CALL"
115 !      if (.not.allocated(remd_ene)) allocate(remd_ene(0:n_ene+4,Nprocs))
116       nres2=2*nres
117       time001=0.0d0
118
119       ntwx_cache=0
120       time00=MPI_WTIME()
121       time01=time00
122       if(me.eq.king.or..not.out1file) then
123        write  (iout,*) 'MREMD',nodes,'time before',time00-walltime
124        write (iout,*) "NREP=",nrep
125       endif
126
127       synflag=.false.
128       if (ilen(tmpdir).gt.0 .and. (me.eq.king)) then
129         call copy_to_tmp(pref_orig(:ilen(pref_orig))//"_mremd.rst")
130       endif
131       mremd_rst_name=prefix(:ilen(prefix))//"_mremd.rst"
132
133       print *,'MREMD',nodes
134       print *,'mmm',me,remd_mlist,(remd_m(i),i=1,nrep)
135       write (iout,*) "Start MREMD: me",me," t_bath",t_bath
136       print *,"NSET=",nset, "MSET=", mset
137       k=0
138       rep2i(k)=-1
139       do il=1,max0(nset,1)
140        do il1=1,max0(mset(il),1)
141         do i=1,nrep
142          iremd_acc(i)=0
143          iremd_acc_usa(i)=0
144          iremd_tot(i)=0
145          do j=1,remd_m(i)
146           i2rep(k)=i
147           i2set(k)=il
148           rep2i(i)=k
149           k=k+1
150           i_index(i,j,il,il1)=k
151          enddo
152         enddo
153        enddo
154       enddo
155
156       if(me.eq.king.or..not.out1file) then
157        write(iout,*) (i2rep(i),i=0,nodes-1)
158        write(iout,*) (i2set(i),i=0,nodes-1)
159        do il=1,nset
160         do il1=1,mset(il)
161          do i=1,nrep
162           do j=1,remd_m(i)
163            write(iout,*) i,j,il,il1,i_index(i,j,il,il1)
164           enddo
165          enddo
166         enddo
167        enddo
168       endif
169
170 !      print *,'i2rep',me,i2rep(me)
171 !      print *,'rep2i',(rep2i(i),i=0,nrep)
172
173 !old       if (i2rep(me).eq.nrep) then
174 !old        nup(0)=0
175 !old       else
176 !old        nup(0)=remd_m(i2rep(me)+1)
177 !old        k=rep2i(int(i2rep(me)))+1
178 !old        do i=1,nup(0)
179 !old         nup(i)=k
180 !old         k=k+1
181 !old        enddo
182 !old       endif
183
184 !d       print '(i4,a4,100i4)',me,' nup',(nup(i),i=0,nup(0))
185
186 !old       if (i2rep(me).eq.1) then
187 !old        ndown(0)=0
188 !old       else
189 !old        ndown(0)=remd_m(i2rep(me)-1)
190 !old        k=rep2i(i2rep(me)-2)+1
191 !old        do i=1,ndown(0)
192 !old         ndown(i)=k
193 !old         k=k+1
194 !old        enddo
195 !old       endif
196
197 !d       print '(i4,a6,100i4)',me,' ndown',(ndown(i),i=0,ndown(0))
198
199 !el  common /przechowalnia/
200        if(.not.allocated(d_restart1)) allocate(d_restart1(3,0:(nres2+1)*nodes))
201        if(.not.allocated(d_restart2)) allocate(d_restart2(3,0:(nres2+1)*nodes))
202        if(.not.allocated(p_c)) allocate(p_c(3,(nres2+2)*nodes))
203 !el-------------
204
205        write (*,*) "Processor",me," rest",rest,&
206          "restart1fie",restart1file
207        if(rest.and.restart1file) then 
208            if (me.eq.king) &
209            inquire(file=mremd_rst_name,exist=file_exist)
210 !d           write (*,*) me," Before broadcast: file_exist",file_exist
211            call MPI_Bcast(file_exist,1,MPI_LOGICAL,king,CG_COMM,&
212                 IERR)
213 !d           write (*,*) me," After broadcast: file_exist",file_exist
214            if(file_exist) then 
215              if(me.eq.king.or..not.out1file) &
216                   write  (iout,*) 'Master is reading restart1file'
217              call read1restart(i_index)
218            else
219              if(me.eq.king.or..not.out1file) &
220                   write  (iout,*) 'WARNING : no restart1file'
221            endif
222
223            if(me.eq.king.or..not.out1file) then
224               write(iout,*) "i2set",(i2set(i),i=0,nodes-1)
225               write(iout,*) "i_index"
226               do il=1,nset
227                do il1=1,mset(il)
228                 do i=1,nrep
229                  do j=1,remd_m(i)
230                   write(iout,*) i,j,il,il1,i_index(i,j,il,il1)
231                  enddo
232                 enddo
233                enddo
234               enddo
235            endif 
236        endif
237
238        if(me.eq.king) then
239         if (rest.and..not.restart1file) &
240                 inquire(file=mremd_rst_name,exist=file_exist)
241         if(.not.file_exist.and.rest.and..not.restart1file) &
242              write(iout,*) 'WARNING : no restart file',mremd_rst_name
243         IF (rest.and.file_exist.and..not.restart1file) THEN
244              write  (iout,*) 'Master is reading restart file',&
245                               mremd_rst_name
246              open(irest2,file=mremd_rst_name,status='unknown')
247              read (irest2,*) 
248              read (irest2,*) (i2rep(i),i=0,nodes-1)
249              read (irest2,*) 
250              read (irest2,*) (ifirst(i),i=1,remd_m(1))
251              do il=1,nodes
252               read (irest2,*) 
253               read (irest2,*) nupa(0,il),(nupa(i,il),i=1,nupa(0,il))
254               read (irest2,*) 
255               read (irest2,*) ndowna(0,il),&
256                           (ndowna(i,il),i=1,ndowna(0,il))
257              enddo
258              if(usampl) then
259               read (irest2,*)
260               read (irest2,*) nset
261               read (irest2,*) 
262               read (irest2,*) (mset(i),i=1,nset)
263               read (irest2,*) 
264               read (irest2,*) (i2set(i),i=0,nodes-1)
265               read (irest2,*) 
266               do il=1,nset
267                do il1=1,mset(il)
268                 do i=1,nrep
269                   read(irest2,*) (i_index(i,j,il,il1),j=1,remd_m(i))
270                 enddo
271                enddo
272               enddo
273
274               write(iout,*) "i2set",(i2set(i),i=0,nodes-1)
275               write(iout,*) "i_index"
276               do il=1,nset
277                do il1=1,mset(il)
278                 do i=1,nrep
279                  do j=1,remd_m(i)
280                   write(iout,*) i,j,il,il1,i_index(i,j,il,il1)
281                  enddo
282                 enddo
283                enddo
284               enddo
285              endif
286
287              close(irest2)
288
289              write (iout,'(a6,1000i5)') "i2rep",(i2rep(i),i=0,nodes-1)
290              write (iout,'(a6,1000i5)') "ifirst",&
291                           (ifirst(i),i=1,remd_m(1))
292              do il=1,nodes
293               write (iout,'(a6,i4,a1,100i4)') "nupa",il,":",&
294                           (nupa(i,il),i=1,nupa(0,il))
295               write (iout,'(a6,i4,a1,100i4)') "ndowna",il,":",&
296                           (ndowna(i,il),i=1,ndowna(0,il))
297              enddo
298         ELSE IF (.not.(rest.and.file_exist)) THEN
299          do il=1,remd_m(1)
300           ifirst(il)=il
301          enddo
302
303          do il=1,nodes
304           if (i2rep(il-1).eq.nrep) then
305            nupa(0,il)=0
306           else
307            nupa(0,il)=remd_m(i2rep(il-1)+1)
308            k=rep2i(int(i2rep(il-1)))+1
309            do i=1,nupa(0,il)
310             nupa(i,il)=k+1
311             k=k+1
312            enddo
313           endif
314           if (i2rep(il-1).eq.1) then
315            ndowna(0,il)=0
316           else
317            ndowna(0,il)=remd_m(i2rep(il-1)-1)
318            k=rep2i(i2rep(il-1)-2)+1
319            do i=1,ndowna(0,il)
320             ndowna(i,il)=k+1
321             k=k+1
322            enddo
323           endif
324          enddo
325         
326         write (iout,'(a6,100i4)') "ifirst",&
327                           (ifirst(i),i=1,remd_m(1))
328         do il=1,nodes
329          write (iout,'(a6,i4,a1,100i4)') "nupa",il,":",&
330                           (nupa(i,il),i=1,nupa(0,il))
331          write (iout,'(a6,i4,a1,100i4)') "ndowna",il,":",&
332                           (ndowna(i,il),i=1,ndowna(0,il))
333         enddo
334         
335         ENDIF
336        endif
337 !
338 !      t_bath=retmin+(retmax-retmin)*me/(nodes-1)
339        if(.not.(rest.and.file_exist.and.restart1file)) then
340          if (me .eq. king) then 
341             t_bath=retmin
342          else 
343             t_bath=retmin+(retmax-retmin)*exp(float(i2rep(me)-nrep))
344          endif
345 !d       print *,'ttt',me,remd_tlist,(remd_t(i),i=1,nrep)
346          if (remd_tlist) t_bath=remd_t(int(i2rep(me)))
347
348        endif
349        if(usampl) then
350           iset=i2set(me)
351           if(me.eq.king.or..not.out1file) &
352            write(iout,*) me,"iset=",iset,"t_bath=",t_bath
353        endif        
354 !    
355        do i=1,5
356        stdfp(i)=dsqrt(2*Rb*t_bath/d_time)
357        enddo
358        do j=1,5
359        do i=1,ntyp
360           stdfsc(i,j)=dsqrt(2*Rb*t_bath/d_time)
361        enddo 
362        enddo
363
364 !      print *,'irep',me,t_bath
365        if (.not.rest) then  
366         if (me.eq.king .or. .not. out1file) &
367          write (iout,'(a60,f10.5)') "REMD Temperature:",t_bath
368         call rescale_weights(t_bath)
369        endif
370
371
372 !------copy MD--------------
373 !  The driver for molecular dynamics subroutines
374 !------------------------------------------------
375       t_MDsetup=0.0d0
376       t_langsetup=0.0d0
377       t_MD=0.0d0
378       t_enegrad=0.0d0
379       t_sdsetup=0.0d0
380       if(me.eq.king.or..not.out1file) &
381        write (iout,'(20(1h=),a20,20(1h=))') "MD calculation started"
382 #ifdef MPI
383       tt0 = MPI_Wtime()
384 #else
385       tt0 = tcpu()
386 #endif
387 ! Determine the inverse of the inertia matrix.
388       call setup_MD_matrices
389 ! Initialize MD
390       call init_MD
391       if (rest) then  
392        if (me.eq.king .or. .not. out1file) &
393         write (iout,'(a60,f10.5)') "REMD restart Temperature:",t_bath
394        do j=1,5
395        stdfp(j)=dsqrt(2*Rb*t_bath/d_time)
396        do i=1,ntyp
397           stdfsc(i,j)=dsqrt(2*Rb*t_bath/d_time)
398        enddo 
399        enddo
400        call rescale_weights(t_bath)
401       endif
402
403 #ifdef MPI
404       t_MDsetup = MPI_Wtime()-tt0
405 #else
406       t_MDsetup = tcpu()-tt0
407 #endif
408       rstcount=0 
409 !   Entering the MD loop       
410 #ifdef MPI
411       tt0 = MPI_Wtime()
412 #else
413       tt0 = tcpu()
414 #endif
415       if (lang.eq.2 .or. lang.eq.3) then
416 #ifndef   LANG0
417         call setup_fricmat
418         if (lang.eq.2) then
419           call sd_verlet_p_setup
420         else
421           call sd_verlet_ciccotti_setup
422         endif
423         do i=1,dimen
424           do j=1,dimen
425             pfric0_mat(i,j,0)=pfric_mat(i,j)
426             afric0_mat(i,j,0)=afric_mat(i,j)
427             vfric0_mat(i,j,0)=vfric_mat(i,j)
428             prand0_mat(i,j,0)=prand_mat(i,j)
429             vrand0_mat1(i,j,0)=vrand_mat1(i,j)
430             vrand0_mat2(i,j,0)=vrand_mat2(i,j)
431           enddo
432         enddo
433         flag_stoch(0)=.true.
434         do i=1,maxflag_stoch
435           flag_stoch(i)=.false.
436         enddo
437 #else
438         write (iout,*) &
439          "LANG=2 or 3 NOT SUPPORTED. Recompile without -DLANG0"
440 #ifdef MPI
441         call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE)
442 #endif
443         stop
444 #endif
445       else if (lang.eq.1 .or. lang.eq.4) then
446         call setup_fricmat
447       endif
448       time00=MPI_WTIME()
449       if (me.eq.king .or. .not. out1file) &
450        write(iout,*) 'Setup time',time00-walltime
451       call flush(iout)
452 #ifdef MPI
453       t_langsetup=MPI_Wtime()-tt0
454       tt0=MPI_Wtime()
455 #else
456       t_langsetup=tcpu()-tt0
457       tt0=tcpu()
458 #endif
459       itime=0
460       end_of_run=.false.
461
462       do while(.not.end_of_run)
463         itime=itime+1
464         if(itime.eq.n_timestep.and.me.eq.king) end_of_run=.true.
465         if(mremdsync.and.itime.eq.n_timestep) end_of_run=.true.
466         rstcount=rstcount+1
467         if (lang.gt.0 .and. surfarea .and. &
468             mod(itime,reset_fricmat).eq.0) then
469           if (lang.eq.2 .or. lang.eq.3) then
470 #ifndef   LANG0
471             call setup_fricmat
472             if (lang.eq.2) then
473               call sd_verlet_p_setup
474             else
475               call sd_verlet_ciccotti_setup
476             endif
477             do i=1,dimen
478               do j=1,dimen
479                 pfric0_mat(i,j,0)=pfric_mat(i,j)
480                 afric0_mat(i,j,0)=afric_mat(i,j)
481                 vfric0_mat(i,j,0)=vfric_mat(i,j)
482                 prand0_mat(i,j,0)=prand_mat(i,j)
483                 vrand0_mat1(i,j,0)=vrand_mat1(i,j)
484                 vrand0_mat2(i,j,0)=vrand_mat2(i,j)
485               enddo
486             enddo
487             flag_stoch(0)=.true.
488             do i=1,maxflag_stoch
489               flag_stoch(i)=.false.
490             enddo   
491 #endif
492           else if (lang.eq.1 .or. lang.eq.4) then
493             call setup_fricmat
494           endif
495           write (iout,'(a,i10)') &
496             "Friction matrix reset based on surface area, itime",itime
497         endif
498         if (reset_vel .and. tbf .and. lang.eq.0 &
499             .and. mod(itime,count_reset_vel).eq.0) then
500           call random_vel
501           if (me.eq.king .or. .not. out1file) &
502            write(iout,'(a,f20.2)') &
503            "Velocities reset to random values, time",totT       
504           do i=0,2*nres
505             do j=1,3
506               d_t_old(j,i)=d_t(j,i)
507             enddo
508           enddo
509         endif
510         if (reset_moment .and. mod(itime,count_reset_moment).eq.0) then
511           call inertia_tensor  
512           call vcm_vel(vcm)
513           do j=1,3
514              d_t(j,0)=d_t(j,0)-vcm(j)
515           enddo
516           call kinetic(EK)
517           kinetic_T=2.0d0/(dimen3*Rb)*EK
518           scalfac=dsqrt(T_bath/kinetic_T)
519 !d          write(iout,'(a,f20.2)') "Momenta zeroed out, time",totT     
520           do i=0,2*nres
521             do j=1,3
522               d_t_old(j,i)=scalfac*d_t(j,i)
523             enddo
524           enddo
525         endif  
526         if (lang.ne.4) then
527           if (RESPA) then
528 ! Time-reversible RESPA algorithm 
529 ! (Tuckerman et al., J. Chem. Phys., 97, 1990, 1992)
530             call RESPA_step(itime)
531           else
532 ! Variable time step algorithm.
533             call velverlet_step(itime)
534           endif
535         else
536 #ifdef BROWN
537           call brown_step(itime)
538 #else
539           print *,"Brown dynamics not here!"
540 #ifdef MPI
541           call MPI_Abort(MPI_COMM_WORLD,IERROR,ERRCODE)
542 #endif
543           stop
544 #endif
545         endif
546         if(ntwe.ne.0) then
547           if (mod(itime,ntwe).eq.0) then
548               call statout(itime)
549               call enerprint(potEcomp)
550           endif
551         endif
552         if (mod(itime,ntwx).eq.0.and..not.traj1file) then
553           write (tytul,'("time",f8.2," temp",f8.1)') totT,t_bath
554           if(mdpdb) then
555              call hairpin(.true.,nharp,iharp)
556              call secondary2(.true.)
557              call pdbout(potE,tytul,ipdb)
558           else 
559              call cartout(totT)
560           endif
561         endif
562         if (mod(itime,ntwx).eq.0.and.traj1file) then
563           if(ntwx_cache.lt.max_cache_traj_use) then
564             ntwx_cache=ntwx_cache+1
565           else
566            if (max_cache_traj_use.ne.1) &
567             print *,itime,"processor ",me," over cache ",ntwx_cache
568            do i=1,ntwx_cache-1
569             call returnbox
570             totT_cache(i)=totT_cache(i+1)
571             EK_cache(i)=EK_cache(i+1)
572             potE_cache(i)=potE_cache(i+1)
573             t_bath_cache(i)=t_bath_cache(i+1)
574             Uconst_cache(i)=Uconst_cache(i+1)
575             iset_cache(i)=iset_cache(i+1)
576
577             do ii=1,nfrag
578              qfrag_cache(ii,i)=qfrag_cache(ii,i+1)
579             enddo
580             do ii=1,npair
581              qpair_cache(ii,i)=qpair_cache(ii,i+1)
582             enddo
583             do ii=1,nfrag_back
584               utheta_cache(ii,i)=utheta_cache(ii,i+1)
585               ugamma_cache(ii,i)=ugamma_cache(ii,i+1)
586               uscdiff_cache(ii,i)=uscdiff_cache(ii,i+1)
587             enddo
588
589
590             do ii=1,nres*2
591              do j=1,3
592               c_cache(j,ii,i)=c_cache(j,ii,i+1)
593              enddo
594             enddo
595            enddo
596           endif
597
598             totT_cache(ntwx_cache)=totT
599             EK_cache(ntwx_cache)=EK
600             potE_cache(ntwx_cache)=potE
601             t_bath_cache(ntwx_cache)=t_bath
602             Uconst_cache(ntwx_cache)=Uconst
603             iset_cache(ntwx_cache)=iset
604
605             do i=1,nfrag
606              qfrag_cache(i,ntwx_cache)=qfrag(i)
607             enddo
608             do i=1,npair
609              qpair_cache(i,ntwx_cache)=qpair(i)
610             enddo
611             do i=1,nfrag_back
612               utheta_cache(i,ntwx_cache)=utheta(i)
613               ugamma_cache(i,ntwx_cache)=ugamma(i)
614               uscdiff_cache(i,ntwx_cache)=uscdiff(i)
615             enddo
616 !            call returnbox
617             do i=1,nres*2
618              do j=1,3
619               c_cache(j,i,ntwx_cache)=c(j,i)
620              enddo
621             enddo
622
623         endif
624         if ((rstcount.eq.1000.or.itime.eq.n_timestep) &
625                                .and..not.restart1file) then
626
627            if(me.eq.king) then
628              open(irest1,file=mremd_rst_name,status='unknown')
629              write (irest1,*) "i2rep"
630              write (irest1,*) (i2rep(i),i=0,nodes-1)
631              write (irest1,*) "ifirst"
632              write (irest1,*) (ifirst(i),i=1,remd_m(1))
633              do il=1,nodes
634               write (irest1,*) "nupa",il
635               write (irest1,*) nupa(0,il),(nupa(i,il),i=1,nupa(0,il))
636               write (irest1,*) "ndowna",il
637               write (irest1,*) ndowna(0,il),&
638                          (ndowna(i,il),i=1,ndowna(0,il))
639              enddo
640              if(usampl) then
641               write (irest1,*) "nset"
642               write (irest1,*) nset
643               write (irest1,*) "mset"
644               write (irest1,*) (mset(i),i=1,nset)
645               write (irest1,*) "i2set"
646               write (irest1,*) (i2set(i),i=0,nodes-1)
647               write (irest1,*) "i_index"
648               do il=1,nset
649                do il1=1,mset(il)
650                 do i=1,nrep
651                   write(irest1,*) (i_index(i,j,il,il1),j=1,remd_m(i))
652                 enddo
653                enddo
654               enddo
655
656              endif
657              close(irest1)
658            endif
659            open(irest2,file=rest2name,status='unknown')
660            write(irest2,*) totT,EK,potE,totE,t_bath
661            do i=1,2*nres
662             write (irest2,'(3e15.5)') (d_t(j,i),j=1,3)
663            enddo
664            do i=1,2*nres
665             write (irest2,'(3e15.5)') (dc(j,i),j=1,3)
666            enddo
667            if(usampl) then
668              write (irest2,*) iset
669            endif
670           close(irest2)
671           rstcount=0
672         endif 
673
674 ! REMD - exchange
675 ! forced synchronization
676         if (mod(itime,i_sync_step).eq.0 .and. me.ne.king &
677                                       .and. .not. mremdsync) then 
678             synflag=.false.
679             call mpi_iprobe(0,101,CG_COMM,synflag,status,ierr)
680             if (synflag) then 
681                call mpi_recv(itime_master, 1, MPI_INTEGER,&
682                                    0,101,CG_COMM, status, ierr)
683                call mpi_barrier(CG_COMM, ierr)
684 !deb               if (out1file.or.traj1file) then
685 !deb                call mpi_gather(itime,1,mpi_integer,
686 !deb     &             icache_all,1,mpi_integer,king,
687 !deb     &             CG_COMM,ierr)                 
688                if(traj1file) &
689                 call mpi_gather(ntwx_cache,1,mpi_integer,&
690                    icache_all,1,mpi_integer,king,&
691                    CG_COMM,ierr)
692                if (.not.out1file) &
693                      write(iout,*) 'REMD synchro at3',itime_master,itime
694                if (itime_master.ge.n_timestep .or. ovrtim()) &
695                   end_of_run=.true.
696 !time               call flush(iout)
697             endif
698         endif
699
700 ! REMD - exchange
701         if ((mod(itime,nstex).eq.0.and.me.eq.king &
702                         .or.end_of_run.and.me.eq.king ) &
703              .and. .not. mremdsync ) then
704            synflag=.true.
705            do i=1,nodes-1
706               call mpi_isend(itime,1,MPI_INTEGER,i,101, &
707                                       CG_COMM, ireqi(i), ierr)
708 !d            write(iout,*) 'REMD synchro with',i
709 !d            call flush(iout)
710            enddo
711            call mpi_waitall(nodes-1,ireqi,statusi,ierr)
712            call mpi_barrier(CG_COMM, ierr)
713            time01=MPI_WTIME()
714            write(iout,*) 'REMD synchro at2',itime,'time=',time01-time00
715            if (out1file.or.traj1file) then
716 !deb            call mpi_gather(itime,1,mpi_integer,
717 !deb     &             itime_all,1,mpi_integer,king,
718 !deb     &             CG_COMM,ierr)
719 !deb            write(iout,'(a19,8000i8)') ' REMD synchro itime',
720 !deb     &                    (itime_all(i),i=1,nodes)
721             if(traj1file) then
722 !deb             imin_itime=itime_all(1)
723 !deb             do i=2,nodes
724 !deb               if(itime_all(i).lt.imin_itime) imin_itime=itime_all(i)
725 !deb             enddo
726 !deb             ii_write=(imin_itime-imin_itime_old)/ntwx
727 !deb             imin_itime_old=int(imin_itime/ntwx)*ntwx
728 !deb             write(iout,*) imin_itime,imin_itime_old,ii_write
729              call mpi_gather(ntwx_cache,1,mpi_integer,&
730                    icache_all,1,mpi_integer,king,&
731                    CG_COMM,ierr)
732 !             write(iout,'(a19,8000i8)') '     ntwx_cache',
733 !     &                    (icache_all(i),i=1,nodes)
734              ii_write=icache_all(1)
735              do i=2,nodes
736                if(icache_all(i).lt.ii_write) ii_write=icache_all(i)
737              enddo
738 !             write(iout,*) "MIN ii_write=",ii_write
739             endif
740            endif
741 !time           call flush(iout)
742         endif
743         if(mremdsync .and. mod(itime,nstex).eq.0) then
744            synflag=.true.
745            if (me.eq.king .or. .not. out1file) &
746             write(iout,*) 'REMD synchro at1',itime,ntwx_cache,Nprocs,nodes
747 !!!!!! TRIAL OF MINIM SYNC
748             if (me.eq.king) then
749            do i=1,nodes-1
750               call mpi_isend(itime,1,MPI_INTEGER,i,101, &
751                                       CG_COMM, ireqi(i), ierr)
752 !d            write(iout,*) 'REMD synchro with',i
753 !d            call flush(iout)
754            enddo
755            call mpi_waitall(nodes-1,ireqi,statusi,ierr)
756            call mpi_barrier(CG_COMM, ierr)
757            time01=MPI_WTIME()
758             endif
759 !!!!!!!!!!!!!!!!!
760             write(iout,*) icache_all
761             if(traj1file) then
762               write(iout,*) "before mpi_gather ntwx_cache"
763              call mpi_gather(ntwx_cache,1,mpi_integer,&
764                    icache_all(1),1,mpi_integer,king,&  ! CONSULT WITH ADAM
765                    CG_COMM,ierr)
766               write(iout,*) "after mpi_gather ntwx_cache"
767
768              if (me.eq.king) then
769                write(iout,'(a19,8000i8)') '     ntwx_cache',&
770                           (icache_all(i),i=1,nodes)
771                ii_write=icache_all(1)
772                do i=2,nodes
773                  if(icache_all(i).lt.ii_write) ii_write=icache_all(i)
774                enddo
775                write(iout,*) "MIN ii_write=",ii_write
776              endif
777             endif
778            call flush(iout)
779         endif
780         if (synflag) then
781 ! Update the time safety limiy
782           if (time001-time00.gt.safety) then
783             safety=time001-time00+600
784             write (iout,*) "****** SAFETY increased to",safety," s"
785           endif
786           if (ovrtim()) end_of_run=.true.
787         endif
788         if(synflag.and..not.end_of_run) then
789            time02=MPI_WTIME()
790            synflag=.false.
791
792            write(iout,*) 'REMD before',me,t_bath
793
794 !           call mpi_gather(t_bath,1,mpi_double_precision,
795 !     &             remd_t_bath,1,mpi_double_precision,king,
796 !     &             CG_COMM,ierr)
797            potEcomp(n_ene+1)=t_bath
798            if (usampl) then
799              potEcomp(n_ene+2)=iset
800              if (iset.lt.nset) then
801                i_set_temp=iset
802                iset=iset+1
803                call EconstrQ
804                potEcomp(n_ene+3)=Uconst
805                iset=i_set_temp
806              endif
807              if (iset.gt.1) then
808                i_set_temp=iset
809                iset=iset-1
810                call EconstrQ
811                potEcomp(n_ene+4)=Uconst 
812                iset=i_set_temp
813              endif
814            endif
815            call mpi_gather(potEcomp(0),n_ene+5,mpi_double_precision,&
816                    remd_ene(0,1),n_ene+5,mpi_double_precision,king,&
817                    CG_COMM,ierr)
818            if(lmuca) then 
819             call mpi_gather(elow,1,mpi_double_precision,&
820                    elowi,1,mpi_double_precision,king,&
821                    CG_COMM,ierr)
822             call mpi_gather(ehigh,1,mpi_double_precision,&
823                    ehighi,1,mpi_double_precision,king,&
824                    CG_COMM,ierr)
825            endif
826
827           time03=MPI_WTIME()
828           if (me.eq.king .or. .not. out1file) then
829             write(iout,*) 'REMD gather times=',time03-time01 &
830                                               ,time03-time02
831           endif
832
833           if (restart1file) call write1rst(i_index)
834
835           time04=MPI_WTIME()
836           if (me.eq.king .or. .not. out1file) then
837             write(iout,*) 'REMD writing rst time=',time04-time03
838           endif
839
840           if (traj1file) call write1traj
841 !d debugging
842 !deb            call mpi_gather(ntwx_cache,1,mpi_integer,
843 !deb     &             icache_all,1,mpi_integer,king,
844 !deb     &             CG_COMM,ierr)
845 !deb            write(iout,'(a19,8000i8)') '  ntwx_cache after traj1file',
846 !deb     &                    (icache_all(i),i=1,nodes)
847 !d end
848
849
850           time05=MPI_WTIME()
851           if (me.eq.king .or. .not. out1file) then
852             write(iout,*) 'REMD writing traj time=',time05-time04
853             call flush(iout)
854           endif
855
856
857           if (me.eq.king) then
858             do i=1,nodes
859                remd_t_bath(i)=remd_ene(n_ene+1,i)
860                iremd_iset(i)=remd_ene(n_ene+2,i)
861             enddo
862 #ifdef DEBUG
863             if(lmuca) then
864 !o             write(iout,*) 'REMD exchange temp,ene,elow,ehigh'
865              do i=1,nodes
866                write(iout,'(i4,4f12.5)') i,remd_t_bath(i),remd_ene(0,i),&
867                   elowi(i),ehighi(i)       
868              enddo
869             else
870               write(iout,*) 'REMD exchange temp,ene'
871               do i=1,nodes
872                 write(iout,'(i4,2f12.5)') i,remd_t_bath(i),remd_ene(0,i)
873                 write(iout,'(6f12.5)') (remd_ene(j,i),j=1,n_ene)
874               enddo
875             endif
876 #endif
877 !-------------------------------------           
878            IF(.not.usampl) THEN
879             write (iout,*) "Enter exchnge, remd_m",remd_m(1),&
880               " nodes",nodes
881             call flush(iout)
882             write (iout,*) "remd_m(1)",remd_m(1)
883             do irr=1,remd_m(1)
884                i=ifirst(iran_num(1,remd_m(1)))
885              write (iout,*) "i",i
886              call flush(iout)
887
888              do ii=1,nodes-1
889
890 #ifdef DEBUG
891               write (iout,*) "i",i," nupa(0,i)",int(nupa(0,i))
892 #endif
893              if(i.gt.0.and.nupa(0,i).gt.0) then
894               iex=i
895 !              if (i.eq.1 .and. int(nupa(0,i)).eq.1) then
896 !                write (iout,*) 
897 !     &  "CHUJ ABSOLUTNY!!! No way to sample a distinct replica in MREMD"
898 !                call flush(iout)
899 !                call MPI_Abort(MPI_COMM_WORLD,ERRCODE,ierr)
900 !              endif
901 !              do while (iex.eq.i)
902 !                write (iout,*) "upper",nupa(int(nupa(0,i)),i)
903                 iex=nupa(iran_num(1,int(nupa(0,i))),i)
904 !              enddo
905 !              write (iout,*) "nupa(0,i)",nupa(0,i)," iex",iex
906               if (lmuca) then
907                call muca_delta(remd_t_bath,remd_ene,i,iex,delta)
908               else
909 ! Swap temperatures between conformations i and iex with recalculating the free energies
910 ! following temperature changes.
911                ene_iex_iex=remd_ene(0,iex)
912                ene_i_i=remd_ene(0,i)
913 !               write (iout,*) "i",i," ene_i_i",ene_i_i,
914 !     &          " iex",iex," ene_iex_iex",ene_iex_iex
915 !               write (iout,*) "rescaling weights with temperature",
916 !     &          remd_t_bath(i)
917 !               call flush(iout)
918                call rescale_weights(remd_t_bath(i))
919
920 !               write (iout,*) "0,iex",remd_t_bath(i)
921 !               call enerprint(remd_ene(0,iex))
922
923                call sum_energy(remd_ene(0,iex),.false.)
924                ene_iex_i=remd_ene(0,iex)
925 !               write (iout,*) "ene_iex_i",remd_ene(0,iex)
926
927 !               write (iout,*) "0,i",remd_t_bath(i)
928 !               call enerprint(remd_ene(0,i))
929
930                call sum_energy(remd_ene(0,i),.false.)
931 !               write (iout,*) "ene_i_i",remd_ene(0,i)
932 !               call flush(iout)
933 !               write (iout,*) "rescaling weights with temperature",
934 !     &          remd_t_bath(iex)
935                if (real(ene_i_i).ne.real(remd_ene(0,i))) then
936                 write (iout,*) "ERROR: inconsistent energies:",i,&
937                   ene_i_i,remd_ene(0,i)
938                endif
939                call rescale_weights(remd_t_bath(iex))
940
941 !               write (iout,*) "0,i",remd_t_bath(iex)
942 !               call enerprint(remd_ene(0,i))
943
944                call sum_energy(remd_ene(0,i),.false.)
945 !               write (iout,*) "ene_i_iex",remd_ene(0,i)
946 !               call flush(iout)
947                ene_i_iex=remd_ene(0,i)
948
949 !               write (iout,*) "0,iex",remd_t_bath(iex)
950                call enerprint(remd_ene(0,iex))
951
952                call sum_energy(remd_ene(0,iex),.false.)
953                if (real(ene_iex_iex).ne.real(remd_ene(0,iex))) then
954                 write (iout,*) "ERROR: inconsistent energies:",iex,&
955                   ene_iex_iex,remd_ene(0,iex)
956                endif
957 !               write (iout,*) "ene_iex_iex",remd_ene(0,iex)
958 !               write (iout,*) "i",i," iex",iex
959 !               write (iout,'(4(a,e15.5))') "ene_i_i",ene_i_i,
960 !     &           " ene_i_iex",ene_i_iex,
961 !     &           " ene_iex_i",ene_iex_i," ene_iex_iex",ene_iex_iex
962 !               call flush(iout)
963                delta=(ene_iex_iex-ene_i_iex)/(Rb*remd_t_bath(iex))- &
964                     (ene_iex_i-ene_i_i)/(Rb*remd_t_bath(i))
965                delta=-delta
966 !               write(iout,*) 'delta',delta
967 !              delta=(remd_t_bath(i)-remd_t_bath(iex))*
968 !     &              (remd_ene(i)-remd_ene(iex))/Rb/
969 !     &              (remd_t_bath(i)*remd_t_bath(iex))
970               endif
971               if (delta .gt. 50.0d0) then
972                 delta=0.0d0
973               else
974 #ifdef OSF 
975                 if(isnan(delta))then
976                   delta=0.0d0
977                 else if (delta.lt.-50.0d0) then
978                   delta=dexp(50.0d0)
979                 else
980                   delta=dexp(-delta)
981                 endif
982 #else
983                 delta=dexp(-delta)
984 #endif
985               endif
986               iremd_tot(int(i2rep(i-1)))=iremd_tot(int(i2rep(i-1)))+1
987               xxx=ran_number(0.0d0,1.0d0)
988 !              write(iout,'(2i4,a6,2f12.5)') i,iex,' delta',delta,xxx
989 !              call flush(iout)
990               if (delta .gt. xxx) then
991                 tmp=remd_t_bath(i)       
992                 remd_t_bath(i)=remd_t_bath(iex)
993                 remd_t_bath(iex)=tmp
994                 remd_ene(0,i)=ene_i_iex
995                 remd_ene(0,iex)=ene_iex_i
996                 if(lmuca) then
997                   tmp=elowi(i)
998                   elowi(i)=elowi(iex)
999                   elowi(iex)=tmp  
1000                   tmp=ehighi(i)
1001                   ehighi(i)=ehighi(iex)
1002                   ehighi(iex)=tmp  
1003                 endif
1004
1005
1006                 do k=0,nodes
1007                   itmp=nupa(k,i)
1008                   nupa(k,i)=nupa(k,iex)
1009                   nupa(k,iex)=itmp
1010                   itmp=ndowna(k,i)
1011                   ndowna(k,i)=ndowna(k,iex)
1012                   ndowna(k,iex)=itmp
1013                 enddo
1014                 do il=1,nodes
1015                  if (ifirst(il).eq.i) ifirst(il)=iex
1016                  do k=1,nupa(0,il)
1017                   if (nupa(k,il).eq.i) then 
1018                      nupa(k,il)=iex
1019                   elseif (nupa(k,il).eq.iex) then 
1020                      nupa(k,il)=i
1021                   endif
1022                  enddo
1023                  do k=1,ndowna(0,il)
1024                   if (ndowna(k,il).eq.i) then 
1025                      ndowna(k,il)=iex
1026                   elseif (ndowna(k,il).eq.iex) then 
1027                      ndowna(k,il)=i
1028                   endif
1029                  enddo
1030                 enddo
1031
1032                 iremd_acc(int(i2rep(i-1)))=iremd_acc(int(i2rep(i-1)))+1
1033                 itmp=i2rep(i-1)
1034                 i2rep(i-1)=i2rep(iex-1)
1035                 i2rep(iex-1)=itmp
1036
1037 !                write(iout,*) 'exchange',i,iex
1038 !                write (iout,'(a8,100i4)') "@ ifirst",
1039 !     &                    (ifirst(k),k=1,remd_m(1))
1040 !                do il=1,nodes
1041 !                 write (iout,'(a8,i4,a1,100i4)') "@ nupa",il,":",
1042 !     &                    (nupa(k,il),k=1,nupa(0,il))
1043 !                 write (iout,'(a8,i4,a1,100i4)') "@ ndowna",il,":",
1044 !     &                    (ndowna(k,il),k=1,ndowna(0,il))
1045 !                enddo
1046 !                call flush(iout) 
1047
1048               else
1049                remd_ene(0,iex)=ene_iex_iex
1050                remd_ene(0,i)=ene_i_i
1051                i=iex
1052               endif 
1053             endif
1054            enddo
1055            enddo
1056 !d           write (iout,*) "exchange completed"
1057 !d           call flush(iout) 
1058         ELSE
1059           do ii=1,nodes  
1060 !d            write(iout,*) "########",ii
1061
1062             i_temp=iran_num(1,nrep)
1063             i_mult=iran_num(1,remd_m(i_temp))
1064             i_iset=iran_num(1,nset)
1065             i_mset=iran_num(1,mset(i_iset))
1066             i=i_index(i_temp,i_mult,i_iset,i_mset)
1067
1068 !d            write(iout,*) "i=",i,i_temp,i_mult,i_iset,i_mset
1069
1070             i_dir=iran_num(1,3)
1071 !d            write(iout,*) "i_dir=",i_dir
1072
1073             if(i_dir.eq.1 .and. remd_m(i_temp+1).gt.0 )then            
1074                
1075                i_temp1=i_temp+1
1076                i_mult1=iran_num(1,remd_m(i_temp1))
1077                i_iset1=i_iset
1078                i_mset1=iran_num(1,mset(i_iset1))
1079                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1080
1081             elseif(i_dir.eq.2 .and. mset(i_iset+1).gt.0)then
1082
1083                i_temp1=i_temp
1084                i_mult1=iran_num(1,remd_m(i_temp1))
1085                i_iset1=i_iset+1
1086                i_mset1=iran_num(1,mset(i_iset1))
1087                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1088                econstr_temp_i=remd_ene(20,i)
1089                econstr_temp_iex=remd_ene(20,iex)
1090                remd_ene(20,i)=remd_ene(n_ene+3,i)
1091                remd_ene(20,iex)=remd_ene(n_ene+4,iex)
1092
1093             elseif(remd_m(i_temp+1).gt.0.and.mset(i_iset+1).gt.0)then
1094
1095                i_temp1=i_temp+1
1096                i_mult1=iran_num(1,remd_m(i_temp1))
1097                i_iset1=i_iset+1
1098                i_mset1=iran_num(1,mset(i_iset1))
1099                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1100                econstr_temp_i=remd_ene(20,i)
1101                econstr_temp_iex=remd_ene(20,iex)
1102                remd_ene(20,i)=remd_ene(n_ene+3,i)
1103                remd_ene(20,iex)=remd_ene(n_ene+4,iex)
1104
1105             else
1106                goto 444 
1107             endif
1108  
1109 !d            write(iout,*) "iex=",iex,i_temp1,i_mult1,i_iset1,i_mset1
1110             call flush(iout)
1111
1112 ! Swap temperatures between conformations i and iex with recalculating the free energies
1113 ! following temperature changes.
1114               ene_iex_iex=remd_ene(0,iex)
1115               ene_i_i=remd_ene(0,i)
1116 !o              write (iout,*) "rescaling weights with temperature",
1117 !o     &          remd_t_bath(i)
1118               call rescale_weights(remd_t_bath(i))
1119               
1120               call sum_energy(remd_ene(0,iex),.false.)
1121               ene_iex_i=remd_ene(0,iex)
1122 !d              write (iout,*) "ene_iex_i",remd_ene(0,iex)
1123 !              call sum_energy(remd_ene(0,i),.false.)
1124 !d              write (iout,*) "ene_i_i",remd_ene(0,i)
1125 !              write (iout,*) "rescaling weights with temperature",
1126 !     &          remd_t_bath(iex)
1127 !              if (real(ene_i_i).ne.real(remd_ene(0,i))) then
1128 !                write (iout,*) "ERROR: inconsistent energies:",i,
1129 !     &            ene_i_i,remd_ene(0,i)
1130 !              endif
1131               call rescale_weights(remd_t_bath(iex))
1132               call sum_energy(remd_ene(0,i),.false.)
1133 !d              write (iout,*) "ene_i_iex",remd_ene(0,i)
1134               ene_i_iex=remd_ene(0,i)
1135 !              call sum_energy(remd_ene(0,iex),.false.)
1136 !              if (real(ene_iex_iex).ne.real(remd_ene(0,iex))) then
1137 !                write (iout,*) "ERROR: inconsistent energies:",iex,
1138 !     &            ene_iex_iex,remd_ene(0,iex)
1139 !              endif
1140 !d              write (iout,*) "ene_iex_iex",remd_ene(0,iex)
1141 !              write (iout,*) "i",i," iex",iex
1142 !d              write (iout,'(4(a,e15.5))') "ene_i_i",ene_i_i,
1143 !d     &           " ene_i_iex",ene_i_iex,
1144 !d     &           " ene_iex_i",ene_iex_i," ene_iex_iex",ene_iex_iex
1145               delta=(ene_iex_iex-ene_i_iex)/(Rb*remd_t_bath(iex))- &
1146                     (ene_iex_i-ene_i_i)/(Rb*remd_t_bath(i))
1147               delta=-delta
1148 !d              write(iout,*) 'delta',delta
1149 !              delta=(remd_t_bath(i)-remd_t_bath(iex))*
1150 !     &              (remd_ene(i)-remd_ene(iex))/Rb/
1151 !     &              (remd_t_bath(i)*remd_t_bath(iex))
1152               if (delta .gt. 50.0d0) then
1153                 delta=0.0d0
1154               else
1155                 delta=dexp(-delta)
1156               endif
1157               if (i_dir.eq.1.or.i_dir.eq.3) &
1158                iremd_tot(int(i2rep(i-1)))=iremd_tot(int(i2rep(i-1)))+1
1159               if (i_dir.eq.2.or.i_dir.eq.3) &
1160                 iremd_tot_usa(int(i2set(i-1)))= &
1161                        iremd_tot_usa(int(i2set(i-1)))+1
1162               xxx=ran_number(0.0d0,1.0d0)
1163 !d              write(iout,'(2i4,a6,2f12.5)') i,iex,' delta',delta,xxx
1164               if (delta .gt. xxx) then
1165                 tmp=remd_t_bath(i)       
1166                 remd_t_bath(i)=remd_t_bath(iex)
1167                 remd_t_bath(iex)=tmp
1168
1169                 itmp=iremd_iset(i)       
1170                 iremd_iset(i)=iremd_iset(iex)
1171                 iremd_iset(iex)=itmp
1172
1173                 remd_ene(0,i)=ene_i_iex
1174                 remd_ene(0,iex)=ene_iex_i
1175
1176                 if (i_dir.eq.1.or.i_dir.eq.3) &
1177                  iremd_acc(int(i2rep(i-1)))=iremd_acc(int(i2rep(i-1)))+1
1178
1179                 itmp=i2rep(i-1)
1180                 i2rep(i-1)=i2rep(iex-1)
1181                 i2rep(iex-1)=itmp
1182
1183                 if (i_dir.eq.2.or.i_dir.eq.3) &
1184                  iremd_acc_usa(int(i2set(i-1)))= &
1185                        iremd_acc_usa(int(i2set(i-1)))+1
1186
1187                 itmp=i2set(i-1)
1188                 i2set(i-1)=i2set(iex-1)
1189                 i2set(iex-1)=itmp
1190         
1191                 itmp=i_index(i_temp,i_mult,i_iset,i_mset)
1192                 i_index(i_temp,i_mult,i_iset,i_mset)= &
1193                       i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1194                 i_index(i_temp1,i_mult1,i_iset1,i_mset1)=itmp
1195
1196               else
1197                remd_ene(0,iex)=ene_iex_iex
1198                remd_ene(0,i)=ene_i_i
1199                remd_ene(20,iex)=econstr_temp_iex
1200                remd_ene(20,i)=econstr_temp_i
1201               endif
1202
1203 !d      do il=1,nset
1204 !d       do il1=1,mset(il)
1205 !d        do i=1,nrep
1206 !d         do j=1,remd_m(i)
1207 !d          write(iout,*) i,j,il,il1,i_index(i,j,il,il1)
1208 !d         enddo
1209 !d        enddo
1210 !d       enddo
1211 !d      enddo
1212
1213  444      continue           
1214
1215           enddo
1216
1217
1218         ENDIF
1219
1220 !-------------------------------------
1221              write (iout,*) "NREP",nrep
1222              do i=1,nrep
1223               if(iremd_tot(i).ne.0) &
1224                 write(iout,'(a3,i4,2f12.5,i5)') 'ACC',i,remd_t(i) &
1225                  ,iremd_acc(i)/(1.0*iremd_tot(i)),iremd_tot(i)
1226              enddo
1227
1228              if(usampl) then
1229               do i=1,nset
1230                if(iremd_tot_usa(i).ne.0) &
1231                  write(iout,'(a10,i4,f12.5,i8)') 'ACC_usampl',i,&
1232                iremd_acc_usa(i)/(1.0*iremd_tot_usa(i)),iremd_tot_usa(i)
1233               enddo
1234              endif
1235
1236              call flush(iout)
1237
1238 !d              write (iout,'(a6,100i4)') "ifirst",
1239 !d     &                    (ifirst(i),i=1,remd_m(1))
1240 !d              do il=1,nodes
1241 !d               write (iout,'(a5,i4,a1,100i4)') "nup",il,":",
1242 !d     &                    (nupa(i,il),i=1,nupa(0,il))
1243 !d               write (iout,'(a5,i4,a1,100i4)') "ndown",il,":",
1244 !d     &                    (ndowna(i,il),i=1,ndowna(0,il))
1245 !d              enddo
1246             endif
1247
1248          time06=MPI_WTIME()
1249 !d         write (iout,*) "Before scatter"
1250 !d         call flush(iout)
1251 #ifdef DEBUG
1252          if (me.eq.king) then
1253            write (iout,*) "t_bath before scatter",remd_t_bath
1254            call flush(iout)
1255          endif
1256 #endif
1257          call mpi_scatter(remd_t_bath,1,mpi_double_precision,&
1258                  t_bath,1,mpi_double_precision,king,&
1259                  CG_COMM,ierr) 
1260 !d         write (iout,*) "After scatter"
1261 !d         call flush(iout)
1262          if(usampl) &
1263           call mpi_scatter(iremd_iset,1,mpi_integer,&
1264                  iset,1,mpi_integer,king,&
1265                  CG_COMM,ierr) 
1266
1267          time07=MPI_WTIME()
1268           if (me.eq.king .or. .not. out1file) then
1269             write(iout,*) 'REMD scatter time=',time07-time06
1270           endif
1271
1272          if(lmuca) then
1273            call mpi_scatter(elowi,1,mpi_double_precision,&
1274                  elow,1,mpi_double_precision,king,&
1275                  CG_COMM,ierr) 
1276            call mpi_scatter(ehighi,1,mpi_double_precision,&
1277                  ehigh,1,mpi_double_precision,king,&
1278                  CG_COMM,ierr) 
1279          endif
1280          call rescale_weights(t_bath)
1281 !o         write (iout,*) "Processor",me,
1282 !o     &    " rescaling weights with temperature",t_bath
1283          do j=1,5
1284          stdfp(j)=dsqrt(2*Rb*t_bath/d_time)
1285          do i=1,ntyp
1286            stdfsc(i,j)=dsqrt(2*Rb*t_bath/d_time)
1287          enddo 
1288          enddo
1289
1290 !de         write(iout,*) 'REMD after',me,t_bath
1291            time08=MPI_WTIME()
1292            if (me.eq.king .or. .not. out1file) then
1293             write(iout,*) 'REMD exchange time 8-0=',time08-time00
1294             write(iout,*) 'REMD exchange time 8-7=',time08-time07
1295             write(iout,*) 'REMD exchange time 7-6=',time07-time06
1296             write(iout,*) 'REMD exchange time 6-5=',time06-time05
1297             write(iout,*) 'REMD exchange time 5-4=',time05-time04
1298             write(iout,*) 'REMD exchange time 4-3=',time04-time03
1299             write(iout,*) 'REMD exchange time 3-2=',time03-time02
1300             write(iout,*) 'REMD exchange time 2-1=',time02-time01
1301             write(iout,*) 'REMD exchange time 1-0=',time01-time00
1302             call flush(iout)
1303            endif
1304         endif
1305       enddo
1306
1307       if (restart1file) then 
1308           if (me.eq.king .or. .not. out1file) &
1309             write(iout,*) 'writing restart at the end of run'
1310            call write1rst(i_index)
1311       endif
1312
1313       if (traj1file) call write1traj
1314 !d debugging
1315 !deb            call mpi_gather(ntwx_cache,1,mpi_integer,
1316 !deb     &             icache_all,1,mpi_integer,king,
1317 !deb     &             CG_COMM,ierr)
1318 !deb            write(iout,'(a40,8000i8)') 
1319 !deb     &             '  ntwx_cache after traj1file at the end',
1320 !deb     &             (icache_all(i),i=1,nodes)
1321 !d end
1322
1323
1324 #ifdef MPI
1325       t_MD=MPI_Wtime()-tt0
1326 #else
1327       t_MD=tcpu()-tt0
1328 #endif
1329       if (me.eq.king .or. .not. out1file) then
1330        write (iout,'(//35(1h=),a10,35(1h=)/10(/a40,1pe15.5))') &
1331         '  Timing  ',&
1332        'MD calculations setup:',t_MDsetup,&
1333        'Energy & gradient evaluation:',t_enegrad,&
1334        'Stochastic MD setup:',t_langsetup,&
1335        'Stochastic MD step setup:',t_sdsetup,&
1336        'MD steps:',t_MD
1337        write (iout,'(/28(1h=),a25,27(1h=))') &
1338        '  End of MD calculation  '
1339       endif
1340 !el  common /przechowalnia/
1341 !      deallocate(d_restart1)
1342 !      deallocate(d_restart2)
1343 !      deallocate(p_c)
1344 !el--------------
1345       return
1346       end subroutine MREMD
1347 !-----------------------------------------------------------------------------
1348       subroutine write1rst(i_index)
1349
1350       use control_data
1351 !      implicit real*8 (a-h,o-z)
1352 !      include 'DIMENSIONS'
1353       include 'mpif.h'
1354 !      include 'COMMON.MD'
1355 !      include 'COMMON.IOUNITS'
1356 !      include 'COMMON.REMD'
1357 !      include 'COMMON.SETUP'
1358 !      include 'COMMON.CHAIN'
1359 !      include 'COMMON.SBRIDGE'
1360 !      include 'COMMON.INTERACT'
1361                
1362 !el      real(kind=4) :: d_restart1(3,2*nres*maxprocs),&
1363 !el           d_restart2(3,2*nres*maxprocs)
1364       real(kind=4) :: r_d(3,0:2*nres)
1365       real(kind=4) :: t5_restart1(5)
1366       integer :: iret,itmp
1367 !      integer(kind=2) :: i_index(Nprocs/4,Nprocs/20,Nprocs/200,Nprocs/200)
1368       integer(kind=2) :: i_index(Nprocs,Nprocs,Nprocs,Nprocs)
1369
1370       !(maxprocs/4,maxprocs/20,maxprocs/200,maxprocs/200)
1371 !el       common /przechowalnia/ d_restart1,d_restart2
1372       integer :: i,j,il,il1,ierr,ixdrf
1373
1374        t5_restart1(1)=totT
1375        t5_restart1(2)=EK
1376        t5_restart1(3)=potE
1377        t5_restart1(4)=t_bath
1378        t5_restart1(5)=Uconst
1379        
1380        call mpi_gather(t5_restart1,5,mpi_real,&
1381             t_restart1,5,mpi_real,king,CG_COMM,ierr)
1382
1383
1384        do i=0,2*nres
1385          do j=1,3
1386            r_d(j,i)=d_t(j,i)
1387          enddo
1388        enddo
1389        call mpi_gather(r_d,3*2*nres+3,mpi_real,&
1390                  d_restart1,3*2*nres+3,mpi_real,king,&
1391                  CG_COMM,ierr)
1392        do j=1,3
1393        dc(j,0)=c(j,1)
1394        enddo
1395
1396
1397        do i=0,2*nres
1398          do j=1,3
1399            r_d(j,i)=dc(j,i)
1400          enddo
1401        enddo
1402        call mpi_gather(r_d,3*2*nres+3,mpi_real,&
1403                  d_restart2,3*2*nres+3,mpi_real,king,&
1404                  CG_COMM,ierr)
1405
1406        if(me.eq.king) then
1407 #ifdef AIX
1408          call xdrfopen_(ixdrf,mremd_rst_name, "w", iret)
1409          do i=0,nodes-1
1410           call xdrfint_(ixdrf, i2rep(i), iret)
1411          enddo
1412          do i=1,remd_m(1)
1413           call xdrfint_(ixdrf, ifirst(i), iret)
1414          enddo
1415          do il=1,nodes
1416               do i=0,nupa(0,il)
1417                call xdrfint_(ixdrf, nupa(i,il), iret)
1418               enddo
1419
1420               do i=0,ndowna(0,il)
1421                call xdrfint_(ixdrf, ndowna(i,il), iret)
1422               enddo
1423          enddo
1424
1425          do il=1,nodes
1426            do j=1,4
1427             call xdrffloat_(ixdrf, t_restart1(j,il), iret)
1428            enddo
1429          enddo
1430
1431          do il=0,nodes-1
1432            do i=0,2*nres
1433             do j=1,3
1434              call xdrffloat_(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1435             enddo
1436            enddo
1437          enddo
1438          do il=0,nodes-1
1439            do i=0,2*nres
1440             do j=1,3
1441              call xdrffloat_(ixdrf, d_restart2(j,i+(2*nres+1)*il), iret)
1442             enddo
1443            enddo
1444          enddo
1445
1446          if(usampl) then
1447            call xdrfint_(ixdrf, nset, iret)
1448            do i=1,nset
1449              call xdrfint_(ixdrf,mset(i), iret)
1450            enddo
1451            do i=0,nodes-1
1452              call xdrfint_(ixdrf,i2set(i), iret)
1453            enddo
1454            do il=1,nset
1455              do il1=1,mset(il)
1456                do i=1,nrep
1457                  do j=1,remd_m(i)
1458                    itmp=i_index(i,j,il,il1)
1459                    call xdrfint_(ixdrf,itmp, iret)
1460                  enddo
1461                enddo
1462              enddo
1463            enddo
1464            
1465          endif
1466          call xdrfclose_(ixdrf, iret)
1467 #else
1468          call xdrfopen(ixdrf,mremd_rst_name, "w", iret)
1469          do i=0,nodes-1
1470           call xdrfint(ixdrf, i2rep(i), iret)
1471          enddo
1472          do i=1,remd_m(1)
1473           call xdrfint(ixdrf, ifirst(i), iret)
1474          enddo
1475          do il=1,nodes
1476               do i=0,nupa(0,il)
1477                call xdrfint(ixdrf, nupa(i,il), iret)
1478               enddo
1479
1480               do i=0,ndowna(0,il)
1481                call xdrfint(ixdrf, ndowna(i,il), iret)
1482               enddo
1483          enddo
1484
1485          do il=1,nodes
1486            do j=1,4
1487             call xdrffloat(ixdrf, t_restart1(j,il), iret)
1488            enddo
1489          enddo
1490
1491          do il=0,nodes-1
1492            do i=0,2*nres
1493             do j=1,3
1494              call xdrffloat(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1495             enddo
1496            enddo
1497          enddo
1498          do il=0,nodes-1
1499            do i=0,2*nres
1500             do j=1,3
1501              call xdrffloat(ixdrf, d_restart2(j,i+(2*nres+1)*il), iret)
1502             enddo
1503            enddo
1504          enddo
1505
1506
1507              if(usampl) then
1508               call xdrfint(ixdrf, nset, iret)
1509               do i=1,nset
1510                 call xdrfint(ixdrf,mset(i), iret)
1511               enddo
1512               do i=0,nodes-1
1513                 call xdrfint(ixdrf,i2set(i), iret)
1514               enddo
1515               do il=1,nset
1516                do il1=1,mset(il)
1517                 do i=1,nrep
1518                  do j=1,remd_m(i)
1519                    itmp=i_index(i,j,il,il1)
1520                    call xdrfint(ixdrf,itmp, iret)
1521                  enddo
1522                 enddo
1523                enddo
1524               enddo
1525            
1526              endif
1527          call xdrfclose(ixdrf, iret)
1528 #endif
1529        endif
1530       return
1531       end subroutine  write1rst
1532 !-----------------------------------------------------------------------------
1533       subroutine write1traj
1534
1535 !      implicit real*8 (a-h,o-z)
1536 !      include 'DIMENSIONS'
1537       include 'mpif.h'
1538 !      include 'COMMON.MD'
1539 !      include 'COMMON.IOUNITS'
1540 !      include 'COMMON.REMD'
1541 !      include 'COMMON.SETUP'
1542 !      include 'COMMON.CHAIN'
1543 !      include 'COMMON.SBRIDGE'
1544 !      include 'COMMON.INTERACT'
1545                
1546       real(kind=4) :: t5_restart1(5)
1547       integer :: iret,itmp
1548       real(kind=4) :: xcoord(3,2*nres+2),prec
1549       real(kind=4) :: r_qfrag(50),r_qpair(100)
1550       real(kind=4) :: r_utheta(50),r_ugamma(100),r_uscdiff(100)
1551       real(kind=4) :: p_qfrag(50*Nprocs),p_qpair(100*Nprocs) !(100*maxprocs)
1552       real(kind=4) :: p_utheta(50*Nprocs),p_ugamma(100*Nprocs),&
1553            p_uscdiff(100*Nprocs) !(100*maxprocs)
1554 !el      real(kind=4) :: p_c(3,(nres2+2)*maxprocs)
1555       real(kind=4) :: r_c(3,2*nres+2)
1556 !el      common /przechowalnia/ p_c
1557
1558       integer :: i,j,il,ierr,ii,ixdrf
1559
1560       call mpi_bcast(ii_write,1,mpi_integer,&
1561                  king,CG_COMM,ierr)
1562
1563 ! debugging
1564 !      print *,'traj1file',me,ii_write,ntwx_cache
1565 ! end debugging
1566
1567 #ifdef AIX
1568       if(me.eq.king) call xdrfopen_(ixdrf,cartname, "a", iret)
1569 #else
1570       if(me.eq.king) call xdrfopen(ixdrf,cartname, "a", iret)
1571 #endif
1572       do ii=1,ii_write
1573 !       write (iout,*) "before gather write1traj: from node",ii
1574 !       call flush(iout)
1575 !       write (iout,*) totT_cache(ii),EK_cache(ii),potE_cache(ii),t_bath_cache(ii),Uconst_cache(ii)
1576 !       call flush(iout)
1577        t5_restart1(1)=totT_cache(ii)
1578        t5_restart1(2)=EK_cache(ii)
1579        t5_restart1(3)=potE_cache(ii)
1580        t5_restart1(4)=t_bath_cache(ii)
1581        t5_restart1(5)=Uconst_cache(ii)
1582 !       write (iout,*) "before gather write1traj: from node",ii,t5_restart1(1),t5_restart1(3),t5_restart1(5),t5_restart1(4)
1583        call flush(iout)
1584        call mpi_gather(t5_restart1,5,mpi_real,&
1585             t_restart1,5,mpi_real,king,CG_COMM,ierr)
1586 !       do il=1,nodes
1587 !       write (iout,*) "after gather write1traj: from node",il,t_restart1(1,il),t_restart1(3,il),t_restart1(5,il),t_restart1(4,il)
1588 !       enddo
1589
1590        call mpi_gather(iset_cache(ii),1,mpi_integer,&
1591             iset_restart1,1,mpi_integer,king,CG_COMM,ierr)
1592
1593           do i=1,nfrag
1594            r_qfrag(i)=qfrag_cache(i,ii)
1595           enddo
1596           do i=1,npair
1597            r_qpair(i)=qpair_cache(i,ii)
1598           enddo
1599           do i=1,nfrag_back
1600            r_utheta(i)=utheta_cache(i,ii)
1601            r_ugamma(i)=ugamma_cache(i,ii)
1602            r_uscdiff(i)=uscdiff_cache(i,ii)
1603           enddo
1604
1605         call mpi_gather(r_qfrag,nfrag,mpi_real,&
1606                  p_qfrag,nfrag,mpi_real,king,&
1607                  CG_COMM,ierr)
1608          call mpi_gather(r_qpair,npair,mpi_real,&
1609                 p_qpair,npair,mpi_real,king,&
1610                  CG_COMM,ierr)
1611          call mpi_gather(r_utheta,nfrag_back,mpi_real,&
1612                 p_utheta,nfrag_back,mpi_real,king,&
1613                  CG_COMM,ierr)
1614         call mpi_gather(r_ugamma,nfrag_back,mpi_real,&
1615                  p_ugamma,nfrag_back,mpi_real,king,&
1616                  CG_COMM,ierr)
1617         call mpi_gather(r_uscdiff,nfrag_back,mpi_real,&
1618                  p_uscdiff,nfrag_back,mpi_real,king,&
1619                  CG_COMM,ierr)
1620
1621 #ifdef DEBUG
1622         write (iout,*) "p_qfrag"
1623         do i=1,nodes
1624           write (iout,*) i,(p_qfrag((i-1)*nfrag+j),j=1,nfrag)
1625         enddo
1626         write (iout,*) "p_qpair"
1627         do i=1,nodes
1628           write (iout,*) i,(p_qpair((i-1)*npair+j),j=1,npair)
1629         enddo
1630         call flush(iout)
1631 #endif
1632         do i=1,nres*2
1633          do j=1,3
1634           r_c(j,i)=c_cache(j,i,ii)
1635          enddo
1636         enddo
1637
1638         call mpi_gather(r_c,3*2*nres,mpi_real,&
1639                  p_c,3*2*nres,mpi_real,king,&
1640                  CG_COMM,ierr)
1641
1642        if(me.eq.king) then
1643 #ifdef AIX
1644          do il=1,nodes
1645           call xdrffloat_(ixdrf, real(t_restart1(1,il)), iret)
1646           call xdrffloat_(ixdrf, real(t_restart1(3,il)), iret)
1647           call xdrffloat_(ixdrf, real(t_restart1(5,il)), iret)
1648           call xdrffloat_(ixdrf, real(t_restart1(4,il)), iret)
1649           call xdrfint_(ixdrf, nss, iret) 
1650           do j=1,nss
1651            if (dyn_ss) then
1652             call xdrfint(ixdrf, idssb(j)+nres, iret)
1653             call xdrfint(ixdrf, jdssb(j)+nres, iret)
1654            else
1655             call xdrfint_(ixdrf, ihpb(j), iret)
1656             call xdrfint_(ixdrf, jhpb(j), iret)
1657            endif
1658           enddo
1659           call xdrfint_(ixdrf, nfrag+npair+3*nfrag_back, iret)
1660           call xdrfint_(ixdrf, iset_restart1(il), iret)
1661           do i=1,nfrag
1662            call xdrffloat_(ixdrf, p_qfrag(i+(il-1)*nfrag), iret)
1663           enddo
1664           do i=1,npair
1665            call xdrffloat_(ixdrf, p_qpair(i+(il-1)*npair), iret)
1666           enddo
1667           do i=1,nfrag_back
1668            call xdrffloat_(ixdrf, p_utheta(i+(il-1)*nfrag_back), iret)
1669            call xdrffloat_(ixdrf, p_ugamma(i+(il-1)*nfrag_back), iret)
1670            call xdrffloat_(ixdrf, p_uscdiff(i+(il-1)*nfrag_back), iret)
1671           enddo
1672           prec=10000.0
1673           do i=1,nres
1674            do j=1,3
1675             xcoord(j,i)=p_c(j,i+(il-1)*nres*2)
1676            enddo
1677           enddo
1678           do i=nnt,nct
1679            do j=1,3
1680             xcoord(j,nres+i-nnt+1)=p_c(j,i+nres+(il-1)*nres*2)
1681            enddo
1682           enddo
1683           itmp=nres+nct-nnt+1
1684           call xdrf3dfcoord_(ixdrf, xcoord, itmp, prec, iret)
1685          enddo
1686 #else
1687          do il=1,nodes
1688           call xdrffloat(ixdrf, real(t_restart1(1,il)), iret)
1689           call xdrffloat(ixdrf, real(t_restart1(3,il)), iret)
1690           call xdrffloat(ixdrf, real(t_restart1(5,il)), iret)
1691           call xdrffloat(ixdrf, real(t_restart1(4,il)), iret)
1692 !          write (iout,*) "write1traj: from node",ii,t_restart1(1,il),t_restart1(3,il),t_restart1(5,il),t_restart1(4,il)
1693           call xdrfint(ixdrf, nss, iret) 
1694           do j=1,nss
1695            if (dyn_ss) then
1696             call xdrfint(ixdrf, idssb(j)+nres, iret)
1697             call xdrfint(ixdrf, jdssb(j)+nres, iret)
1698            else
1699             call xdrfint(ixdrf, ihpb(j), iret)
1700             call xdrfint(ixdrf, jhpb(j), iret)
1701            endif
1702           enddo
1703           call xdrfint(ixdrf, nfrag+npair+3*nfrag_back, iret)
1704           call xdrfint(ixdrf, iset_restart1(il), iret)
1705           do i=1,nfrag
1706            call xdrffloat(ixdrf, p_qfrag(i+(il-1)*nfrag), iret)
1707           enddo
1708           do i=1,npair
1709            call xdrffloat(ixdrf, p_qpair(i+(il-1)*npair), iret)
1710           enddo
1711           do i=1,nfrag_back
1712            call xdrffloat(ixdrf, p_utheta(i+(il-1)*nfrag_back), iret)
1713            call xdrffloat(ixdrf, p_ugamma(i+(il-1)*nfrag_back), iret)
1714            call xdrffloat(ixdrf, p_uscdiff(i+(il-1)*nfrag_back), iret)
1715           enddo
1716           prec=10000.0
1717           do i=1,nres
1718            do j=1,3
1719             xcoord(j,i)=p_c(j,i+(il-1)*nres*2)
1720            enddo
1721           enddo
1722           do i=nnt,nct
1723            do j=1,3
1724             xcoord(j,nres+i-nnt+1)=p_c(j,i+nres+(il-1)*nres*2)
1725            enddo
1726           enddo
1727           itmp=nres+nct-nnt+1
1728           call xdrf3dfcoord(ixdrf, xcoord, itmp, prec, iret)
1729          enddo
1730 #endif
1731        endif
1732       enddo
1733 #ifdef AIX
1734       if(me.eq.king) call xdrfclose_(ixdrf, iret)
1735 #else
1736       if(me.eq.king) call xdrfclose(ixdrf, iret)
1737 #endif
1738       do i=1,ntwx_cache-ii_write
1739
1740             totT_cache(i)=totT_cache(ii_write+i)
1741             EK_cache(i)=EK_cache(ii_write+i)
1742             potE_cache(i)=potE_cache(ii_write+i)
1743             t_bath_cache(i)=t_bath_cache(ii_write+i)
1744             Uconst_cache(i)=Uconst_cache(ii_write+i)
1745             iset_cache(i)=iset_cache(ii_write+i)
1746
1747             do ii=1,nfrag
1748              qfrag_cache(ii,i)=qfrag_cache(ii,ii_write+i)
1749             enddo
1750             do ii=1,npair
1751              qpair_cache(ii,i)=qpair_cache(ii,ii_write+i)
1752             enddo
1753             do ii=1,nfrag_back
1754               utheta_cache(ii,i)=utheta_cache(ii,ii_write+i)
1755               ugamma_cache(ii,i)=ugamma_cache(ii,ii_write+i)
1756               uscdiff_cache(ii,i)=uscdiff_cache(ii,ii_write+i)
1757             enddo
1758
1759             do ii=1,nres*2
1760              do j=1,3
1761               c_cache(j,ii,i)=c_cache(j,ii,ii_write+i)
1762              enddo
1763             enddo
1764       enddo
1765       ntwx_cache=ntwx_cache-ii_write
1766       return
1767       end subroutine write1traj
1768 !-----------------------------------------------------------------------------
1769       subroutine read1restart(i_index)
1770
1771 !      implicit real*8 (a-h,o-z)
1772 !      include 'DIMENSIONS'
1773       include 'mpif.h'
1774 !      include 'COMMON.MD'
1775 !      include 'COMMON.IOUNITS'
1776 !      include 'COMMON.REMD'
1777 !      include 'COMMON.SETUP'
1778 !      include 'COMMON.CHAIN'
1779 !      include 'COMMON.SBRIDGE'
1780 !      include 'COMMON.INTERACT'
1781 !el      real(kind=4) :: d_restart1(3,2*nres*maxprocs)
1782       real(kind=4) :: r_d(3,0:2*nres),t5_restart1(5)
1783 !      integer(kind=2) :: i_index(Nprocs/4,Nprocs/20,Nprocs/200,Nprocs/200)
1784       integer(kind=2) :: i_index(Nprocs,Nprocs,Nprocs,Nprocs)
1785
1786       !(maxprocs/4,maxprocs/20,maxprocs/200,maxprocs/200)
1787 !el      common /przechowalnia/ d_restart1
1788       integer :: i,j,il,il1,ierr,itmp,iret,ixdrf
1789
1790       write (*,*) "Processor",me," called read1restart"
1791
1792          if(me.eq.king)then
1793               open(irest2,file=mremd_rst_name,status='unknown')
1794               read(irest2,*,err=334) i
1795               write(iout,*) "Reading old rst in ASCI format"
1796               close(irest2)
1797                call read1restart_old
1798                return
1799  334          continue
1800 #ifdef AIX
1801               call xdrfopen_(ixdrf,mremd_rst_name, "r", iret)
1802
1803               do i=0,nodes-1
1804                call xdrfint_(ixdrf, i2rep(i), iret)
1805               enddo
1806               do i=1,remd_m(1)
1807                call xdrfint_(ixdrf, ifirst(i), iret)
1808               enddo
1809              do il=1,nodes
1810               call xdrfint_(ixdrf, nupa(0,il), iret)
1811               do i=1,nupa(0,il)
1812                call xdrfint_(ixdrf, nupa(i,il), iret)
1813               enddo
1814
1815               call xdrfint_(ixdrf, ndowna(0,il), iret)
1816               do i=1,ndowna(0,il)
1817                call xdrfint_(ixdrf, ndowna(i,il), iret)
1818               enddo
1819              enddo
1820              do il=1,nodes
1821                do j=1,4
1822                 call xdrffloat_(ixdrf, t_restart1(j,il), iret)
1823                enddo
1824              enddo
1825 #else
1826               call xdrfopen(ixdrf,mremd_rst_name, "r", iret)
1827
1828               do i=0,nodes-1
1829                call xdrfint(ixdrf, i2rep(i), iret)
1830               enddo
1831               do i=1,remd_m(1)
1832                call xdrfint(ixdrf, ifirst(i), iret)
1833               enddo
1834              do il=1,nodes
1835               call xdrfint(ixdrf, nupa(0,il), iret)
1836               do i=1,nupa(0,il)
1837                call xdrfint(ixdrf, nupa(i,il), iret)
1838               enddo
1839
1840               call xdrfint(ixdrf, ndowna(0,il), iret)
1841               do i=1,ndowna(0,il)
1842                call xdrfint(ixdrf, ndowna(i,il), iret)
1843               enddo
1844              enddo
1845              do il=1,nodes
1846                do j=1,4
1847                 call xdrffloat(ixdrf, t_restart1(j,il), iret)
1848                enddo
1849              enddo
1850 #endif
1851          endif
1852          call mpi_scatter(t_restart1,5,mpi_real,&
1853                  t5_restart1,5,mpi_real,king,CG_COMM,ierr)
1854          totT=t5_restart1(1)              
1855          EK=t5_restart1(2)
1856          potE=t5_restart1(3)
1857          t_bath=t5_restart1(4)
1858
1859          if(me.eq.king)then
1860               do il=0,nodes-1
1861                do i=0,2*nres
1862 !                read(irest2,'(3e15.5)') 
1863 !     &                (d_restart1(j,i+2*nres*il),j=1,3)
1864             do j=1,3
1865 #ifdef AIX
1866              call xdrffloat_(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1867 #else
1868              call xdrffloat(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1869 #endif
1870             enddo
1871                enddo
1872               enddo
1873          endif
1874          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
1875                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
1876
1877          do i=0,2*nres
1878            do j=1,3
1879             d_t(j,i)=r_d(j,i)
1880            enddo
1881          enddo
1882          if(me.eq.king)then 
1883               do il=0,nodes-1
1884                do i=0,2*nres
1885 !                read(irest2,'(3e15.5)') 
1886 !     &                (d_restart1(j,i+2*nres*il),j=1,3)
1887             do j=1,3
1888 #ifdef AIX
1889              call xdrffloat_(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1890 #else
1891              call xdrffloat(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1892 #endif
1893             enddo
1894                enddo
1895               enddo
1896          endif
1897          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
1898                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
1899          do i=0,2*nres
1900            do j=1,3
1901             dc(j,i)=r_d(j,i)
1902            enddo
1903          enddo
1904        
1905
1906            if(usampl) then
1907 #ifdef AIX
1908              if(me.eq.king)then
1909               call xdrfint_(ixdrf, nset, iret)
1910               do i=1,nset
1911                 call xdrfint_(ixdrf,mset(i), iret)
1912               enddo
1913               do i=0,nodes-1
1914                 call xdrfint_(ixdrf,i2set(i), iret)
1915               enddo
1916               do il=1,nset
1917                do il1=1,mset(il)
1918                 do i=1,nrep
1919                  do j=1,remd_m(i)
1920                    call xdrfint_(ixdrf,itmp, iret)
1921                    i_index(i,j,il,il1)=itmp
1922                  enddo
1923                 enddo
1924                enddo
1925               enddo
1926              endif
1927 #else
1928              if(me.eq.king)then
1929               call xdrfint(ixdrf, nset, iret)
1930               do i=1,nset
1931                 call xdrfint(ixdrf,mset(i), iret)
1932               enddo
1933               do i=0,nodes-1
1934                 call xdrfint(ixdrf,i2set(i), iret)
1935               enddo
1936               do il=1,nset
1937                do il1=1,mset(il)
1938                 do i=1,nrep
1939                  do j=1,remd_m(i)
1940                    call xdrfint(ixdrf,itmp, iret)
1941                    i_index(i,j,il,il1)=itmp
1942                  enddo
1943                 enddo
1944                enddo
1945               enddo
1946              endif
1947 #endif
1948               call mpi_scatter(i2set,1,mpi_integer,&
1949                  iset,1,mpi_integer,king,&
1950                  CG_COMM,ierr) 
1951
1952            endif
1953
1954         if(me.eq.king) close(irest2)
1955       return
1956       end subroutine read1restart
1957 !-----------------------------------------------------------------------------
1958       subroutine read1restart_old
1959
1960 !      implicit real*8 (a-h,o-z)
1961 !      include 'DIMENSIONS'
1962       include 'mpif.h'
1963 !      include 'COMMON.MD'
1964 !      include 'COMMON.IOUNITS'
1965 !      include 'COMMON.REMD'
1966 !      include 'COMMON.SETUP'
1967 !      include 'COMMON.CHAIN'
1968 !      include 'COMMON.SBRIDGE'
1969 !      include 'COMMON.INTERACT'
1970 !el      real(kind=4) :: d_restart1(3,2*nres*maxprocs)
1971       real(kind=4) :: r_d(3,0:2*nres),t5_restart1(5)
1972 !el      common /przechowalnia/ d_restart1
1973
1974       integer :: i,j,il,ierr
1975
1976          if(me.eq.king)then
1977              open(irest2,file=mremd_rst_name,status='unknown')
1978              read (irest2,*) (i2rep(i),i=0,nodes-1)
1979              read (irest2,*) (ifirst(i),i=1,remd_m(1))
1980              do il=1,nodes
1981               read (irest2,*) nupa(0,il),(nupa(i,il),i=1,nupa(0,il))
1982               read (irest2,*) ndowna(0,il),&
1983                           (ndowna(i,il),i=1,ndowna(0,il))
1984              enddo
1985              do il=1,nodes
1986                read(irest2,*) (t_restart1(j,il),j=1,4)
1987              enddo
1988          endif
1989          call mpi_scatter(t_restart1,5,mpi_real,&
1990                  t5_restart1,5,mpi_real,king,CG_COMM,ierr)
1991          totT=t5_restart1(1)              
1992          EK=t5_restart1(2)
1993          potE=t5_restart1(3)
1994          t_bath=t5_restart1(4)
1995
1996          if(me.eq.king)then
1997               do il=0,nodes-1
1998                do i=0,2*nres
1999                 read(irest2,'(3e15.5)')  &
2000                       (d_restart1(j,i+(2*nres+1)*il),j=1,3)
2001                enddo
2002               enddo
2003          endif
2004          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
2005                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
2006
2007          do i=0,2*nres
2008            do j=1,3
2009             d_t(j,i)=r_d(j,i)
2010            enddo
2011          enddo
2012          if(me.eq.king)then 
2013               do il=0,nodes-1
2014                do i=1,2*nres
2015                 read(irest2,'(3e15.5)') &
2016                       (d_restart1(j,i+(2*nres+1)*il),j=1,3)
2017                enddo
2018               enddo
2019          endif
2020          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
2021                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
2022          do i=0,2*nres
2023            do j=1,3
2024             dc(j,i)=r_d(j,i)
2025            enddo
2026          enddo
2027         if(me.eq.king) close(irest2)
2028       return
2029       end subroutine read1restart_old
2030 !----------------------------------------------------------------
2031       subroutine alloc_MREMD_arrays
2032
2033 !      if(.not.allocated(mset)) allocate(mset(max0(nset,1)))
2034       if(.not.allocated(stdfsc)) allocate(stdfsc(ntyp1,5)) !(ntyp1))
2035 ! commom.remd
2036 !      common /remdcommon/ in io: read_REMDpar
2037 !      real(kind=8),dimension(:),allocatable :: remd_t !(maxprocs)
2038 !      integer,dimension(:),allocatable :: remd_m !(maxprocs)
2039 !      common /remdrestart/
2040       if(.not.allocated(i2rep)) allocate(i2rep(0:2*nodes))
2041
2042       allocate(i2set(0:2*nodes)) !(0:maxprocs)
2043       allocate(ifirst(0:nodes)) !(maxprocs)
2044       allocate(nupa(0:nodes,0:2*nodes))
2045       allocate(ndowna(0:nodes,0:2*nodes)) !(0:maxprocs/4,0:maxprocs)
2046       allocate(t_restart1(5,nodes)) !(5,maxprocs)
2047       allocate(iset_restart1(nodes)) !(maxprocs)
2048 !      common /traj1cache/
2049       allocate(totT_cache(max_cache_traj),EK_cache(max_cache_traj))
2050       allocate(potE_cache(max_cache_traj),t_bath_cache(max_cache_traj))
2051       allocate(Uconst_cache(max_cache_traj)) !(max_cache_traj)
2052       allocate(qfrag_cache(nfrag,max_cache_traj)) !(50,max_cache_traj)
2053       allocate(qpair_cache(npair,max_cache_traj)) !(100,max_cache_traj)
2054       allocate(ugamma_cache(nfrag_back,max_cache_traj))
2055       allocate(utheta_cache(nfrag_back,max_cache_traj))
2056       allocate(uscdiff_cache(nfrag_back,max_cache_traj)) !(maxfrag_back,max_cache_traj)
2057       allocate(c_cache(3,2*nres+2,max_cache_traj)) !(3,maxres2+2,max_cache_traj)
2058       allocate(iset_cache(max_cache_traj)) !(max_cache_traj)
2059
2060       return
2061       end subroutine alloc_MREMD_arrays
2062 !-----------------------------------------------------------------------------
2063 !-----------------------------------------------------------------------------
2064       end module MREMDyn