remd with homol_nset>1 working with FG parallelization
authorCezary Czaplewski <czarek@chem.univ.gda.pl>
Tue, 17 Mar 2015 16:04:38 +0000 (17:04 +0100)
committerCezary Czaplewski <czarek@chem.univ.gda.pl>
Tue, 17 Mar 2015 16:04:38 +0000 (17:04 +0100)
source/unres/src_MD/MREMD.F
source/unres/src_MD/cinfo.f
source/unres/src_MD/minimize_p.F
source/unres/src_MD/readrtns.F

index ed655a8..bb1353e 100644 (file)
@@ -1,3 +1,4 @@
+#define DEBUG
 #ifdef MPI
       subroutine MREMD
       implicit real*8 (a-h,o-z)
@@ -29,7 +30,7 @@
       integer iremd_iset(maxprocs)
       integer*2 i_index
      &            (maxprocs/4,maxprocs/20,maxprocs/200,maxprocs/200)
-      double precision remd_ene(0:n_ene+4,maxprocs),t_bath_old
+      double precision remd_ene(0:n_ene+4,maxprocs),t_bath_old,e_tmp
       integer iremd_acc(maxprocs),iremd_tot(maxprocs)
       integer iremd_acc_usa(maxprocs),iremd_tot_usa(maxprocs)
       integer ilen,rstcount
@@ -290,6 +291,11 @@ cd       print *,'ttt',me,remd_tlist,(remd_t(i),i=1,nrep)
        endif
        if(usampl.or.hremd.gt.0.or.homol_nset.gt.1) then
           iset=i2set(me)
+c broadcast iset to slaves
+          if (nfgtasks.gt.1) then         
+           call MPI_Bcast(11,1,MPI_INTEGER,king,FG_COMM,IERROR)
+           call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+          endif
           if (hremd.gt.0) call set_hweights(iset)
           if(me.eq.king.or..not.out1file) 
      &     write(iout,*) me,"iset=",iset,"t_bath=",t_bath
@@ -734,25 +740,61 @@ c             write(iout,*) potEcomp(n_ene+1),potEcomp(n_ene+2),iset,nset
                i_set_temp=iset
                iset=iset+1
                if (homol_nset.gt.1) then
-                call e_modeller(potEcomp(n_ene+3))
+c broadcast iset to slaves and reduce energy
+                if (nfgtasks.gt.1) then         
+                 call MPI_Bcast(12,1,MPI_INTEGER,king,FG_COMM,IERROR)
+                 call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+                 call e_modeller(e_tmp)
+c                 write(iout,*) "iset+1 before reduce",e_tmp
+                 call MPI_Barrier(FG_COMM,IERR)
+                 call MPI_Reduce(e_tmp,potEcomp(n_ene+3),1,
+     &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
+                else
+                 call e_modeller(potEcomp(n_ene+3))
+                endif
 c                write(iout,*) "iset+1",potEcomp(n_ene+3)
                else
                 call EconstrQ
                 potEcomp(n_ene+3)=Uconst
                endif
                iset=i_set_temp
+c broadcast iset to slaves 
+               if (nfgtasks.gt.1) then
+                 call MPI_Bcast(11,1,MPI_INTEGER,king,FG_COMM,IERROR)
+                 call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+               endif
+             else
+              potEcomp(n_ene+3)=0.0
              endif
              if (iset.gt.1) then
                i_set_temp=iset
                iset=iset-1
                if (homol_nset.gt.1) then
-                call e_modeller(potEcomp(n_ene+4))
+c broadcast iset to slaves and reduce energy
+                if (nfgtasks.gt.1) then
+                 call MPI_Bcast(12,1,MPI_INTEGER,king,FG_COMM,IERROR)
+                 call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+                 call e_modeller(e_tmp)
+c                 write(iout,*) "iset-1 before reduce",e_tmp
+                 call MPI_Barrier(FG_COMM,IERR)
+                 call MPI_Reduce(e_tmp,potEcomp(n_ene+4),1,
+     &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)  
+                else
+                 call e_modeller(potEcomp(n_ene+4))
+                endif
 c                write(iout,*) "iset-1",potEcomp(n_ene+4)
                else
                 call EconstrQ
                 potEcomp(n_ene+4)=Uconst
                endif
                iset=i_set_temp
+c broadcast iset to slaves 
+               if (nfgtasks.gt.1) then
+                 call MPI_Bcast(11,1,MPI_INTEGER,king,FG_COMM,IERROR)
+                 call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+               endif
+             else
+               potEcomp(n_ene+4)=0.0
              endif
            endif
            if(hremd.gt.0) potEcomp(n_ene+2)=iset   
@@ -2033,11 +2075,7 @@ c     &                (d_restart1(j,i+2*nres*il),j=1,3)
               enddo
              endif
 #endif
-c<<<<<<< HEAD
 c Corrected AL 8/19/2014: each processor needs whole iset array not only its
-cv=======
-Corrected AL 8/19/2014: each processor needs whole iset array not only its
-c>>>>>>> 66693b0684c228404e7aadcffe6d2a8c9f489063
 c own element
 c              call mpi_scatter(i2set,1,mpi_integer,
 c     &           iset,1,mpi_integer,king,
@@ -2045,7 +2083,11 @@ c     &           CG_COMM,ierr)
               call mpi_bcast(i2set(0),nodes,mpi_integer,king,
      &         CG_COMM,ierr)
               iset=i2set(me)
-
+c broadcast iset to slaves
+              if (nfgtasks.gt.1) then         
+               call MPI_Bcast(11,1,MPI_INTEGER,king,FG_COMM,IERROR)
+               call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+              endif
            endif
 
 
index a688eb6..7b76974 100644 (file)
@@ -1,11 +1,11 @@
 C DO NOT EDIT THIS FILE - IT HAS BEEN GENERATED BY COMPINFO.C
-C 0 40376 53
+C 0 40376 64
       subroutine cinfo
       include 'COMMON.IOUNITS'
       write(iout,*)'++++ Compile info ++++'
-      write(iout,*)'Version 0.40376 build 53'
-      write(iout,*)'compiled Fri Mar  6 10:04:56 2015'
-      write(iout,*)'compiled by felipe@piasek4'
+      write(iout,*)'Version 0.40376 build 64'
+      write(iout,*)'compiled Tue Mar 17 12:39:20 2015'
+      write(iout,*)'compiled by czarek@piasek4'
       write(iout,*)'OS name:    Linux '
       write(iout,*)'OS release: 3.2.0-70-generic '
       write(iout,*)'OS version:',
index 9192bc1..14c8bd8 100644 (file)
@@ -134,6 +134,7 @@ c----------------------------------------------------------------------------
       double precision z(maxres6),d_a_tmp(maxres6)
       double precision edum(0:n_ene),time_order(0:10)
       double precision Gcopy(maxres2,maxres2)
+      double precision e_tmp
       common /przechowalnia/ Gcopy
       integer icall /0/
 C Workers wait for variables and NF, and NFL from the boss 
@@ -229,12 +230,24 @@ c           call MATOUT2(my_ng_count,dimen3,maxres2,maxers2,ginv_block)
           call fricmat_mult(z,d_a_tmp)
         else if (iorder.eq.10) then
           call setup_fricmat
+        else if (iorder.eq.11) then
+c broadcast new value of iset
+          call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+        else if (iorder.eq.12) then
+c broadcast new value of iset and
+c calculate restraint homology energy and
+c sum it over FG_COMM
+          call MPI_Bcast(iset,1,MPI_INTEGER,king,FG_COMM,IERROR)
+          call e_modeller(e_tmp)
+          call MPI_Barrier(FG_COMM,IERR)
+          call MPI_Reduce(e_tmp,ehomology_constr,1,
+     &    MPI_DOUBLE_PRECISION,MPI_SUM,king,FG_COMM,IERR)
         endif
       enddo
       write (*,*) 'Processor',fg_rank,' CG group',kolor,
      &  ' absolute rank',myrank,' leves ERGASTULUM.'
       write(*,*)'Processor',fg_rank,' wait times for respective orders',
-     & (' order[',i,']',time_order(i),i=0,10)
+     & (' order[',i,']',time_order(i),i=0,12)
       return
       end
 #endif
index b731cf2..69c7fe8 100644 (file)
@@ -2508,33 +2508,36 @@ c-------------------------------------------------------------------------------
       include 'COMMON.IOUNITS'
       include 'COMMON.MD'
       include 'COMMON.CONTROL'
+      integer iset1
       read(inp,*) nset,nfrag,npair,nfrag_back
       if(me.eq.king.or..not.out1file)
      & write(iout,*) "nset",nset," nfrag",nfrag," npair",npair,
      &  " nfrag_back",nfrag_back
-      do iset=1,nset
-         read(inp,*) mset(iset)
+      do iset1=1,nset
+         read(inp,*) mset(iset1)
        do i=1,nfrag
-         read(inp,*) wfrag(i,iset),ifrag(1,i,iset),ifrag(2,i,iset), 
-     &     qinfrag(i,iset)
+         read(inp,*) wfrag(i,iset1),ifrag(1,i,iset1),ifrag(2,i,iset1), 
+     &     qinfrag(i,iset1)
          if(me.eq.king.or..not.out1file)
-     &    write(iout,*) "R ",i,wfrag(i,iset),ifrag(1,i,iset),
-     &     ifrag(2,i,iset), qinfrag(i,iset)
+     &    write(iout,*) "R ",i,wfrag(i,iset1),ifrag(1,i,iset1),
+     &     ifrag(2,i,iset1), qinfrag(i,iset1)
        enddo
        do i=1,npair
-        read(inp,*) wpair(i,iset),ipair(1,i,iset),ipair(2,i,iset), 
-     &    qinpair(i,iset)
+        read(inp,*) wpair(i,iset1),ipair(1,i,iset1),ipair(2,i,iset1), 
+     &    qinpair(i,iset1)
         if(me.eq.king.or..not.out1file)
-     &   write(iout,*) "R ",i,wpair(i,iset),ipair(1,i,iset),
-     &    ipair(2,i,iset), qinpair(i,iset)
+     &   write(iout,*) "R ",i,wpair(i,iset1),ipair(1,i,iset1),
+     &    ipair(2,i,iset1), qinpair(i,iset1)
        enddo 
        do i=1,nfrag_back
-        read(inp,*) wfrag_back(1,i,iset),wfrag_back(2,i,iset),
-     &     wfrag_back(3,i,iset),
-     &     ifrag_back(1,i,iset),ifrag_back(2,i,iset)
+        read(inp,*) wfrag_back(1,i,iset1),wfrag_back(2,i,iset1),
+     &     wfrag_back(3,i,iset1),
+     &     ifrag_back(1,i,iset1),ifrag_back(2,i,iset1)
         if(me.eq.king.or..not.out1file)
-     &   write(iout,*) "A",i,wfrag_back(1,i,iset),wfrag_back(2,i,iset),
-     &   wfrag_back(3,i,iset),ifrag_back(1,i,iset),ifrag_back(2,i,iset)
+     &   write(iout,*) "A",i,wfrag_back(1,i,iset1),
+     &   wfrag_back(2,i,iset1),
+     &   wfrag_back(3,i,iset1),ifrag_back(1,i,iset1),
+     &   ifrag_back(2,i,iset1)
        enddo 
       enddo
       return