|
RInside Version 0.2.6
|
00001 // -*- mode: C++; c-indent-level: 4; c-basic-offset: 4; tab-width: 8; -*- 00002 // 00003 // Simple mpi example: simulate sampling/averaging on multiple nodes and gathering the results. 00004 // 00005 // MPI C++ API version of file contributed by Jianping Hua 00006 // 00007 // Copyright (C) 2010 - 2011 Dirk Eddelbuettel and Romain Francois 00008 // 00009 // GPL'ed 00010 00011 #include <mpi.h> // mpi header file 00012 #include <RInside.h> // for the embedded R via RInside 00013 00014 int main(int argc, char *argv[]) { 00015 00016 MPI::Init(argc, argv); // mpi initialization 00017 int myrank = MPI::COMM_WORLD.Get_rank(); // obtain current node rank 00018 int nodesize = MPI::COMM_WORLD.Get_size(); // obtain total nodes running. 00019 00020 int sndcnt = 1, rcvcnt = 1; // # of elements in send/recv buffer 00021 double sendValue; // value to be collected in current node 00022 double *allvalues = new double[nodesize]; // to save all results 00023 00024 // simulation info 00025 // to sample from a uniform distribution 00026 int rangeMin = 0, rangeMax = 10; // range of uniform distribution 00027 int sampleSize = 2; // points in each sample 00028 00029 try { 00030 RInside R(argc, argv); // create an embedded R instance 00031 00032 std::stringstream txt; 00033 txt << "x <- " << rangeMin << std::endl; 00034 R.parseEvalQ( txt.str() ); // assign x with lower range of uniform distribution 00035 00036 txt << "y <- " << rangeMax << std::endl; 00037 R.parseEvalQ( txt.str() ); // assign y with upper range of uniform distribution 00038 00039 txt << "n <- " << sampleSize << std::endl; 00040 R.parseEvalQ( txt.str() ); // assign n with the size of sample 00041 00042 std::string evalstr = "mean(runif(n,x,y))"; // sampling, compute the mean 00043 Rcpp::NumericVector m = R.parseEval(evalstr); // eval str, convert result to NumericVector 00044 sendValue = m( 0 ); // assign the return value to the variable to be gathered 00045 00046 //gather together values from all processes to allvalues 00047 MPI::COMM_WORLD.Gather((const void*)&sendValue, sndcnt, MPI::DOUBLE, (void*)allvalues, rcvcnt, MPI::DOUBLE, 0); 00048 00049 // show what inidividual node's contribution 00050 std::cout << "node " << myrank << " has mean " << m(0) << std::endl; 00051 00052 } catch(std::exception& ex) { 00053 std::cerr << "Exception caught: " << ex.what() << std::endl; 00054 } catch(...) { 00055 std::cerr << "Unknown exception caught" << std::endl; 00056 } 00057 00058 // show gathered results in node 0 00059 if ( myrank == 0 ) { 00060 std::cout << "values of all " << nodesize << " trials: " << std::endl; 00061 for ( int i = 0; i < nodesize; i++ ) 00062 std::cout << allvalues[ i ] << ", "; 00063 std::cout << std::endl; 00064 } 00065 00066 // clean up 00067 delete[] allvalues; 00068 00069 MPI::Finalize(); // mpi finalization 00070 00071 exit(0); 00072 } 00073