changes in wham and cluser + cutoff corr
[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             write(iout,*) icache_all
748             if(traj1file) then
749               write(iout,*) "before mpi_gather ntwx_cache"
750              call mpi_gather(ntwx_cache,1,mpi_integer,&
751                    icache_all(1),1,mpi_integer,king,&  ! CONSULT WITH ADAM
752                    CG_COMM,ierr)
753               write(iout,*) "after mpi_gather ntwx_cache"
754
755              if (me.eq.king) then
756                write(iout,'(a19,8000i8)') '     ntwx_cache',&
757                           (icache_all(i),i=1,nodes)
758                ii_write=icache_all(1)
759                do i=2,nodes
760                  if(icache_all(i).lt.ii_write) ii_write=icache_all(i)
761                enddo
762                write(iout,*) "MIN ii_write=",ii_write
763              endif
764             endif
765            call flush(iout)
766         endif
767         if (synflag) then
768 ! Update the time safety limiy
769           if (time001-time00.gt.safety) then
770             safety=time001-time00+600
771             write (iout,*) "****** SAFETY increased to",safety," s"
772           endif
773           if (ovrtim()) end_of_run=.true.
774         endif
775         if(synflag.and..not.end_of_run) then
776            time02=MPI_WTIME()
777            synflag=.false.
778
779            write(iout,*) 'REMD before',me,t_bath
780
781 !           call mpi_gather(t_bath,1,mpi_double_precision,
782 !     &             remd_t_bath,1,mpi_double_precision,king,
783 !     &             CG_COMM,ierr)
784            potEcomp(n_ene+1)=t_bath
785            if (usampl) then
786              potEcomp(n_ene+2)=iset
787              if (iset.lt.nset) then
788                i_set_temp=iset
789                iset=iset+1
790                call EconstrQ
791                potEcomp(n_ene+3)=Uconst
792                iset=i_set_temp
793              endif
794              if (iset.gt.1) then
795                i_set_temp=iset
796                iset=iset-1
797                call EconstrQ
798                potEcomp(n_ene+4)=Uconst 
799                iset=i_set_temp
800              endif
801            endif
802            call mpi_gather(potEcomp(0),n_ene+5,mpi_double_precision,&
803                    remd_ene(0,1),n_ene+5,mpi_double_precision,king,&
804                    CG_COMM,ierr)
805            if(lmuca) then 
806             call mpi_gather(elow,1,mpi_double_precision,&
807                    elowi,1,mpi_double_precision,king,&
808                    CG_COMM,ierr)
809             call mpi_gather(ehigh,1,mpi_double_precision,&
810                    ehighi,1,mpi_double_precision,king,&
811                    CG_COMM,ierr)
812            endif
813
814           time03=MPI_WTIME()
815           if (me.eq.king .or. .not. out1file) then
816             write(iout,*) 'REMD gather times=',time03-time01 &
817                                               ,time03-time02
818           endif
819
820           if (restart1file) call write1rst(i_index)
821
822           time04=MPI_WTIME()
823           if (me.eq.king .or. .not. out1file) then
824             write(iout,*) 'REMD writing rst time=',time04-time03
825           endif
826
827           if (traj1file) call write1traj
828 !d debugging
829 !deb            call mpi_gather(ntwx_cache,1,mpi_integer,
830 !deb     &             icache_all,1,mpi_integer,king,
831 !deb     &             CG_COMM,ierr)
832 !deb            write(iout,'(a19,8000i8)') '  ntwx_cache after traj1file',
833 !deb     &                    (icache_all(i),i=1,nodes)
834 !d end
835
836
837           time05=MPI_WTIME()
838           if (me.eq.king .or. .not. out1file) then
839             write(iout,*) 'REMD writing traj time=',time05-time04
840             call flush(iout)
841           endif
842
843
844           if (me.eq.king) then
845             do i=1,nodes
846                remd_t_bath(i)=remd_ene(n_ene+1,i)
847                iremd_iset(i)=remd_ene(n_ene+2,i)
848             enddo
849 #ifdef DEBUG
850             if(lmuca) then
851 !o             write(iout,*) 'REMD exchange temp,ene,elow,ehigh'
852              do i=1,nodes
853                write(iout,'(i4,4f12.5)') i,remd_t_bath(i),remd_ene(0,i),&
854                   elowi(i),ehighi(i)       
855              enddo
856             else
857               write(iout,*) 'REMD exchange temp,ene'
858               do i=1,nodes
859                 write(iout,'(i4,2f12.5)') i,remd_t_bath(i),remd_ene(0,i)
860                 write(iout,'(6f12.5)') (remd_ene(j,i),j=1,n_ene)
861               enddo
862             endif
863 #endif
864 !-------------------------------------           
865            IF(.not.usampl) THEN
866             write (iout,*) "Enter exchnge, remd_m",remd_m(1),&
867               " nodes",nodes
868             call flush(iout)
869             write (iout,*) "remd_m(1)",remd_m(1)
870             do irr=1,remd_m(1)
871                i=ifirst(iran_num(1,remd_m(1)))
872              write (iout,*) "i",i
873              call flush(iout)
874
875              do ii=1,nodes-1
876
877 #ifdef DEBUG
878               write (iout,*) "i",i," nupa(0,i)",int(nupa(0,i))
879 #endif
880              if(i.gt.0.and.nupa(0,i).gt.0) then
881               iex=i
882 !              if (i.eq.1 .and. int(nupa(0,i)).eq.1) then
883 !                write (iout,*) 
884 !     &  "CHUJ ABSOLUTNY!!! No way to sample a distinct replica in MREMD"
885 !                call flush(iout)
886 !                call MPI_Abort(MPI_COMM_WORLD,ERRCODE,ierr)
887 !              endif
888 !              do while (iex.eq.i)
889 !                write (iout,*) "upper",nupa(int(nupa(0,i)),i)
890                 iex=nupa(iran_num(1,int(nupa(0,i))),i)
891 !              enddo
892 !              write (iout,*) "nupa(0,i)",nupa(0,i)," iex",iex
893               if (lmuca) then
894                call muca_delta(remd_t_bath,remd_ene,i,iex,delta)
895               else
896 ! Swap temperatures between conformations i and iex with recalculating the free energies
897 ! following temperature changes.
898                ene_iex_iex=remd_ene(0,iex)
899                ene_i_i=remd_ene(0,i)
900 !               write (iout,*) "i",i," ene_i_i",ene_i_i,
901 !     &          " iex",iex," ene_iex_iex",ene_iex_iex
902 !               write (iout,*) "rescaling weights with temperature",
903 !     &          remd_t_bath(i)
904 !               call flush(iout)
905                call rescale_weights(remd_t_bath(i))
906
907 !               write (iout,*) "0,iex",remd_t_bath(i)
908 !               call enerprint(remd_ene(0,iex))
909
910                call sum_energy(remd_ene(0,iex),.false.)
911                ene_iex_i=remd_ene(0,iex)
912 !               write (iout,*) "ene_iex_i",remd_ene(0,iex)
913
914 !               write (iout,*) "0,i",remd_t_bath(i)
915 !               call enerprint(remd_ene(0,i))
916
917                call sum_energy(remd_ene(0,i),.false.)
918 !               write (iout,*) "ene_i_i",remd_ene(0,i)
919 !               call flush(iout)
920 !               write (iout,*) "rescaling weights with temperature",
921 !     &          remd_t_bath(iex)
922                if (real(ene_i_i).ne.real(remd_ene(0,i))) then
923                 write (iout,*) "ERROR: inconsistent energies:",i,&
924                   ene_i_i,remd_ene(0,i)
925                endif
926                call rescale_weights(remd_t_bath(iex))
927
928 !               write (iout,*) "0,i",remd_t_bath(iex)
929 !               call enerprint(remd_ene(0,i))
930
931                call sum_energy(remd_ene(0,i),.false.)
932 !               write (iout,*) "ene_i_iex",remd_ene(0,i)
933 !               call flush(iout)
934                ene_i_iex=remd_ene(0,i)
935
936 !               write (iout,*) "0,iex",remd_t_bath(iex)
937                call enerprint(remd_ene(0,iex))
938
939                call sum_energy(remd_ene(0,iex),.false.)
940                if (real(ene_iex_iex).ne.real(remd_ene(0,iex))) then
941                 write (iout,*) "ERROR: inconsistent energies:",iex,&
942                   ene_iex_iex,remd_ene(0,iex)
943                endif
944 !               write (iout,*) "ene_iex_iex",remd_ene(0,iex)
945 !               write (iout,*) "i",i," iex",iex
946 !               write (iout,'(4(a,e15.5))') "ene_i_i",ene_i_i,
947 !     &           " ene_i_iex",ene_i_iex,
948 !     &           " ene_iex_i",ene_iex_i," ene_iex_iex",ene_iex_iex
949 !               call flush(iout)
950                delta=(ene_iex_iex-ene_i_iex)/(Rb*remd_t_bath(iex))- &
951                     (ene_iex_i-ene_i_i)/(Rb*remd_t_bath(i))
952                delta=-delta
953 !               write(iout,*) 'delta',delta
954 !              delta=(remd_t_bath(i)-remd_t_bath(iex))*
955 !     &              (remd_ene(i)-remd_ene(iex))/Rb/
956 !     &              (remd_t_bath(i)*remd_t_bath(iex))
957               endif
958               if (delta .gt. 50.0d0) then
959                 delta=0.0d0
960               else
961 #ifdef OSF 
962                 if(isnan(delta))then
963                   delta=0.0d0
964                 else if (delta.lt.-50.0d0) then
965                   delta=dexp(50.0d0)
966                 else
967                   delta=dexp(-delta)
968                 endif
969 #else
970                 delta=dexp(-delta)
971 #endif
972               endif
973               iremd_tot(int(i2rep(i-1)))=iremd_tot(int(i2rep(i-1)))+1
974               xxx=ran_number(0.0d0,1.0d0)
975 !              write(iout,'(2i4,a6,2f12.5)') i,iex,' delta',delta,xxx
976 !              call flush(iout)
977               if (delta .gt. xxx) then
978                 tmp=remd_t_bath(i)       
979                 remd_t_bath(i)=remd_t_bath(iex)
980                 remd_t_bath(iex)=tmp
981                 remd_ene(0,i)=ene_i_iex
982                 remd_ene(0,iex)=ene_iex_i
983                 if(lmuca) then
984                   tmp=elowi(i)
985                   elowi(i)=elowi(iex)
986                   elowi(iex)=tmp  
987                   tmp=ehighi(i)
988                   ehighi(i)=ehighi(iex)
989                   ehighi(iex)=tmp  
990                 endif
991
992
993                 do k=0,nodes
994                   itmp=nupa(k,i)
995                   nupa(k,i)=nupa(k,iex)
996                   nupa(k,iex)=itmp
997                   itmp=ndowna(k,i)
998                   ndowna(k,i)=ndowna(k,iex)
999                   ndowna(k,iex)=itmp
1000                 enddo
1001                 do il=1,nodes
1002                  if (ifirst(il).eq.i) ifirst(il)=iex
1003                  do k=1,nupa(0,il)
1004                   if (nupa(k,il).eq.i) then 
1005                      nupa(k,il)=iex
1006                   elseif (nupa(k,il).eq.iex) then 
1007                      nupa(k,il)=i
1008                   endif
1009                  enddo
1010                  do k=1,ndowna(0,il)
1011                   if (ndowna(k,il).eq.i) then 
1012                      ndowna(k,il)=iex
1013                   elseif (ndowna(k,il).eq.iex) then 
1014                      ndowna(k,il)=i
1015                   endif
1016                  enddo
1017                 enddo
1018
1019                 iremd_acc(int(i2rep(i-1)))=iremd_acc(int(i2rep(i-1)))+1
1020                 itmp=i2rep(i-1)
1021                 i2rep(i-1)=i2rep(iex-1)
1022                 i2rep(iex-1)=itmp
1023
1024 !                write(iout,*) 'exchange',i,iex
1025 !                write (iout,'(a8,100i4)') "@ ifirst",
1026 !     &                    (ifirst(k),k=1,remd_m(1))
1027 !                do il=1,nodes
1028 !                 write (iout,'(a8,i4,a1,100i4)') "@ nupa",il,":",
1029 !     &                    (nupa(k,il),k=1,nupa(0,il))
1030 !                 write (iout,'(a8,i4,a1,100i4)') "@ ndowna",il,":",
1031 !     &                    (ndowna(k,il),k=1,ndowna(0,il))
1032 !                enddo
1033 !                call flush(iout) 
1034
1035               else
1036                remd_ene(0,iex)=ene_iex_iex
1037                remd_ene(0,i)=ene_i_i
1038                i=iex
1039               endif 
1040             endif
1041            enddo
1042            enddo
1043 !d           write (iout,*) "exchange completed"
1044 !d           call flush(iout) 
1045         ELSE
1046           do ii=1,nodes  
1047 !d            write(iout,*) "########",ii
1048
1049             i_temp=iran_num(1,nrep)
1050             i_mult=iran_num(1,remd_m(i_temp))
1051             i_iset=iran_num(1,nset)
1052             i_mset=iran_num(1,mset(i_iset))
1053             i=i_index(i_temp,i_mult,i_iset,i_mset)
1054
1055 !d            write(iout,*) "i=",i,i_temp,i_mult,i_iset,i_mset
1056
1057             i_dir=iran_num(1,3)
1058 !d            write(iout,*) "i_dir=",i_dir
1059
1060             if(i_dir.eq.1 .and. remd_m(i_temp+1).gt.0 )then            
1061                
1062                i_temp1=i_temp+1
1063                i_mult1=iran_num(1,remd_m(i_temp1))
1064                i_iset1=i_iset
1065                i_mset1=iran_num(1,mset(i_iset1))
1066                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1067
1068             elseif(i_dir.eq.2 .and. mset(i_iset+1).gt.0)then
1069
1070                i_temp1=i_temp
1071                i_mult1=iran_num(1,remd_m(i_temp1))
1072                i_iset1=i_iset+1
1073                i_mset1=iran_num(1,mset(i_iset1))
1074                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1075                econstr_temp_i=remd_ene(20,i)
1076                econstr_temp_iex=remd_ene(20,iex)
1077                remd_ene(20,i)=remd_ene(n_ene+3,i)
1078                remd_ene(20,iex)=remd_ene(n_ene+4,iex)
1079
1080             elseif(remd_m(i_temp+1).gt.0.and.mset(i_iset+1).gt.0)then
1081
1082                i_temp1=i_temp+1
1083                i_mult1=iran_num(1,remd_m(i_temp1))
1084                i_iset1=i_iset+1
1085                i_mset1=iran_num(1,mset(i_iset1))
1086                iex=i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1087                econstr_temp_i=remd_ene(20,i)
1088                econstr_temp_iex=remd_ene(20,iex)
1089                remd_ene(20,i)=remd_ene(n_ene+3,i)
1090                remd_ene(20,iex)=remd_ene(n_ene+4,iex)
1091
1092             else
1093                goto 444 
1094             endif
1095  
1096 !d            write(iout,*) "iex=",iex,i_temp1,i_mult1,i_iset1,i_mset1
1097             call flush(iout)
1098
1099 ! Swap temperatures between conformations i and iex with recalculating the free energies
1100 ! following temperature changes.
1101               ene_iex_iex=remd_ene(0,iex)
1102               ene_i_i=remd_ene(0,i)
1103 !o              write (iout,*) "rescaling weights with temperature",
1104 !o     &          remd_t_bath(i)
1105               call rescale_weights(remd_t_bath(i))
1106               
1107               call sum_energy(remd_ene(0,iex),.false.)
1108               ene_iex_i=remd_ene(0,iex)
1109 !d              write (iout,*) "ene_iex_i",remd_ene(0,iex)
1110 !              call sum_energy(remd_ene(0,i),.false.)
1111 !d              write (iout,*) "ene_i_i",remd_ene(0,i)
1112 !              write (iout,*) "rescaling weights with temperature",
1113 !     &          remd_t_bath(iex)
1114 !              if (real(ene_i_i).ne.real(remd_ene(0,i))) then
1115 !                write (iout,*) "ERROR: inconsistent energies:",i,
1116 !     &            ene_i_i,remd_ene(0,i)
1117 !              endif
1118               call rescale_weights(remd_t_bath(iex))
1119               call sum_energy(remd_ene(0,i),.false.)
1120 !d              write (iout,*) "ene_i_iex",remd_ene(0,i)
1121               ene_i_iex=remd_ene(0,i)
1122 !              call sum_energy(remd_ene(0,iex),.false.)
1123 !              if (real(ene_iex_iex).ne.real(remd_ene(0,iex))) then
1124 !                write (iout,*) "ERROR: inconsistent energies:",iex,
1125 !     &            ene_iex_iex,remd_ene(0,iex)
1126 !              endif
1127 !d              write (iout,*) "ene_iex_iex",remd_ene(0,iex)
1128 !              write (iout,*) "i",i," iex",iex
1129 !d              write (iout,'(4(a,e15.5))') "ene_i_i",ene_i_i,
1130 !d     &           " ene_i_iex",ene_i_iex,
1131 !d     &           " ene_iex_i",ene_iex_i," ene_iex_iex",ene_iex_iex
1132               delta=(ene_iex_iex-ene_i_iex)/(Rb*remd_t_bath(iex))- &
1133                     (ene_iex_i-ene_i_i)/(Rb*remd_t_bath(i))
1134               delta=-delta
1135 !d              write(iout,*) 'delta',delta
1136 !              delta=(remd_t_bath(i)-remd_t_bath(iex))*
1137 !     &              (remd_ene(i)-remd_ene(iex))/Rb/
1138 !     &              (remd_t_bath(i)*remd_t_bath(iex))
1139               if (delta .gt. 50.0d0) then
1140                 delta=0.0d0
1141               else
1142                 delta=dexp(-delta)
1143               endif
1144               if (i_dir.eq.1.or.i_dir.eq.3) &
1145                iremd_tot(int(i2rep(i-1)))=iremd_tot(int(i2rep(i-1)))+1
1146               if (i_dir.eq.2.or.i_dir.eq.3) &
1147                 iremd_tot_usa(int(i2set(i-1)))= &
1148                        iremd_tot_usa(int(i2set(i-1)))+1
1149               xxx=ran_number(0.0d0,1.0d0)
1150 !d              write(iout,'(2i4,a6,2f12.5)') i,iex,' delta',delta,xxx
1151               if (delta .gt. xxx) then
1152                 tmp=remd_t_bath(i)       
1153                 remd_t_bath(i)=remd_t_bath(iex)
1154                 remd_t_bath(iex)=tmp
1155
1156                 itmp=iremd_iset(i)       
1157                 iremd_iset(i)=iremd_iset(iex)
1158                 iremd_iset(iex)=itmp
1159
1160                 remd_ene(0,i)=ene_i_iex
1161                 remd_ene(0,iex)=ene_iex_i
1162
1163                 if (i_dir.eq.1.or.i_dir.eq.3) &
1164                  iremd_acc(int(i2rep(i-1)))=iremd_acc(int(i2rep(i-1)))+1
1165
1166                 itmp=i2rep(i-1)
1167                 i2rep(i-1)=i2rep(iex-1)
1168                 i2rep(iex-1)=itmp
1169
1170                 if (i_dir.eq.2.or.i_dir.eq.3) &
1171                  iremd_acc_usa(int(i2set(i-1)))= &
1172                        iremd_acc_usa(int(i2set(i-1)))+1
1173
1174                 itmp=i2set(i-1)
1175                 i2set(i-1)=i2set(iex-1)
1176                 i2set(iex-1)=itmp
1177         
1178                 itmp=i_index(i_temp,i_mult,i_iset,i_mset)
1179                 i_index(i_temp,i_mult,i_iset,i_mset)= &
1180                       i_index(i_temp1,i_mult1,i_iset1,i_mset1)
1181                 i_index(i_temp1,i_mult1,i_iset1,i_mset1)=itmp
1182
1183               else
1184                remd_ene(0,iex)=ene_iex_iex
1185                remd_ene(0,i)=ene_i_i
1186                remd_ene(20,iex)=econstr_temp_iex
1187                remd_ene(20,i)=econstr_temp_i
1188               endif
1189
1190 !d      do il=1,nset
1191 !d       do il1=1,mset(il)
1192 !d        do i=1,nrep
1193 !d         do j=1,remd_m(i)
1194 !d          write(iout,*) i,j,il,il1,i_index(i,j,il,il1)
1195 !d         enddo
1196 !d        enddo
1197 !d       enddo
1198 !d      enddo
1199
1200  444      continue           
1201
1202           enddo
1203
1204
1205         ENDIF
1206
1207 !-------------------------------------
1208              write (iout,*) "NREP",nrep
1209              do i=1,nrep
1210               if(iremd_tot(i).ne.0) &
1211                 write(iout,'(a3,i4,2f12.5,i5)') 'ACC',i,remd_t(i) &
1212                  ,iremd_acc(i)/(1.0*iremd_tot(i)),iremd_tot(i)
1213              enddo
1214
1215              if(usampl) then
1216               do i=1,nset
1217                if(iremd_tot_usa(i).ne.0) &
1218                  write(iout,'(a10,i4,f12.5,i8)') 'ACC_usampl',i,&
1219                iremd_acc_usa(i)/(1.0*iremd_tot_usa(i)),iremd_tot_usa(i)
1220               enddo
1221              endif
1222
1223              call flush(iout)
1224
1225 !d              write (iout,'(a6,100i4)') "ifirst",
1226 !d     &                    (ifirst(i),i=1,remd_m(1))
1227 !d              do il=1,nodes
1228 !d               write (iout,'(a5,i4,a1,100i4)') "nup",il,":",
1229 !d     &                    (nupa(i,il),i=1,nupa(0,il))
1230 !d               write (iout,'(a5,i4,a1,100i4)') "ndown",il,":",
1231 !d     &                    (ndowna(i,il),i=1,ndowna(0,il))
1232 !d              enddo
1233             endif
1234
1235          time06=MPI_WTIME()
1236 !d         write (iout,*) "Before scatter"
1237 !d         call flush(iout)
1238 #ifdef DEBUG
1239          if (me.eq.king) then
1240            write (iout,*) "t_bath before scatter",remd_t_bath
1241            call flush(iout)
1242          endif
1243 #endif
1244          call mpi_scatter(remd_t_bath,1,mpi_double_precision,&
1245                  t_bath,1,mpi_double_precision,king,&
1246                  CG_COMM,ierr) 
1247 !d         write (iout,*) "After scatter"
1248 !d         call flush(iout)
1249          if(usampl) &
1250           call mpi_scatter(iremd_iset,1,mpi_integer,&
1251                  iset,1,mpi_integer,king,&
1252                  CG_COMM,ierr) 
1253
1254          time07=MPI_WTIME()
1255           if (me.eq.king .or. .not. out1file) then
1256             write(iout,*) 'REMD scatter time=',time07-time06
1257           endif
1258
1259          if(lmuca) then
1260            call mpi_scatter(elowi,1,mpi_double_precision,&
1261                  elow,1,mpi_double_precision,king,&
1262                  CG_COMM,ierr) 
1263            call mpi_scatter(ehighi,1,mpi_double_precision,&
1264                  ehigh,1,mpi_double_precision,king,&
1265                  CG_COMM,ierr) 
1266          endif
1267          call rescale_weights(t_bath)
1268 !o         write (iout,*) "Processor",me,
1269 !o     &    " rescaling weights with temperature",t_bath
1270          do j=1,5
1271          stdfp(j)=dsqrt(2*Rb*t_bath/d_time)
1272          do i=1,ntyp
1273            stdfsc(i,j)=dsqrt(2*Rb*t_bath/d_time)
1274          enddo 
1275          enddo
1276
1277 !de         write(iout,*) 'REMD after',me,t_bath
1278            time08=MPI_WTIME()
1279            if (me.eq.king .or. .not. out1file) then
1280             write(iout,*) 'REMD exchange time 8-0=',time08-time00
1281             write(iout,*) 'REMD exchange time 8-7=',time08-time07
1282             write(iout,*) 'REMD exchange time 7-6=',time07-time06
1283             write(iout,*) 'REMD exchange time 6-5=',time06-time05
1284             write(iout,*) 'REMD exchange time 5-4=',time05-time04
1285             write(iout,*) 'REMD exchange time 4-3=',time04-time03
1286             write(iout,*) 'REMD exchange time 3-2=',time03-time02
1287             write(iout,*) 'REMD exchange time 2-1=',time02-time01
1288             write(iout,*) 'REMD exchange time 1-0=',time01-time00
1289             call flush(iout)
1290            endif
1291         endif
1292       enddo
1293
1294       if (restart1file) then 
1295           if (me.eq.king .or. .not. out1file) &
1296             write(iout,*) 'writing restart at the end of run'
1297            call write1rst(i_index)
1298       endif
1299
1300       if (traj1file) call write1traj
1301 !d debugging
1302 !deb            call mpi_gather(ntwx_cache,1,mpi_integer,
1303 !deb     &             icache_all,1,mpi_integer,king,
1304 !deb     &             CG_COMM,ierr)
1305 !deb            write(iout,'(a40,8000i8)') 
1306 !deb     &             '  ntwx_cache after traj1file at the end',
1307 !deb     &             (icache_all(i),i=1,nodes)
1308 !d end
1309
1310
1311 #ifdef MPI
1312       t_MD=MPI_Wtime()-tt0
1313 #else
1314       t_MD=tcpu()-tt0
1315 #endif
1316       if (me.eq.king .or. .not. out1file) then
1317        write (iout,'(//35(1h=),a10,35(1h=)/10(/a40,1pe15.5))') &
1318         '  Timing  ',&
1319        'MD calculations setup:',t_MDsetup,&
1320        'Energy & gradient evaluation:',t_enegrad,&
1321        'Stochastic MD setup:',t_langsetup,&
1322        'Stochastic MD step setup:',t_sdsetup,&
1323        'MD steps:',t_MD
1324        write (iout,'(/28(1h=),a25,27(1h=))') &
1325        '  End of MD calculation  '
1326       endif
1327 !el  common /przechowalnia/
1328 !      deallocate(d_restart1)
1329 !      deallocate(d_restart2)
1330 !      deallocate(p_c)
1331 !el--------------
1332       return
1333       end subroutine MREMD
1334 !-----------------------------------------------------------------------------
1335       subroutine write1rst(i_index)
1336
1337       use control_data
1338 !      implicit real*8 (a-h,o-z)
1339 !      include 'DIMENSIONS'
1340       include 'mpif.h'
1341 !      include 'COMMON.MD'
1342 !      include 'COMMON.IOUNITS'
1343 !      include 'COMMON.REMD'
1344 !      include 'COMMON.SETUP'
1345 !      include 'COMMON.CHAIN'
1346 !      include 'COMMON.SBRIDGE'
1347 !      include 'COMMON.INTERACT'
1348                
1349 !el      real(kind=4) :: d_restart1(3,2*nres*maxprocs),&
1350 !el           d_restart2(3,2*nres*maxprocs)
1351       real(kind=4) :: r_d(3,0:2*nres)
1352       real(kind=4) :: t5_restart1(5)
1353       integer :: iret,itmp
1354 !      integer(kind=2) :: i_index(Nprocs/4,Nprocs/20,Nprocs/200,Nprocs/200)
1355       integer(kind=2) :: i_index(Nprocs,Nprocs,Nprocs,Nprocs)
1356
1357       !(maxprocs/4,maxprocs/20,maxprocs/200,maxprocs/200)
1358 !el       common /przechowalnia/ d_restart1,d_restart2
1359       integer :: i,j,il,il1,ierr,ixdrf
1360
1361        t5_restart1(1)=totT
1362        t5_restart1(2)=EK
1363        t5_restart1(3)=potE
1364        t5_restart1(4)=t_bath
1365        t5_restart1(5)=Uconst
1366        
1367        call mpi_gather(t5_restart1,5,mpi_real,&
1368             t_restart1,5,mpi_real,king,CG_COMM,ierr)
1369
1370
1371        do i=0,2*nres
1372          do j=1,3
1373            r_d(j,i)=d_t(j,i)
1374          enddo
1375        enddo
1376        call mpi_gather(r_d,3*2*nres+3,mpi_real,&
1377                  d_restart1,3*2*nres+3,mpi_real,king,&
1378                  CG_COMM,ierr)
1379        do j=1,3
1380        dc(j,0)=c(j,1)
1381        enddo
1382
1383
1384        do i=0,2*nres
1385          do j=1,3
1386            r_d(j,i)=dc(j,i)
1387          enddo
1388        enddo
1389        call mpi_gather(r_d,3*2*nres+3,mpi_real,&
1390                  d_restart2,3*2*nres+3,mpi_real,king,&
1391                  CG_COMM,ierr)
1392
1393        if(me.eq.king) then
1394 #ifdef AIX
1395          call xdrfopen_(ixdrf,mremd_rst_name, "w", iret)
1396          do i=0,nodes-1
1397           call xdrfint_(ixdrf, i2rep(i), iret)
1398          enddo
1399          do i=1,remd_m(1)
1400           call xdrfint_(ixdrf, ifirst(i), iret)
1401          enddo
1402          do il=1,nodes
1403               do i=0,nupa(0,il)
1404                call xdrfint_(ixdrf, nupa(i,il), iret)
1405               enddo
1406
1407               do i=0,ndowna(0,il)
1408                call xdrfint_(ixdrf, ndowna(i,il), iret)
1409               enddo
1410          enddo
1411
1412          do il=1,nodes
1413            do j=1,4
1414             call xdrffloat_(ixdrf, t_restart1(j,il), iret)
1415            enddo
1416          enddo
1417
1418          do il=0,nodes-1
1419            do i=0,2*nres
1420             do j=1,3
1421              call xdrffloat_(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1422             enddo
1423            enddo
1424          enddo
1425          do il=0,nodes-1
1426            do i=0,2*nres
1427             do j=1,3
1428              call xdrffloat_(ixdrf, d_restart2(j,i+(2*nres+1)*il), iret)
1429             enddo
1430            enddo
1431          enddo
1432
1433          if(usampl) then
1434            call xdrfint_(ixdrf, nset, iret)
1435            do i=1,nset
1436              call xdrfint_(ixdrf,mset(i), iret)
1437            enddo
1438            do i=0,nodes-1
1439              call xdrfint_(ixdrf,i2set(i), iret)
1440            enddo
1441            do il=1,nset
1442              do il1=1,mset(il)
1443                do i=1,nrep
1444                  do j=1,remd_m(i)
1445                    itmp=i_index(i,j,il,il1)
1446                    call xdrfint_(ixdrf,itmp, iret)
1447                  enddo
1448                enddo
1449              enddo
1450            enddo
1451            
1452          endif
1453          call xdrfclose_(ixdrf, iret)
1454 #else
1455          call xdrfopen(ixdrf,mremd_rst_name, "w", iret)
1456          do i=0,nodes-1
1457           call xdrfint(ixdrf, i2rep(i), iret)
1458          enddo
1459          do i=1,remd_m(1)
1460           call xdrfint(ixdrf, ifirst(i), iret)
1461          enddo
1462          do il=1,nodes
1463               do i=0,nupa(0,il)
1464                call xdrfint(ixdrf, nupa(i,il), iret)
1465               enddo
1466
1467               do i=0,ndowna(0,il)
1468                call xdrfint(ixdrf, ndowna(i,il), iret)
1469               enddo
1470          enddo
1471
1472          do il=1,nodes
1473            do j=1,4
1474             call xdrffloat(ixdrf, t_restart1(j,il), iret)
1475            enddo
1476          enddo
1477
1478          do il=0,nodes-1
1479            do i=0,2*nres
1480             do j=1,3
1481              call xdrffloat(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1482             enddo
1483            enddo
1484          enddo
1485          do il=0,nodes-1
1486            do i=0,2*nres
1487             do j=1,3
1488              call xdrffloat(ixdrf, d_restart2(j,i+(2*nres+1)*il), iret)
1489             enddo
1490            enddo
1491          enddo
1492
1493
1494              if(usampl) then
1495               call xdrfint(ixdrf, nset, iret)
1496               do i=1,nset
1497                 call xdrfint(ixdrf,mset(i), iret)
1498               enddo
1499               do i=0,nodes-1
1500                 call xdrfint(ixdrf,i2set(i), iret)
1501               enddo
1502               do il=1,nset
1503                do il1=1,mset(il)
1504                 do i=1,nrep
1505                  do j=1,remd_m(i)
1506                    itmp=i_index(i,j,il,il1)
1507                    call xdrfint(ixdrf,itmp, iret)
1508                  enddo
1509                 enddo
1510                enddo
1511               enddo
1512            
1513              endif
1514          call xdrfclose(ixdrf, iret)
1515 #endif
1516        endif
1517       return
1518       end subroutine  write1rst
1519 !-----------------------------------------------------------------------------
1520       subroutine write1traj
1521
1522 !      implicit real*8 (a-h,o-z)
1523 !      include 'DIMENSIONS'
1524       include 'mpif.h'
1525 !      include 'COMMON.MD'
1526 !      include 'COMMON.IOUNITS'
1527 !      include 'COMMON.REMD'
1528 !      include 'COMMON.SETUP'
1529 !      include 'COMMON.CHAIN'
1530 !      include 'COMMON.SBRIDGE'
1531 !      include 'COMMON.INTERACT'
1532                
1533       real(kind=4) :: t5_restart1(5)
1534       integer :: iret,itmp
1535       real(kind=4) :: xcoord(3,2*nres+2),prec
1536       real(kind=4) :: r_qfrag(50),r_qpair(100)
1537       real(kind=4) :: r_utheta(50),r_ugamma(100),r_uscdiff(100)
1538       real(kind=4) :: p_qfrag(50*Nprocs),p_qpair(100*Nprocs) !(100*maxprocs)
1539       real(kind=4) :: p_utheta(50*Nprocs),p_ugamma(100*Nprocs),&
1540            p_uscdiff(100*Nprocs) !(100*maxprocs)
1541 !el      real(kind=4) :: p_c(3,(nres2+2)*maxprocs)
1542       real(kind=4) :: r_c(3,2*nres+2)
1543 !el      common /przechowalnia/ p_c
1544
1545       integer :: i,j,il,ierr,ii,ixdrf
1546
1547       call mpi_bcast(ii_write,1,mpi_integer,&
1548                  king,CG_COMM,ierr)
1549
1550 ! debugging
1551       print *,'traj1file',me,ii_write,ntwx_cache
1552 ! end debugging
1553
1554 #ifdef AIX
1555       if(me.eq.king) call xdrfopen_(ixdrf,cartname, "a", iret)
1556 #else
1557       if(me.eq.king) call xdrfopen(ixdrf,cartname, "a", iret)
1558 #endif
1559       do ii=1,ii_write
1560 !       write (iout,*) "before gather write1traj: from node",ii
1561 !       call flush(iout)
1562 !       write (iout,*) totT_cache(ii),EK_cache(ii),potE_cache(ii),t_bath_cache(ii),Uconst_cache(ii)
1563 !       call flush(iout)
1564        t5_restart1(1)=totT_cache(ii)
1565        t5_restart1(2)=EK_cache(ii)
1566        t5_restart1(3)=potE_cache(ii)
1567        t5_restart1(4)=t_bath_cache(ii)
1568        t5_restart1(5)=Uconst_cache(ii)
1569 !       write (iout,*) "before gather write1traj: from node",ii,t5_restart1(1),t5_restart1(3),t5_restart1(5),t5_restart1(4)
1570        call flush(iout)
1571        call mpi_gather(t5_restart1,5,mpi_real,&
1572             t_restart1,5,mpi_real,king,CG_COMM,ierr)
1573 !       do il=1,nodes
1574 !       write (iout,*) "after gather write1traj: from node",il,t_restart1(1,il),t_restart1(3,il),t_restart1(5,il),t_restart1(4,il)
1575 !       enddo
1576
1577        call mpi_gather(iset_cache(ii),1,mpi_integer,&
1578             iset_restart1,1,mpi_integer,king,CG_COMM,ierr)
1579
1580           do i=1,nfrag
1581            r_qfrag(i)=qfrag_cache(i,ii)
1582           enddo
1583           do i=1,npair
1584            r_qpair(i)=qpair_cache(i,ii)
1585           enddo
1586           do i=1,nfrag_back
1587            r_utheta(i)=utheta_cache(i,ii)
1588            r_ugamma(i)=ugamma_cache(i,ii)
1589            r_uscdiff(i)=uscdiff_cache(i,ii)
1590           enddo
1591
1592         call mpi_gather(r_qfrag,nfrag,mpi_real,&
1593                  p_qfrag,nfrag,mpi_real,king,&
1594                  CG_COMM,ierr)
1595          call mpi_gather(r_qpair,npair,mpi_real,&
1596                 p_qpair,npair,mpi_real,king,&
1597                  CG_COMM,ierr)
1598          call mpi_gather(r_utheta,nfrag_back,mpi_real,&
1599                 p_utheta,nfrag_back,mpi_real,king,&
1600                  CG_COMM,ierr)
1601         call mpi_gather(r_ugamma,nfrag_back,mpi_real,&
1602                  p_ugamma,nfrag_back,mpi_real,king,&
1603                  CG_COMM,ierr)
1604         call mpi_gather(r_uscdiff,nfrag_back,mpi_real,&
1605                  p_uscdiff,nfrag_back,mpi_real,king,&
1606                  CG_COMM,ierr)
1607
1608 #ifdef DEBUG
1609         write (iout,*) "p_qfrag"
1610         do i=1,nodes
1611           write (iout,*) i,(p_qfrag((i-1)*nfrag+j),j=1,nfrag)
1612         enddo
1613         write (iout,*) "p_qpair"
1614         do i=1,nodes
1615           write (iout,*) i,(p_qpair((i-1)*npair+j),j=1,npair)
1616         enddo
1617         call flush(iout)
1618 #endif
1619         do i=1,nres*2
1620          do j=1,3
1621           r_c(j,i)=c_cache(j,i,ii)
1622          enddo
1623         enddo
1624
1625         call mpi_gather(r_c,3*2*nres,mpi_real,&
1626                  p_c,3*2*nres,mpi_real,king,&
1627                  CG_COMM,ierr)
1628
1629        if(me.eq.king) then
1630 #ifdef AIX
1631          do il=1,nodes
1632           call xdrffloat_(ixdrf, real(t_restart1(1,il)), iret)
1633           call xdrffloat_(ixdrf, real(t_restart1(3,il)), iret)
1634           call xdrffloat_(ixdrf, real(t_restart1(5,il)), iret)
1635           call xdrffloat_(ixdrf, real(t_restart1(4,il)), iret)
1636           call xdrfint_(ixdrf, nss, iret) 
1637           do j=1,nss
1638            if (dyn_ss) then
1639             call xdrfint(ixdrf, idssb(j)+nres, iret)
1640             call xdrfint(ixdrf, jdssb(j)+nres, iret)
1641            else
1642             call xdrfint_(ixdrf, ihpb(j), iret)
1643             call xdrfint_(ixdrf, jhpb(j), iret)
1644            endif
1645           enddo
1646           call xdrfint_(ixdrf, nfrag+npair+3*nfrag_back, iret)
1647           call xdrfint_(ixdrf, iset_restart1(il), iret)
1648           do i=1,nfrag
1649            call xdrffloat_(ixdrf, p_qfrag(i+(il-1)*nfrag), iret)
1650           enddo
1651           do i=1,npair
1652            call xdrffloat_(ixdrf, p_qpair(i+(il-1)*npair), iret)
1653           enddo
1654           do i=1,nfrag_back
1655            call xdrffloat_(ixdrf, p_utheta(i+(il-1)*nfrag_back), iret)
1656            call xdrffloat_(ixdrf, p_ugamma(i+(il-1)*nfrag_back), iret)
1657            call xdrffloat_(ixdrf, p_uscdiff(i+(il-1)*nfrag_back), iret)
1658           enddo
1659           prec=10000.0
1660           do i=1,nres
1661            do j=1,3
1662             xcoord(j,i)=p_c(j,i+(il-1)*nres*2)
1663            enddo
1664           enddo
1665           do i=nnt,nct
1666            do j=1,3
1667             xcoord(j,nres+i-nnt+1)=p_c(j,i+nres+(il-1)*nres*2)
1668            enddo
1669           enddo
1670           itmp=nres+nct-nnt+1
1671           call xdrf3dfcoord_(ixdrf, xcoord, itmp, prec, iret)
1672          enddo
1673 #else
1674          do il=1,nodes
1675           call xdrffloat(ixdrf, real(t_restart1(1,il)), iret)
1676           call xdrffloat(ixdrf, real(t_restart1(3,il)), iret)
1677           call xdrffloat(ixdrf, real(t_restart1(5,il)), iret)
1678           call xdrffloat(ixdrf, real(t_restart1(4,il)), iret)
1679 !          write (iout,*) "write1traj: from node",ii,t_restart1(1,il),t_restart1(3,il),t_restart1(5,il),t_restart1(4,il)
1680           call xdrfint(ixdrf, nss, iret) 
1681           do j=1,nss
1682            if (dyn_ss) then
1683             call xdrfint(ixdrf, idssb(j)+nres, iret)
1684             call xdrfint(ixdrf, jdssb(j)+nres, iret)
1685            else
1686             call xdrfint(ixdrf, ihpb(j), iret)
1687             call xdrfint(ixdrf, jhpb(j), iret)
1688            endif
1689           enddo
1690           call xdrfint(ixdrf, nfrag+npair+3*nfrag_back, iret)
1691           call xdrfint(ixdrf, iset_restart1(il), iret)
1692           do i=1,nfrag
1693            call xdrffloat(ixdrf, p_qfrag(i+(il-1)*nfrag), iret)
1694           enddo
1695           do i=1,npair
1696            call xdrffloat(ixdrf, p_qpair(i+(il-1)*npair), iret)
1697           enddo
1698           do i=1,nfrag_back
1699            call xdrffloat(ixdrf, p_utheta(i+(il-1)*nfrag_back), iret)
1700            call xdrffloat(ixdrf, p_ugamma(i+(il-1)*nfrag_back), iret)
1701            call xdrffloat(ixdrf, p_uscdiff(i+(il-1)*nfrag_back), iret)
1702           enddo
1703           prec=10000.0
1704           do i=1,nres
1705            do j=1,3
1706             xcoord(j,i)=p_c(j,i+(il-1)*nres*2)
1707            enddo
1708           enddo
1709           do i=nnt,nct
1710            do j=1,3
1711             xcoord(j,nres+i-nnt+1)=p_c(j,i+nres+(il-1)*nres*2)
1712            enddo
1713           enddo
1714           itmp=nres+nct-nnt+1
1715           call xdrf3dfcoord(ixdrf, xcoord, itmp, prec, iret)
1716          enddo
1717 #endif
1718        endif
1719       enddo
1720 #ifdef AIX
1721       if(me.eq.king) call xdrfclose_(ixdrf, iret)
1722 #else
1723       if(me.eq.king) call xdrfclose(ixdrf, iret)
1724 #endif
1725       do i=1,ntwx_cache-ii_write
1726
1727             totT_cache(i)=totT_cache(ii_write+i)
1728             EK_cache(i)=EK_cache(ii_write+i)
1729             potE_cache(i)=potE_cache(ii_write+i)
1730             t_bath_cache(i)=t_bath_cache(ii_write+i)
1731             Uconst_cache(i)=Uconst_cache(ii_write+i)
1732             iset_cache(i)=iset_cache(ii_write+i)
1733
1734             do ii=1,nfrag
1735              qfrag_cache(ii,i)=qfrag_cache(ii,ii_write+i)
1736             enddo
1737             do ii=1,npair
1738              qpair_cache(ii,i)=qpair_cache(ii,ii_write+i)
1739             enddo
1740             do ii=1,nfrag_back
1741               utheta_cache(ii,i)=utheta_cache(ii,ii_write+i)
1742               ugamma_cache(ii,i)=ugamma_cache(ii,ii_write+i)
1743               uscdiff_cache(ii,i)=uscdiff_cache(ii,ii_write+i)
1744             enddo
1745
1746             do ii=1,nres*2
1747              do j=1,3
1748               c_cache(j,ii,i)=c_cache(j,ii,ii_write+i)
1749              enddo
1750             enddo
1751       enddo
1752       ntwx_cache=ntwx_cache-ii_write
1753       return
1754       end subroutine write1traj
1755 !-----------------------------------------------------------------------------
1756       subroutine read1restart(i_index)
1757
1758 !      implicit real*8 (a-h,o-z)
1759 !      include 'DIMENSIONS'
1760       include 'mpif.h'
1761 !      include 'COMMON.MD'
1762 !      include 'COMMON.IOUNITS'
1763 !      include 'COMMON.REMD'
1764 !      include 'COMMON.SETUP'
1765 !      include 'COMMON.CHAIN'
1766 !      include 'COMMON.SBRIDGE'
1767 !      include 'COMMON.INTERACT'
1768 !el      real(kind=4) :: d_restart1(3,2*nres*maxprocs)
1769       real(kind=4) :: r_d(3,0:2*nres),t5_restart1(5)
1770 !      integer(kind=2) :: i_index(Nprocs/4,Nprocs/20,Nprocs/200,Nprocs/200)
1771       integer(kind=2) :: i_index(Nprocs,Nprocs,Nprocs,Nprocs)
1772
1773       !(maxprocs/4,maxprocs/20,maxprocs/200,maxprocs/200)
1774 !el      common /przechowalnia/ d_restart1
1775       integer :: i,j,il,il1,ierr,itmp,iret,ixdrf
1776
1777       write (*,*) "Processor",me," called read1restart"
1778
1779          if(me.eq.king)then
1780               open(irest2,file=mremd_rst_name,status='unknown')
1781               read(irest2,*,err=334) i
1782               write(iout,*) "Reading old rst in ASCI format"
1783               close(irest2)
1784                call read1restart_old
1785                return
1786  334          continue
1787 #ifdef AIX
1788               call xdrfopen_(ixdrf,mremd_rst_name, "r", iret)
1789
1790               do i=0,nodes-1
1791                call xdrfint_(ixdrf, i2rep(i), iret)
1792               enddo
1793               do i=1,remd_m(1)
1794                call xdrfint_(ixdrf, ifirst(i), iret)
1795               enddo
1796              do il=1,nodes
1797               call xdrfint_(ixdrf, nupa(0,il), iret)
1798               do i=1,nupa(0,il)
1799                call xdrfint_(ixdrf, nupa(i,il), iret)
1800               enddo
1801
1802               call xdrfint_(ixdrf, ndowna(0,il), iret)
1803               do i=1,ndowna(0,il)
1804                call xdrfint_(ixdrf, ndowna(i,il), iret)
1805               enddo
1806              enddo
1807              do il=1,nodes
1808                do j=1,4
1809                 call xdrffloat_(ixdrf, t_restart1(j,il), iret)
1810                enddo
1811              enddo
1812 #else
1813               call xdrfopen(ixdrf,mremd_rst_name, "r", iret)
1814
1815               do i=0,nodes-1
1816                call xdrfint(ixdrf, i2rep(i), iret)
1817               enddo
1818               do i=1,remd_m(1)
1819                call xdrfint(ixdrf, ifirst(i), iret)
1820               enddo
1821              do il=1,nodes
1822               call xdrfint(ixdrf, nupa(0,il), iret)
1823               do i=1,nupa(0,il)
1824                call xdrfint(ixdrf, nupa(i,il), iret)
1825               enddo
1826
1827               call xdrfint(ixdrf, ndowna(0,il), iret)
1828               do i=1,ndowna(0,il)
1829                call xdrfint(ixdrf, ndowna(i,il), iret)
1830               enddo
1831              enddo
1832              do il=1,nodes
1833                do j=1,4
1834                 call xdrffloat(ixdrf, t_restart1(j,il), iret)
1835                enddo
1836              enddo
1837 #endif
1838          endif
1839          call mpi_scatter(t_restart1,5,mpi_real,&
1840                  t5_restart1,5,mpi_real,king,CG_COMM,ierr)
1841          totT=t5_restart1(1)              
1842          EK=t5_restart1(2)
1843          potE=t5_restart1(3)
1844          t_bath=t5_restart1(4)
1845
1846          if(me.eq.king)then
1847               do il=0,nodes-1
1848                do i=0,2*nres
1849 !                read(irest2,'(3e15.5)') 
1850 !     &                (d_restart1(j,i+2*nres*il),j=1,3)
1851             do j=1,3
1852 #ifdef AIX
1853              call xdrffloat_(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1854 #else
1855              call xdrffloat(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1856 #endif
1857             enddo
1858                enddo
1859               enddo
1860          endif
1861          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
1862                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
1863
1864          do i=0,2*nres
1865            do j=1,3
1866             d_t(j,i)=r_d(j,i)
1867            enddo
1868          enddo
1869          if(me.eq.king)then 
1870               do il=0,nodes-1
1871                do i=0,2*nres
1872 !                read(irest2,'(3e15.5)') 
1873 !     &                (d_restart1(j,i+2*nres*il),j=1,3)
1874             do j=1,3
1875 #ifdef AIX
1876              call xdrffloat_(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1877 #else
1878              call xdrffloat(ixdrf, d_restart1(j,i+(2*nres+1)*il), iret)
1879 #endif
1880             enddo
1881                enddo
1882               enddo
1883          endif
1884          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
1885                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
1886          do i=0,2*nres
1887            do j=1,3
1888             dc(j,i)=r_d(j,i)
1889            enddo
1890          enddo
1891        
1892
1893            if(usampl) then
1894 #ifdef AIX
1895              if(me.eq.king)then
1896               call xdrfint_(ixdrf, nset, iret)
1897               do i=1,nset
1898                 call xdrfint_(ixdrf,mset(i), iret)
1899               enddo
1900               do i=0,nodes-1
1901                 call xdrfint_(ixdrf,i2set(i), iret)
1902               enddo
1903               do il=1,nset
1904                do il1=1,mset(il)
1905                 do i=1,nrep
1906                  do j=1,remd_m(i)
1907                    call xdrfint_(ixdrf,itmp, iret)
1908                    i_index(i,j,il,il1)=itmp
1909                  enddo
1910                 enddo
1911                enddo
1912               enddo
1913              endif
1914 #else
1915              if(me.eq.king)then
1916               call xdrfint(ixdrf, nset, iret)
1917               do i=1,nset
1918                 call xdrfint(ixdrf,mset(i), iret)
1919               enddo
1920               do i=0,nodes-1
1921                 call xdrfint(ixdrf,i2set(i), iret)
1922               enddo
1923               do il=1,nset
1924                do il1=1,mset(il)
1925                 do i=1,nrep
1926                  do j=1,remd_m(i)
1927                    call xdrfint(ixdrf,itmp, iret)
1928                    i_index(i,j,il,il1)=itmp
1929                  enddo
1930                 enddo
1931                enddo
1932               enddo
1933              endif
1934 #endif
1935               call mpi_scatter(i2set,1,mpi_integer,&
1936                  iset,1,mpi_integer,king,&
1937                  CG_COMM,ierr) 
1938
1939            endif
1940
1941         if(me.eq.king) close(irest2)
1942       return
1943       end subroutine read1restart
1944 !-----------------------------------------------------------------------------
1945       subroutine read1restart_old
1946
1947 !      implicit real*8 (a-h,o-z)
1948 !      include 'DIMENSIONS'
1949       include 'mpif.h'
1950 !      include 'COMMON.MD'
1951 !      include 'COMMON.IOUNITS'
1952 !      include 'COMMON.REMD'
1953 !      include 'COMMON.SETUP'
1954 !      include 'COMMON.CHAIN'
1955 !      include 'COMMON.SBRIDGE'
1956 !      include 'COMMON.INTERACT'
1957 !el      real(kind=4) :: d_restart1(3,2*nres*maxprocs)
1958       real(kind=4) :: r_d(3,0:2*nres),t5_restart1(5)
1959 !el      common /przechowalnia/ d_restart1
1960
1961       integer :: i,j,il,ierr
1962
1963          if(me.eq.king)then
1964              open(irest2,file=mremd_rst_name,status='unknown')
1965              read (irest2,*) (i2rep(i),i=0,nodes-1)
1966              read (irest2,*) (ifirst(i),i=1,remd_m(1))
1967              do il=1,nodes
1968               read (irest2,*) nupa(0,il),(nupa(i,il),i=1,nupa(0,il))
1969               read (irest2,*) ndowna(0,il),&
1970                           (ndowna(i,il),i=1,ndowna(0,il))
1971              enddo
1972              do il=1,nodes
1973                read(irest2,*) (t_restart1(j,il),j=1,4)
1974              enddo
1975          endif
1976          call mpi_scatter(t_restart1,5,mpi_real,&
1977                  t5_restart1,5,mpi_real,king,CG_COMM,ierr)
1978          totT=t5_restart1(1)              
1979          EK=t5_restart1(2)
1980          potE=t5_restart1(3)
1981          t_bath=t5_restart1(4)
1982
1983          if(me.eq.king)then
1984               do il=0,nodes-1
1985                do i=0,2*nres
1986                 read(irest2,'(3e15.5)')  &
1987                       (d_restart1(j,i+(2*nres+1)*il),j=1,3)
1988                enddo
1989               enddo
1990          endif
1991          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
1992                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
1993
1994          do i=0,2*nres
1995            do j=1,3
1996             d_t(j,i)=r_d(j,i)
1997            enddo
1998          enddo
1999          if(me.eq.king)then 
2000               do il=0,nodes-1
2001                do i=1,2*nres
2002                 read(irest2,'(3e15.5)') &
2003                       (d_restart1(j,i+(2*nres+1)*il),j=1,3)
2004                enddo
2005               enddo
2006          endif
2007          call mpi_scatter(d_restart1,3*2*nres+3,mpi_real,&
2008                  r_d,3*2*nres+3,mpi_real,king,CG_COMM,ierr)
2009          do i=0,2*nres
2010            do j=1,3
2011             dc(j,i)=r_d(j,i)
2012            enddo
2013          enddo
2014         if(me.eq.king) close(irest2)
2015       return
2016       end subroutine read1restart_old
2017 !----------------------------------------------------------------
2018       subroutine alloc_MREMD_arrays
2019
2020 !      if(.not.allocated(mset)) allocate(mset(max0(nset,1)))
2021       if(.not.allocated(stdfsc)) allocate(stdfsc(ntyp1,5)) !(ntyp1))
2022 ! commom.remd
2023 !      common /remdcommon/ in io: read_REMDpar
2024 !      real(kind=8),dimension(:),allocatable :: remd_t !(maxprocs)
2025 !      integer,dimension(:),allocatable :: remd_m !(maxprocs)
2026 !      common /remdrestart/
2027       if(.not.allocated(i2rep)) allocate(i2rep(0:2*nodes))
2028
2029       allocate(i2set(0:2*nodes)) !(0:maxprocs)
2030       allocate(ifirst(0:nodes)) !(maxprocs)
2031       allocate(nupa(0:nodes,0:2*nodes))
2032       allocate(ndowna(0:nodes,0:2*nodes)) !(0:maxprocs/4,0:maxprocs)
2033       allocate(t_restart1(5,nodes)) !(5,maxprocs)
2034       allocate(iset_restart1(nodes)) !(maxprocs)
2035 !      common /traj1cache/
2036       allocate(totT_cache(max_cache_traj),EK_cache(max_cache_traj))
2037       allocate(potE_cache(max_cache_traj),t_bath_cache(max_cache_traj))
2038       allocate(Uconst_cache(max_cache_traj)) !(max_cache_traj)
2039       allocate(qfrag_cache(nfrag,max_cache_traj)) !(50,max_cache_traj)
2040       allocate(qpair_cache(npair,max_cache_traj)) !(100,max_cache_traj)
2041       allocate(ugamma_cache(nfrag_back,max_cache_traj))
2042       allocate(utheta_cache(nfrag_back,max_cache_traj))
2043       allocate(uscdiff_cache(nfrag_back,max_cache_traj)) !(maxfrag_back,max_cache_traj)
2044       allocate(c_cache(3,2*nres+2,max_cache_traj)) !(3,maxres2+2,max_cache_traj)
2045       allocate(iset_cache(max_cache_traj)) !(max_cache_traj)
2046
2047       return
2048       end subroutine alloc_MREMD_arrays
2049 !-----------------------------------------------------------------------------
2050 !-----------------------------------------------------------------------------
2051       end module MREMDyn