Epetra Package Browser (Single Doxygen Collection) Development
example/MapColoring/cxx_main.cpp
Go to the documentation of this file.
00001 //@HEADER
00002 // ************************************************************************
00003 //
00004 //               Epetra: Linear Algebra Services Package
00005 //                 Copyright 2011 Sandia Corporation
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ************************************************************************
00040 //@HEADER
00041 
00042 #include <assert.h>
00043 #include <stdio.h>
00044 #include <set>
00045 #include <algorithm>
00046 #include <vector>
00047 #include "Epetra_Map.h"
00048 #include "Epetra_BlockMap.h"
00049 #include "Epetra_MapColoring.h"
00050 #include "Epetra_CrsMatrix.h"
00051 #ifdef EPETRA_MPI
00052 #include "mpi.h"
00053 #include "Epetra_MpiComm.h"
00054 #else
00055 #include "Epetra_SerialComm.h"
00056 #endif
00057 
00058 using namespace std;
00059 
00060 int main(int argc, char * argv[]) {
00061 
00062 // Set up the Epetra communicator
00063 #ifdef EPETRA_MPI
00064   MPI_Init(&argc, &argv);
00065   int size;      // Number of MPI processes
00066   int rank;      // My process ID
00067   MPI_Comm_size(MPI_COMM_WORLD, &size);
00068   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00069   Epetra_MpiComm Comm(MPI_COMM_WORLD);
00070 #else
00071   int size = 1;  // Serial case (not using MPI)
00072   int rank = 0;  // Processor ID
00073   Epetra_SerialComm Comm;
00074 #endif
00075   //  cout << Comm << endl;
00076 
00077   int MyPID    = Comm.MyPID();
00078   int NumProc  = Comm.NumProc();
00079   bool verbose = (MyPID == 0);
00080   cout << "MyPID   = " << MyPID   << endl
00081        << "NumProc = " << NumProc << endl;
00082 
00083   // Get the problem size from the command line argument
00084   if (argc < 2 || argc > 3) {
00085     if (verbose)
00086       cout << "Usage: " << argv[0] << " nx [ny]" << endl;
00087     exit(1);
00088   } // end if
00089   int nx = atoi(argv[1]);            // Get the dimensions for a 1D or 2D
00090   int ny = 1;                        // central difference problem
00091   if (argc == 3) ny = atoi(argv[2]);
00092   int NumGlobalElements = nx * ny;
00093   if (NumGlobalElements < NumProc) {
00094     if (verbose)
00095       cout << "numGlobalElements = " << NumGlobalElements
00096      << " cannot be < number of processors = " << NumProc << endl;
00097     exit(1);
00098   } // end if
00099 
00100   // Epetra distribution map
00101   int IndexBase = 0;       // Zero-based indices
00102   Epetra_Map Map(NumGlobalElements, IndexBase, Comm);
00103   //  if (verbose) cout << Map << endl;
00104 
00105   // Extract the global indices of the elements local to this processor
00106   int NumMyElements = Map.NumMyElements();
00107   int * MyGlobalElements = new int[NumMyElements];
00108   Map.MyGlobalElements(MyGlobalElements);
00109   for (int p = 0; p < NumProc; p++)
00110     if (p == MyPID) {
00111       cout << endl << "Processor " << MyPID << ": Global elements = ";
00112       for (int i = 0; i < NumMyElements; i++)
00113   cout << MyGlobalElements[i] << " ";
00114       cout << endl;
00115       Comm.Barrier();
00116     } // end if
00117 
00118   // Create the number of non-zeros for a tridiagonal (1D problem) or banded
00119   // (2D problem) matrix
00120   int * NumNz = new int[NumMyElements];
00121   int global_i;
00122   int global_j;
00123   for (int i = 0; i < NumMyElements; i++) {
00124     NumNz[i] = 5;
00125     global_j = MyGlobalElements[i] / nx;
00126     global_i = MyGlobalElements[i] - global_j * nx;
00127     if (global_i == 0)    NumNz[i] -= 1;  // By having separate statements,
00128     if (global_i == nx-1) NumNz[i] -= 1;  // this works for 2D as well as 1D
00129     if (global_j == 0)    NumNz[i] -= 1;  // systems (i.e. nx x 1 or 1 x ny)
00130     if (global_j == ny-1) NumNz[i] -= 1;  // or even a 1 x 1 system
00131   }
00132   if (verbose) {
00133     cout << endl << "NumNz: ";
00134     for (int i = 0; i < NumMyElements; i++) cout << NumNz[i] << " ";
00135     cout << endl;
00136   } // end if
00137 
00138   // Create the Epetra Compressed Row Sparse matrix
00139   // Note: the actual values for the matrix entries are meaningless for
00140   // this exercise, but I assign them anyway.
00141   Epetra_CrsMatrix A(Copy, Map, NumNz);
00142 
00143   double * Values = new double[4];
00144   for (int i = 0; i < 4; i++) Values[i] = -1.0;
00145   int * Indices = new int[4];
00146   double diag = 2.0;
00147   if (ny > 1) diag = 4.0;
00148   int NumEntries;
00149 
00150   for (int i = 0; i < NumMyElements; i++) {
00151     global_j = MyGlobalElements[i] / nx;
00152     global_i = MyGlobalElements[i] - global_j * nx;
00153     NumEntries = 0;
00154     if (global_j > 0 && ny > 1)
00155       Indices[NumEntries++] = global_i   + (global_j-1)*nx;
00156     if (global_i > 0)
00157       Indices[NumEntries++] = global_i-1 +  global_j   *nx;
00158     if (global_i < nx-1)
00159       Indices[NumEntries++] = global_i+1 +  global_j   *nx;
00160     if (global_j < ny-1 && ny > 1)
00161       Indices[NumEntries++] = global_i   + (global_j+1)*nx;
00162 
00163     // Put in the off-diagonal terms
00164     assert(A.InsertGlobalValues(MyGlobalElements[i], NumEntries, Values,
00165         Indices) == 0);
00166     // Put in the diagonal entry
00167     assert(A.InsertGlobalValues(MyGlobalElements[i], 1, &diag,
00168         MyGlobalElements+i) == 0);
00169   } // end i loop
00170 
00171   // Finish up matrix construction
00172   delete [] Values;
00173   delete [] Indices;
00174   assert(A.FillComplete() == 0);
00175   // cout << endl << A << endl;
00176 
00177   // Create the local distance-1 adjancency graph
00178   // This is essentially a transpose of the Epetra_CrsGraph, where off-
00179   // processor couplings are ignored and global indexes are converted to
00180   // local.  We use the C++ standard libraries vector and set, since we
00181   // don't know how many nonzeroes we will end up with for each column.
00182 
00183   vector< set<int> > adj1(NumMyElements);
00184   for (int lr = 0; lr < adj1.size(); lr++) {
00185     int lrid;   // Local row ID
00186     double * Values = new double[NumNz[lr]];
00187     int * Indices = new int[NumNz[lr]];
00188     assert(A.ExtractMyRowCopy(lr, NumNz[lr], NumNz[lr], Values, Indices) == 0);
00189     for (int i = 0; i < NumNz[lr]; i++) {
00190       lrid = A.LRID(Indices[i]);
00191       if (lrid >= 0) adj1[lrid].insert(lr);
00192     } // end i loop
00193     delete [] Values;
00194     delete [] Indices;
00195   } // end lr loop
00196 
00197   if (verbose) {
00198     cout << endl;
00199     for (int lr = 0; lr < NumMyElements; lr++) {
00200       cout << "adj1[" << lr << "] = { ";
00201       for (set<int>::const_iterator p = adj1[lr].begin(); p != adj1[lr].end();
00202      p++) cout << *p << " ";
00203       cout << "}" << endl;
00204     } // end lr loop
00205   } // end if
00206 
00207   // Create the local distance-2 adjancency graph
00208   // This is built from the distance-1 adjancency graph.  We use the C++
00209   // standard libraries vector and set, since we don't know how many
00210   // nonzeroes we will end up with for each column.
00211 
00212   vector< set<int> > adj2(NumMyElements);
00213   for (int lc = 0; lc < NumMyElements; lc++) {
00214     for (set<int>::const_iterator p = adj1[lc].begin(); p != adj1[lc].end();
00215    p++) {
00216       int lrid;    // Local row ID
00217       double * Values = new double[NumNz[*p]];
00218       int * Indices = new int[NumNz[*p]];
00219       assert(A.ExtractMyRowCopy(*p, NumNz[*p], NumNz[*p], Values, Indices)
00220        == 0);
00221       for (int i = 0; i < NumNz[*p]; i++) {
00222   lrid = A.LRID(Indices[i]);
00223   if (lrid >= 0) adj2[lc].insert(lrid);
00224       } // end i loop
00225       delete [] Values;
00226       delete [] Indices;
00227     } // end p loop
00228   } // end lc loop
00229 
00230   cout << endl;
00231   for (int lc = 0; lc < NumMyElements; lc++) {
00232     cout << "adj2[" << lc << "] = { ";
00233     for (set<int>::const_iterator p = adj2[lc].begin(); p != adj2[lc].end();
00234    p++) cout << *p << " ";
00235     cout << "}" << endl;
00236   } // end lc loop
00237 
00238   // Now that we have the local distance-2 adjacency graph, we can compute a
00239   // color map using a greedy algorithm.  The first step is to compute Delta,
00240   // the maximum size (degree) of adj1.
00241   size_t Delta = 0;
00242   for (int i = 0; i < NumMyElements; i++)
00243     Delta = max(Delta, adj1[i].size());
00244   cout << endl << "Delta = " << Delta << endl << endl;
00245 
00246   // Now create a color map and initialize all values to 0, which
00247   // indicates that none of the columns have yet been colored.
00248   int * color_map = new int[NumMyElements];
00249   for (int i = 0; i < NumMyElements; i++) color_map[i] = 0;
00250 
00251   // Apply the distance-2 greedy coloring algorithm
00252   for (int column = 0; column < NumMyElements; column++) {
00253     set<int> allowedColors;    // Create the set of allowed colors
00254     for (int i = 1; i < Delta*Delta+1; i++) allowedColors.insert(i);
00255     for (set<int>::const_iterator p = adj2[column].begin();
00256    p != adj2[column].end(); p++) if (color_map[*p] > 0)
00257      allowedColors.erase(color_map[*p]);
00258     color_map[column] = *(allowedColors.begin());
00259     cout << "color_map[" << column << "] = " << color_map[column] << endl;
00260   } // end col loop
00261 
00262   // New section to Epetra_MapColoring
00263   Epetra_MapColoring C1(Map, color_map);
00264 
00265   cout << C1;
00266 
00267   // Clean up
00268 
00269   delete [] MyGlobalElements;
00270   delete [] NumNz;
00271   delete [] color_map;
00272 
00273   cout << endl << argv[0] << " done." << endl;
00274 
00275 } // end main
00276 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines