example/MapColoring/cxx_main.cpp

Go to the documentation of this file.
00001 //@HEADER
00002 // ************************************************************************
00003 // 
00004 //               Epetra: Linear Algebra Services Package 
00005 //                 Copyright (2001) Sandia Corporation
00006 // 
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 // 
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //  
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //  
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 
00025 // 
00026 // ************************************************************************
00027 //@HEADER
00028 
00029 #include <assert.h>
00030 #include <stdio.h>
00031 #include <set>
00032 #include <algorithm>
00033 #include <vector>
00034 #include "Epetra_Map.h"
00035 #include "Epetra_BlockMap.h"
00036 #include "Epetra_MapColoring.h"
00037 #include "Epetra_CrsMatrix.h"
00038 #ifdef EPETRA_MPI
00039 #include "mpi.h"
00040 #include "Epetra_MpiComm.h"
00041 #else
00042 #include "Epetra_SerialComm.h"
00043 #endif
00044 
00045 using namespace std;
00046 
00047 int main(int argc, char * argv[]) {
00048 
00049 // Set up the Epetra communicator
00050 #ifdef EPETRA_MPI
00051   MPI_Init(&argc, &argv);
00052   int size;      // Number of MPI processes
00053   int rank;      // My process ID
00054   MPI_Comm_size(MPI_COMM_WORLD, &size);
00055   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00056   Epetra_MpiComm Comm(MPI_COMM_WORLD);
00057 #else
00058   int size = 1;  // Serial case (not using MPI)
00059   int rank = 0;  // Processor ID
00060   Epetra_SerialComm Comm;
00061 #endif
00062   //  cout << Comm << endl;
00063 
00064   int MyPID    = Comm.MyPID();
00065   int NumProc  = Comm.NumProc();
00066   bool verbose = (MyPID == 0);
00067   cout << "MyPID   = " << MyPID   << endl
00068        << "NumProc = " << NumProc << endl;
00069 
00070   // Get the problem size from the command line argument
00071   if (argc < 2 || argc > 3) {
00072     if (verbose)
00073       cout << "Usage: " << argv[0] << " nx [ny]" << endl;
00074     exit(1);
00075   } // end if
00076   int nx = atoi(argv[1]);            // Get the dimensions for a 1D or 2D
00077   int ny = 1;                        // central difference problem
00078   if (argc == 3) ny = atoi(argv[2]);
00079   int NumGlobalElements = nx * ny;
00080   if (NumGlobalElements < NumProc) {
00081     if (verbose)
00082       cout << "numGlobalElements = " << NumGlobalElements 
00083      << " cannot be < number of processors = " << NumProc << endl;
00084     exit(1);
00085   } // end if
00086 
00087   // Epetra distribution map
00088   int IndexBase = 0;       // Zero-based indices
00089   Epetra_Map Map(NumGlobalElements, IndexBase, Comm);
00090   //  if (verbose) cout << Map << endl;
00091 
00092   // Extract the global indices of the elements local to this processor
00093   int NumMyElements = Map.NumMyElements();
00094   int * MyGlobalElements = new int[NumMyElements];
00095   Map.MyGlobalElements(MyGlobalElements);
00096   for (int p = 0; p < NumProc; p++) 
00097     if (p == MyPID) {
00098       cout << endl << "Processor " << MyPID << ": Global elements = ";
00099       for (int i = 0; i < NumMyElements; i++)
00100   cout << MyGlobalElements[i] << " ";
00101       cout << endl;
00102       Comm.Barrier();
00103     } // end if
00104 
00105   // Create the number of non-zeros for a tridiagonal (1D problem) or banded
00106   // (2D problem) matrix
00107   int * NumNz = new int[NumMyElements];
00108   int global_i;
00109   int global_j;
00110   for (int i = 0; i < NumMyElements; i++) {
00111     NumNz[i] = 5;
00112     global_j = MyGlobalElements[i] / nx;
00113     global_i = MyGlobalElements[i] - global_j * nx;
00114     if (global_i == 0)    NumNz[i] -= 1;  // By having separate statements,
00115     if (global_i == nx-1) NumNz[i] -= 1;  // this works for 2D as well as 1D
00116     if (global_j == 0)    NumNz[i] -= 1;  // systems (i.e. nx x 1 or 1 x ny)
00117     if (global_j == ny-1) NumNz[i] -= 1;  // or even a 1 x 1 system
00118   }
00119   if (verbose) {
00120     cout << endl << "NumNz: ";
00121     for (int i = 0; i < NumMyElements; i++) cout << NumNz[i] << " ";
00122     cout << endl;
00123   } // end if
00124 
00125   // Create the Epetra Compressed Row Sparse matrix
00126   // Note: the actual values for the matrix entries are meaningless for
00127   // this exercise, but I assign them anyway.
00128   Epetra_CrsMatrix A(Copy, Map, NumNz);
00129 
00130   double * Values = new double[4];
00131   for (int i = 0; i < 4; i++) Values[i] = -1.0;
00132   int * Indices = new int[4];
00133   double diag = 2.0;
00134   if (ny > 1) diag = 4.0;
00135   int NumEntries;
00136 
00137   for (int i = 0; i < NumMyElements; i++) {
00138     global_j = MyGlobalElements[i] / nx;
00139     global_i = MyGlobalElements[i] - global_j * nx;
00140     NumEntries = 0;
00141     if (global_j > 0 && ny > 1)
00142       Indices[NumEntries++] = global_i   + (global_j-1)*nx;
00143     if (global_i > 0)
00144       Indices[NumEntries++] = global_i-1 +  global_j   *nx;
00145     if (global_i < nx-1)
00146       Indices[NumEntries++] = global_i+1 +  global_j   *nx;
00147     if (global_j < ny-1 && ny > 1)
00148       Indices[NumEntries++] = global_i   + (global_j+1)*nx;
00149 
00150     // Put in the off-diagonal terms
00151     assert(A.InsertGlobalValues(MyGlobalElements[i], NumEntries, Values,
00152         Indices) == 0);
00153     // Put in the diagonal entry
00154     assert(A.InsertGlobalValues(MyGlobalElements[i], 1, &diag,
00155         MyGlobalElements+i) == 0);
00156   } // end i loop
00157 
00158   // Finish up matrix construction
00159   delete [] Values;
00160   delete [] Indices;
00161   assert(A.FillComplete() == 0);
00162   // cout << endl << A << endl;
00163 
00164   // Create the local distance-1 adjancency graph
00165   // This is essentially a transpose of the Epetra_CrsGraph, where off-
00166   // processor couplings are ignored and global indexes are converted to
00167   // local.  We use the C++ standard libraries vector and set, since we
00168   // don't know how many nonzeroes we will end up with for each column.
00169 
00170   vector< set<int> > adj1(NumMyElements);
00171   for (int lr = 0; lr < adj1.size(); lr++) {
00172     int lrid;   // Local row ID
00173     double * Values = new double[NumNz[lr]];
00174     int * Indices = new int[NumNz[lr]];
00175     assert(A.ExtractMyRowCopy(lr, NumNz[lr], NumNz[lr], Values, Indices) == 0);
00176     for (int i = 0; i < NumNz[lr]; i++) {
00177       lrid = A.LRID(Indices[i]);
00178       if (lrid >= 0) adj1[lrid].insert(lr);
00179     } // end i loop
00180     delete [] Values;
00181     delete [] Indices;
00182   } // end lr loop
00183 
00184   if (verbose) {
00185     cout << endl;
00186     for (int lr = 0; lr < NumMyElements; lr++) {
00187       cout << "adj1[" << lr << "] = { ";
00188       for (set<int>::const_iterator p = adj1[lr].begin(); p != adj1[lr].end();
00189      p++) cout << *p << " ";
00190       cout << "}" << endl;
00191     } // end lr loop
00192   } // end if
00193 
00194   // Create the local distance-2 adjancency graph
00195   // This is built from the distance-1 adjancency graph.  We use the C++
00196   // standard libraries vector and set, since we don't know how many
00197   // nonzeroes we will end up with for each column.
00198 
00199   vector< set<int> > adj2(NumMyElements);
00200   for (int lc = 0; lc < NumMyElements; lc++) {
00201     for (set<int>::const_iterator p = adj1[lc].begin(); p != adj1[lc].end();
00202    p++) {
00203       int lrid;    // Local row ID
00204       double * Values = new double[NumNz[*p]];
00205       int * Indices = new int[NumNz[*p]];
00206       assert(A.ExtractMyRowCopy(*p, NumNz[*p], NumNz[*p], Values, Indices)
00207        == 0);
00208       for (int i = 0; i < NumNz[*p]; i++) {
00209   lrid = A.LRID(Indices[i]);
00210   if (lrid >= 0) adj2[lc].insert(lrid);
00211       } // end i loop
00212       delete [] Values;
00213       delete [] Indices;
00214     } // end p loop
00215   } // end lc loop
00216 
00217   cout << endl;
00218   for (int lc = 0; lc < NumMyElements; lc++) {
00219     cout << "adj2[" << lc << "] = { ";
00220     for (set<int>::const_iterator p = adj2[lc].begin(); p != adj2[lc].end();
00221    p++) cout << *p << " ";
00222     cout << "}" << endl;
00223   } // end lc loop
00224 
00225   // Now that we have the local distance-2 adjacency graph, we can compute a
00226   // color map using a greedy algorithm.  The first step is to compute Delta,
00227   // the maximum size (degree) of adj1.
00228   size_t Delta = 0;
00229   for (int i = 0; i < NumMyElements; i++)
00230     Delta = max(Delta, adj1[i].size());
00231   cout << endl << "Delta = " << Delta << endl << endl;
00232   
00233   // Now create a color map and initialize all values to 0, which
00234   // indicates that none of the columns have yet been colored.
00235   int * color_map = new int[NumMyElements];
00236   for (int i = 0; i < NumMyElements; i++) color_map[i] = 0;
00237 
00238   // Apply the distance-2 greedy coloring algorithm
00239   for (int column = 0; column < NumMyElements; column++) {
00240     set<int> allowedColors;    // Create the set of allowed colors
00241     for (int i = 1; i < Delta*Delta+1; i++) allowedColors.insert(i);
00242     for (set<int>::const_iterator p = adj2[column].begin();
00243    p != adj2[column].end(); p++) if (color_map[*p] > 0)
00244      allowedColors.erase(color_map[*p]);
00245     color_map[column] = *(allowedColors.begin());
00246     cout << "color_map[" << column << "] = " << color_map[column] << endl;
00247   } // end col loop
00248 
00249   // New section to Epetra_MapColoring
00250   Epetra_MapColoring C1(Map, color_map);
00251 
00252   cout << C1;
00253   
00254   // Clean up
00255 
00256   delete [] MyGlobalElements;
00257   delete [] NumNz;
00258   delete [] color_map;
00259 
00260   cout << endl << argv[0] << " done." << endl;
00261 
00262 } // end main
00263 

Generated on Thu Sep 18 12:37:56 2008 for Epetra Package Browser (Single Doxygen Collection) by doxygen 1.3.9.1