]> git.donarmstrong.com Git - libparallel-mpi-simple-perl.git/commitdiff
add support for intercomm merge, comm_self comm_get_parent and fix a few problems... features/mpi_comm_spawn
authorDon Armstrong <don@donarmstrong.com>
Thu, 7 Jun 2012 23:32:10 +0000 (16:32 -0700)
committerDon Armstrong <don@donarmstrong.com>
Thu, 7 Jun 2012 23:32:10 +0000 (16:32 -0700)
Simple.pm
Simple.xs

index ea709cc1c4081390cc2a54ecb547c766347c5dc5..8fc58b3682c9cd8515f507deb5ee41a735d66ff6 100644 (file)
--- a/Simple.pm
+++ b/Simple.pm
@@ -16,10 +16,11 @@ sub import {
     my $call = (caller())[0];
     no strict 'refs';
     # subs (MPI_ function calls)
-    foreach (qw(Init Finalize COMM_WORLD ANY_SOURCE ANY_TAG Comm_rank Comm_size
+    foreach (qw(Init Finalize COMM_WORLD COMM_SELF ANY_SOURCE ANY_TAG Comm_rank Comm_size
                Recv Send Barrier Bcast Gather
                Scatter Allgather Alltoall Reduce
                Comm_compare Comm_dup Comm_free Comm_split Comm_spawn
+               Intercomm_merge Comm_get_parent
               )) {
        *{$call.'::MPI_'.$_} = \&$_;
     }
@@ -515,12 +516,53 @@ communicator which communicates between the parents and the children.
 
 =cut
 
-
 sub Comm_spawn {
     my ($maxprocs,$comm,$root,$command,@arguments) = @_;
-    return _Comm_spawn($command,@arguments,$maxprocs,$root,$comm);
+    return _Comm_spawn($command,\@arguments,$maxprocs,$root,$comm);
+}
+
+
+=head2 MPI_Comm_spawn
+
+ $parent_comm = MPI_Comm_get_parent();
+
+Returns the parent intercommunicator. This is primarily useful when
+merging intercommunicators into one intracommunicator. [For example,
+if you wanted to communicate between spawned children and a parent
+process.]
+
+=cut
+
+sub Comm_get_parent{
+    return _Comm_get_parent();
 }
 
+=head2 MPI_Intercomm_merge
+
+ my $intracomm = MPI_Intercomm_merge($intercomm,$high);
+
+Returns a merged intracommunicator given an intercomm. If $high is
+true, the rank in the merged intracommunicator of processes will be
+higher than those with high false. If high is the same value in all
+merging processes, the order is undefined.
+
+An example:
+ # in the master process
+ my $intercomm = MPI_Comm_spawn(1,MPI_COMM_SELF,$rank,"test2.pl");
+ my $intracomm = MPI_Intercomm_merge($intercomm,0);
+ print(MPI_Recv(1,0,$intracomm));
+
+ # in test2.pl
+ my $intercomm = MPI_Comm_get_parent();
+ my $intracomm = MPI_Intercomm_merge($intercomm,1);
+ print(MPI_Send("Hello world!",0,0,$intracomm));
+
+=cut
+
+sub Intercomm_merge {
+    my ($comm,$high) = @_;
+    return _Intercomm_merge($comm,$high?1:0);
+}
 
 1; # I am the ANTI-POD!
 
index 8090de30eda9c257a08656aeff6f7a1848d281e7..74d2f8fdb304ef4279b6f9629109938aa7639519 100644 (file)
--- a/Simple.xs
+++ b/Simple.xs
@@ -95,6 +95,21 @@ int _Comm_size (SV* comm) {
   return tsize;
 }
 
+SV* _Comm_get_parent () {
+  MPI_Comm parent;
+  MPI_Comm_get_parent(&parent);
+  return newSViv((IV)parent);
+}
+
+SV* _Intercomm_merge(SV* intercomm, int high) {
+  MPI_Comm newintracomm;
+
+  MPI_Intercomm_merge((MPI_Comm)SvIVX(intercomm),high,&newintracomm);
+  return newSViv((IV)newintracomm);
+}
+
+
+
 /* spawns a command running on other hosts */
 SV* _Comm_spawn(SV* command, AV* argv, int maxprocs, int root, SV* comm) {
   MPI_Comm intercomm;
@@ -106,34 +121,36 @@ SV* _Comm_spawn(SV* command, AV* argv, int maxprocs, int root, SV* comm) {
   int i;
   int size;
   
-  /* create the argv needed to pass to MPI_comm_spawn */
+  /* create the argv needed to pass to MPI_Comm_spawn */
   len = av_len(argv);
   Newx(argv_new,len < 0 ? 1 : len+2,char*);
+  Newx(array_of_errcodes,maxprocs,int);
+  Newx(intercomm,1,MPI_Comm);
   for (key = 0; key <= len; key++) {
     argv_new[key]=SvPV_nolen(*av_fetch(argv,key,0));
   }
-  Newx(argv_new[len < 0 ? 0:len+1],1,char);
-  argv_new[len < 0 ? 0:len+1][0] = 0;
+  argv_new[len < 0 ? 0:len+1] = 0;
 
   /* eventually we should handle MPI_INFO, but since I don't need it
      yet, not bothering. */
-  error = MPI_comm_spawn(SvPV_nolen(command),
+  error = MPI_Comm_spawn(SvPV_nolen(command),
                         argv_new,maxprocs,
                         MPI_INFO_NULL,root,
                         (MPI_Comm)SvIVX(comm),
-                        &intercomm,&array_of_errcodes);
-  Safefree(argv_new[len < 0 ? 0:len+1]);
+                        &intercomm,array_of_errcodes);
   Safefree(argv_new);
   if (error != 0)
     croak("Unable to spawn process");
   /* figure out how many processes were spawned, and check to see if there were errors */
   MPI_Comm_size(intercomm, &size);
   for(i = 0; i < size; i++) {
-    if (array_of_errcodes[i] != 0)
+    if (array_of_errcodes[i] != 0) {
+      Safefree(array_of_errcodes);
       croak("Process did not spawn properly");
+    }
   }
+  Safefree(array_of_errcodes);
   return newSViv((IV)intercomm);
-  
 }
 
 /* returns SV whose IV slot is a cast pointer to the MPI_COMM_WORLD object */
@@ -141,16 +158,25 @@ SV* COMM_WORLD () {
   return newSViv((IV)MPI_COMM_WORLD);
 }
 
+
+/* returns SV whose IV slot is a cast pointer to the MPI_COMM_SELF object */
+SV* COMM_SELF () {
+  return newSViv((IV)MPI_COMM_SELF);
+}
+
+
 /* returns SV whose IV slot is a cast pointer to the MPI_ANY_SOURCE value */
 SV* ANY_SOURCE () {
   return newSViv((IV)MPI_ANY_SOURCE);
 }
 
+
 /* returns SV whose IV slot is a cast pointer to the MPI_ANY_TAG value */
 SV* ANY_TAG () {
   return newSViv((IV)MPI_ANY_TAG);
 }
 
+
 /* calls MPI_Barrier for comm */
 int Barrier (SV*comm) {
   MPI_Barrier((MPI_Comm)SvIVX(comm));
@@ -297,15 +323,22 @@ int
 _Comm_size (comm)
        SV *    comm
 
+
 SV *
 COMM_WORLD ()
 
+SV *
+COMM_SELF ()
+
 SV *
 ANY_SOURCE ()
 
 SV *
 ANY_TAG ()
 
+SV*
+_Comm_get_parent()
+
 int
 Barrier (comm)
        SV *    comm
@@ -324,6 +357,12 @@ _Comm_compare (comm1, comm2)
        SV *    comm1
        SV *    comm2
 
+SV*
+_Intercomm_merge (intercomm, high)
+       SV *    intercomm
+       int     high
+
+
 void
 _Comm_free (comm)
        SV *    comm
@@ -350,3 +389,10 @@ _Comm_split (comm, colour, key)
        int     colour
        int     key
 
+SV *
+_Comm_spawn (command, argv, maxprocs, root, comm)
+        SV* command
+        AV* argv
+        int maxprocs
+        int root
+        SV* comm