|
RInside Version 0.2.10
|
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 // This file was contributed by Jianping Hua 00006 // 00007 // Copyright (C) 2010 - 2011 Jianping Hua, 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 related 00017 int myrank, nodesize; // node information 00018 int sndcnt = 1, rcvcnt = 1; // # of elements in send/recv buffer 00019 MPI_Init(&argc,&argv); // mpi initialization 00020 MPI_Comm_rank(MPI_COMM_WORLD, &myrank); // obtain current node rank 00021 MPI_Comm_size(MPI_COMM_WORLD, &nodesize); // obtain total nodes running. 00022 double sendValue; // value to be collected in current node 00023 double *allvalues = new double[nodesize]; // to save all results 00024 00025 // simulation info 00026 // to sample from a uniform distribution 00027 int rangeMin = 0, rangeMax = 10; // range of uniform distribution 00028 int sampleSize = 2; // points in each sample 00029 00030 try { 00031 RInside R(argc, argv); // create an embedded R instance 00032 00033 std::stringstream txt; 00034 txt << "x <- " << rangeMin << std::endl; 00035 R.parseEvalQ( txt.str() ); // assign x with lower range of uniform distribution 00036 00037 txt << "y <- " << rangeMax << std::endl; 00038 R.parseEvalQ( txt.str() ); // assign y with upper range of uniform distribution 00039 00040 txt << "n <- " << sampleSize << std::endl; 00041 R.parseEvalQ( txt.str() ); // assign n with the size of sample 00042 00043 std::string evalstr = "mean(runif(n,x,y))"; // sampling, compute the mean 00044 Rcpp::NumericVector m = R.parseEval(evalstr); // eval str, convert result to NumericVector 00045 sendValue = m( 0 ); // assign the return value to the variable to be gathered 00046 00047 //gather together values from all processes to allvalues 00048 MPI_Gather(&sendValue, sndcnt, MPI_DOUBLE, allvalues, rcvcnt, MPI_DOUBLE, 0, MPI_COMM_WORLD); 00049 00050 // show what inidividual node's contribution 00051 std::cout << "node " << myrank << " has mean " << m(0) << std::endl; 00052 00053 } catch(std::exception& ex) { 00054 std::cerr << "Exception caught: " << ex.what() << std::endl; 00055 } catch(...) { 00056 std::cerr << "Unknown exception caught" << std::endl; 00057 } 00058 00059 // show gathered results in node 0 00060 if ( myrank == 0 ) { 00061 std::cout << "values of all " << nodesize << " trials: " << std::endl; 00062 for ( int i = 0; i < nodesize; i++ ) 00063 std::cout << allvalues[ i ] << ", "; 00064 std::cout << std::endl; 00065 } 00066 00067 // clean up 00068 delete[] allvalues; 00069 00070 MPI_Finalize(); // mpi finalization 00071 00072 exit(0); 00073 } 00074