[Trilinos-Users] Replicate Serial Vectors
Michael A Heroux
maherou at sandia.gov
Thu Aug 4 22:54:15 MDT 2005
Ammar,
Here is a simple example that illustrates one way to replicate a vector on a
parallel machine. I hope you find it useful.
Mike
#include "Epetra_MpiComm.h"
#include "Epetra_Map.h"
#include "Epetra_LocalMap.h"
#include "Epetra_Vector.h"
#include "Epetra_Import.h"
#include "Epetra_IntSerialDenseVector.h"
int main(int argc, char *argv[])
{
// init mpi
int size, rank;
MPI_Init(&argc,&argv);
MPI_Comm_size(MPI_COMM_WORLD, &size);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
Epetra_MpiComm Comm(MPI_COMM_WORLD);
int lenX = 5; // The length of the vector on each processor.
Epetra_LocalMap localMap(lenX, 0, Comm); // This map has GIDs replicated
on each processor.
// Now build source map. It will have GIDs on PE 0 only.
int lenX0 = 0;
if(Comm.MyPID() == 0) lenX0= lenX;
Epetra_IntSerialDenseVector myGIDs(lenX0);
for (int i=0; i<lenX0; i++) myGIDs[i] = i;
Epetra_Map map0(-1, lenX0, myGIDs.Values(), 0, Comm);
// Now build source vector
Epetra_Vector x0(map0);
if (Comm.MyPID()==0) x0.SetLabel("Source Vector"); else x0.SetLabel("");
for (int i=0; i<lenX0; i++) x0[i] = (double) (i+1)*(i+1);
cout << x0 << endl;
// Build import object will get the values from x0 and bring them to all
processors.
Epetra_Import importer(localMap, map0);
// Build replicated vector, do import and finish.
Epetra_Vector x(localMap); // This is our target vector.
x.Import(x0, importer, Insert);
if (Comm.MyPID()==0) x.SetLabel("Replicated Vector"); else
x.SetLabel("");
cout << x << endl;
// finalize mpi
MPI_Finalize();
return(0);
}
More information about the Trilinos-Users
mailing list