00001
00002 #include "mpi.h"
00003 #include "petsc.h"
00004 #include "sys.h"
00005 #include "parUtils.h"
00006 #include "octUtils.h"
00007 #include "TreeNode.h"
00008 #include "externVars.h"
00009 #include "dendro.h"
00010
00011 #ifdef MPI_WTIME_IS_GLOBAL
00012 #undef MPI_WTIME_IS_GLOBAL
00013 #endif
00014
00015 double gaussian(double mean, double std_deviation);
00016
00017 int main(int argc, char ** argv ) {
00018 int size, rank;
00019 double startTime, endTime;
00020 double localTime, globalTime;
00021 double gSize[3];
00022 unsigned int local_num_pts = 5000;
00023 unsigned int dim = 3;
00024 unsigned int maxDepth = 30;
00025 unsigned int maxNumPts = 1;
00026 bool incCorner = 1;
00027 std::vector<ot::TreeNode> linOct;
00028 std::vector<ot::TreeNode> balOct;
00029 std::vector<ot::TreeNode> coarseOct;
00030 std::vector<ot::TreeNode> tmpNodes;
00031 std::vector<double> pts;
00032 unsigned int ptsLen;
00033 DendroIntL localSz, totalSz;
00034
00035 PetscInitialize(&argc, &argv, "options", NULL);
00036 ot::RegisterEvents();
00037
00038 #ifdef PETSC_USE_LOG
00039 int stages[4];
00040 PetscLogStageRegister(&stages[0], "Main");
00041 PetscLogStageRegister(&stages[1], "P2O");
00042 PetscLogStageRegister(&stages[2], "Bal");
00043 PetscLogStageRegister(&stages[3], "Coarsen");
00044 #endif
00045
00046 #ifdef PETSC_USE_LOG
00047 PetscLogStagePush(stages[0]);
00048 #endif
00049
00050 MPI_Comm_size(MPI_COMM_WORLD, &size);
00051 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00052
00053 if(argc > 1) {
00054 local_num_pts = atoi(argv[1]);
00055 }
00056
00057 if(argc > 2) {
00058 dim = atoi(argv[2]);
00059 }
00060
00061 if(argc > 3) {
00062 maxDepth = atoi(argv[3]);
00063 }
00064
00065 if(argc > 4) {
00066 maxNumPts = atoi(argv[4]);
00067 }
00068
00069 if(argc > 5) {
00070 incCorner = (bool)(atoi(argv[5]));
00071 }
00072
00073 MPI_Barrier(MPI_COMM_WORLD);
00074 if(!rank){
00075 std::cout << "Creating points on the fly... "<<std::endl;
00076 }
00077
00078 startTime = MPI_Wtime();
00079
00080 srand48(0x12345678 + 76543*rank);
00081
00082 pts.resize(3*local_num_pts);
00083 for(int i=0; i<3*local_num_pts; i++) {
00084 pts[i]= gaussian(0.5, 0.16);
00085 }
00086
00087 MPI_Barrier(MPI_COMM_WORLD);
00088 if(!rank){
00089 std::cout << " finished creating points "<<std::endl;
00090 }
00091 endTime = MPI_Wtime();
00092 localTime = endTime - startTime;
00093 if(!rank){
00094 std::cout <<"point generation time: "<<localTime << std::endl;
00095 }
00096
00097 ptsLen = pts.size();
00098 for(int i = 0; i < ptsLen; i+=3) {
00099 if( (pts[i] > 0.0) &&
00100 (pts[i+1] > 0.0)
00101 && (pts[i+2] > 0.0) &&
00102 ( ((unsigned int)(pts[i]*((double)(1u << maxDepth)))) < (1u << maxDepth)) &&
00103 ( ((unsigned int)(pts[i+1]*((double)(1u << maxDepth)))) < (1u << maxDepth)) &&
00104 ( ((unsigned int)(pts[i+2]*((double)(1u << maxDepth)))) < (1u << maxDepth)) ) {
00105 tmpNodes.push_back( ot::TreeNode((unsigned int)(pts[i]*(double)(1u << maxDepth)),
00106 (unsigned int)(pts[i+1]*(double)(1u << maxDepth)),
00107 (unsigned int)(pts[i+2]*(double)(1u << maxDepth)),
00108 maxDepth,dim,maxDepth) );
00109 }
00110 }
00111 pts.clear();
00112
00113 par::removeDuplicates<ot::TreeNode>(tmpNodes,false,MPI_COMM_WORLD);
00114 linOct = tmpNodes;
00115 tmpNodes.clear();
00116 par::partitionW<ot::TreeNode>(linOct, NULL,MPI_COMM_WORLD);
00117
00118 localSz = linOct.size();
00119 par::Mpi_Reduce<DendroIntL>(&localSz, &totalSz, 1, MPI_SUM, 0, MPI_COMM_WORLD);
00120 if(rank==0) {
00121 std::cout<<"# pts= " << totalSz<<std::endl;
00122 }
00123
00124 pts.resize(3*(linOct.size()));
00125 ptsLen = (3*(linOct.size()));
00126 for(int i = 0; i < linOct.size(); i++) {
00127 pts[3*i] = (((double)(linOct[i].getX())) + 0.5)/((double)(1u << maxDepth));
00128 pts[(3*i)+1] = (((double)(linOct[i].getY())) +0.5)/((double)(1u << maxDepth));
00129 pts[(3*i)+2] = (((double)(linOct[i].getZ())) +0.5)/((double)(1u << maxDepth));
00130 }
00131 linOct.clear();
00132
00133 gSize[0] = 1.0;
00134 gSize[1] = 1.0;
00135 gSize[2] = 1.0;
00136
00137 #ifdef PETSC_USE_LOG
00138 PetscLogStagePop();
00139 #endif
00140
00141 MPI_Barrier(MPI_COMM_WORLD);
00142
00143 #ifdef PETSC_USE_LOG
00144 PetscLogStagePush(stages[1]);
00145 #endif
00146 startTime = MPI_Wtime();
00147 ot::points2Octree(pts, gSize, linOct, dim, maxDepth, maxNumPts, MPI_COMM_WORLD);
00148 endTime = MPI_Wtime();
00149 localTime = endTime - startTime;
00150 #ifdef PETSC_USE_LOG
00151 PetscLogStagePop();
00152 #endif
00153 par::Mpi_Reduce<double>(&localTime, &globalTime, 1, MPI_MAX, 0, MPI_COMM_WORLD);
00154 if(!rank){
00155 std::cout <<"P2n Time: "<<globalTime << "secs " << std::endl;
00156 }
00157 pts.clear();
00158
00159 localSz = linOct.size();
00160 par::Mpi_Reduce<DendroIntL>(&localSz, &totalSz, 1, MPI_SUM, 0, MPI_COMM_WORLD);
00161 if(rank==0) {
00162 std::cout<<"linOct.size = " << totalSz<<std::endl;
00163 }
00164
00165 MPI_Barrier(MPI_COMM_WORLD);
00166
00167 #ifdef PETSC_USE_LOG
00168 PetscLogStagePush(stages[2]);
00169 #endif
00170 startTime = MPI_Wtime();
00171 ot::balanceOctree (linOct, balOct, dim, maxDepth, incCorner, MPI_COMM_WORLD, NULL, NULL);
00172 endTime = MPI_Wtime();
00173 localTime = endTime - startTime;
00174 #ifdef PETSC_USE_LOG
00175 PetscLogStagePop();
00176 #endif
00177 par::Mpi_Reduce<double>(&localTime, &globalTime, 1, MPI_MAX, 0, MPI_COMM_WORLD);
00178 if(!rank){
00179 std::cout <<"Bal Time: "<<globalTime << "secs " << std::endl;
00180 }
00181 linOct.clear();
00182
00183 localSz = balOct.size();
00184 par::Mpi_Reduce<DendroIntL>(&localSz, &totalSz, 1, MPI_SUM, 0, MPI_COMM_WORLD);
00185 if(rank==0) {
00186 std::cout<<"balOct.size = " << totalSz<<std::endl;
00187 }
00188 #ifdef PETSC_USE_LOG
00189 PetscLogStagePush(stages[3]);
00190 #endif
00191 startTime = MPI_Wtime();
00192 ot::coarsenOctree(balOct, coarseOct, dim, maxDepth, MPI_COMM_WORLD, false, NULL, NULL);
00193 endTime = MPI_Wtime();
00194 localTime = endTime - startTime;
00195 #ifdef PETSC_USE_LOG
00196 PetscLogStagePop();
00197 #endif
00198 par::Mpi_Reduce<double>(&localTime, &globalTime, 1, MPI_MAX, 0, MPI_COMM_WORLD);
00199 if(!rank){
00200 std::cout <<"coarse Time: "<<globalTime << "secs " << std::endl;
00201 }
00202 balOct.clear();
00203
00204 localSz = coarseOct.size();
00205 par::Mpi_Reduce<DendroIntL>(&localSz, &totalSz, 1, MPI_SUM, 0, MPI_COMM_WORLD);
00206 if(rank==0) {
00207 std::cout<<"coarseOct.size = " << totalSz<<std::endl;
00208 }
00209 coarseOct.clear();
00210
00211 PetscFinalize();
00212 }
00213
00214 double gaussian(double mean, double std_deviation) {
00215 static double t1 = 0, t2=0;
00216 double x1, x2, x3, r;
00217
00218 using namespace std;
00219
00220
00221 if(t1) {
00222 const double tmp = t1;
00223 t1 = 0;
00224 return mean + std_deviation * tmp;
00225 }
00226 if(t2) {
00227 const double tmp = t2;
00228 t2 = 0;
00229 return mean + std_deviation * tmp;
00230 }
00231
00232
00233 do {
00234 x1 = 2 * drand48() - 1;
00235 x2 = 2 * drand48() - 1;
00236 x3 = 2 * drand48() - 1;
00237 r = x1 * x1 + x2 * x2 + x3*x3;
00238 } while(r >= 1);
00239
00240
00241 r = sqrt(-2.0 * log(r) / r);
00242
00243
00244 t1 = (r * x2);
00245 t2 = (r * x3);
00246
00247 return mean + (std_deviation * r * x1);
00248 }
00249
00250