Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/MapColoring/cxx_main.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Epetra: Linear Algebra Services Package
5 // Copyright 2011 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 //@HEADER
41 
42 
43 // Epetra_BlockMap Test routine
44 
45 #include "Epetra_BlockMap.h"
46 #include "Epetra_Map.h"
47 #include "Epetra_MapColoring.h"
48 #include "Epetra_IntVector.h"
49 #include "Epetra_Import.h"
50 #ifdef EPETRA_MPI
51 #include "Epetra_MpiComm.h"
52 #include <mpi.h>
53 #endif
54 
55 #include "Epetra_SerialComm.h"
56 #include "Epetra_Time.h"
57 #include "Epetra_Version.h"
58 
59 int main(int argc, char *argv[]) {
60 
61  int i, returnierr=0;
62 
63 #ifdef EPETRA_MPI
64  // Initialize MPI
65  MPI_Init(&argc,&argv);
66  Epetra_MpiComm Comm(MPI_COMM_WORLD);
67 #else
68  Epetra_SerialComm Comm;
69 #endif
70 
71  // Uncomment to debug in parallel int tmp; if (Comm.MyPID()==0) cin >> tmp; Comm.Barrier();
72 
73  bool verbose = false;
74  bool veryVerbose = false;
75 
76  // Check if we should print results to standard out
77  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;
78 
79  // Check if we should print lots of results to standard out
80  if (argc>2) if (argv[2][0]=='-' && argv[2][1]=='v') veryVerbose = true;
81 
82  if (verbose && Comm.MyPID()==0)
83  std::cout << Epetra_Version() << std::endl << std::endl;
84 
85  if (!verbose) Comm.SetTracebackMode(0); // This should shut down any error traceback reporting
86 
87  if (verbose) std::cout << Comm << std::endl << std::flush;
88 
89  bool verbose1 = verbose;
90  if (verbose) verbose = (Comm.MyPID()==0);
91 
92  bool veryVerbose1 = veryVerbose;
93  if (veryVerbose) veryVerbose = (Comm.MyPID()==0);
94 
95  int NumMyElements = 100;
96  if (veryVerbose1) NumMyElements = 10;
97  NumMyElements += Comm.MyPID();
98  int MaxNumMyElements = NumMyElements+Comm.NumProc()-1;
99  int * ElementSizeList = new int[NumMyElements];
100  int * MyGlobalElements = new int[NumMyElements];
101 
102  for (i = 0; i<NumMyElements; i++) {
103  MyGlobalElements[i] = (Comm.MyPID()*MaxNumMyElements+i)*2;
104  ElementSizeList[i] = i%6 + 2; // elementsizes go from 2 to 7
105  }
106 
107  Epetra_BlockMap Map(-1, NumMyElements, MyGlobalElements, ElementSizeList,
108  0, Comm);
109 
110  delete [] ElementSizeList;
111  delete [] MyGlobalElements;
112 
113  Epetra_MapColoring C0(Map);
114 
115  int * elementColors = new int[NumMyElements];
116 
117  int maxcolor = 24;
118  int * colorCount = new int[maxcolor];
119  int ** colorLIDs = new int*[maxcolor];
120  for (i=0; i<maxcolor; i++) colorCount[i] = 0;
121  for (i=0; i<maxcolor; i++) colorLIDs[i] = 0;
122 
123  int defaultColor = C0.DefaultColor();
124  for (i=0; i<Map.NumMyElements(); i++) {
125  assert(C0[i]==defaultColor);
126  assert(C0(Map.GID(i))==defaultColor);
127  if (i%2==0) C0[i] = i%6+5+i%14; // cycle through 5...23 on even elements
128  else C0(Map.GID(i)) = i%5+1; // cycle through 1...5 on odd elements
129  elementColors[i] = C0[i]; // Record color of ith element for use below
130  colorCount[C0[i]]++; // Count how many of each color for checking below
131  }
132 
133  if (veryVerbose)
134  std::cout << "Original Map Coloring using element-by-element definitions" << std::endl;
135  if (veryVerbose1)
136  std::cout << C0 << std::endl;
137 
138  int numColors = 0;
139  for (i=0; i<maxcolor; i++)
140  if (colorCount[i]>0) {
141  numColors++;
142  colorLIDs[i] = new int[colorCount[i]];
143  }
144  for (i=0; i<maxcolor; i++) colorCount[i] = 0;
145  for (i=0; i<Map.NumMyElements(); i++) colorLIDs[C0[i]][colorCount[C0[i]]++] = i;
146 
147 
148 
149  int newDefaultColor = -1;
150  Epetra_MapColoring C1(Map, elementColors, newDefaultColor);
151  if (veryVerbose)
152  std::cout << "Same Map Coloring using one-time construction" << std::endl;
153  if (veryVerbose1)
154  std::cout << C1 << std::endl;
155  assert(C1.DefaultColor()==newDefaultColor);
156  for (i=0; i<Map.NumMyElements(); i++) assert(C1[i]==C0[i]);
157 
158  Epetra_MapColoring C2(C1);
159  if (veryVerbose)
160  std::cout << "Same Map Coloring using copy constructor" << std::endl;
161  if (veryVerbose1)
162  std::cout << C1 << std::endl;
163  for (i=0; i<Map.NumMyElements(); i++) assert(C2[i]==C0[i]);
164  assert(C2.DefaultColor()==newDefaultColor);
165 
166  assert(numColors==C2.NumColors());
167 
168  for (i=0; i<maxcolor; i++) {
169  int curNumElementsWithColor = C2.NumElementsWithColor(i);
170  assert(colorCount[i]==curNumElementsWithColor);
171  int * curColorLIDList = C2.ColorLIDList(i);
172  if (curNumElementsWithColor==0) {
173  assert(curColorLIDList==0);
174  }
175  else
176  for (int j=0; j<curNumElementsWithColor; j++) assert(curColorLIDList[j]==colorLIDs[i][j]);
177  }
178  int curColor = 1;
179  Epetra_Map * Map1 = C2.GenerateMap(curColor);
180  Epetra_BlockMap * Map2 = C2.GenerateBlockMap(curColor);
181 
182  assert(Map1->NumMyElements()==colorCount[curColor]);
183  assert(Map2->NumMyElements()==colorCount[curColor]);
184 
185  for (i=0; i<Map1->NumMyElements(); i++) {
186  assert(Map1->GID(i)==Map.GID(colorLIDs[curColor][i]));
187  assert(Map2->GID(i)==Map.GID(colorLIDs[curColor][i]));
188  assert(Map2->ElementSize(i)==Map.ElementSize(colorLIDs[curColor][i]));
189  }
190 
191  // Now test data redistribution capabilities
192 
193 
194  Epetra_Map ContiguousMap(-1, Map.NumMyElements(), Map.IndexBase(), Comm);
195  // This vector contains the element sizes for the original map.
196  Epetra_IntVector elementSizes(Copy, ContiguousMap, Map.ElementSizeList());
197  Epetra_IntVector elementIDs(Copy, ContiguousMap, Map.MyGlobalElements());
198  Epetra_IntVector elementColorValues(Copy, ContiguousMap, C2.ElementColors());
199 
200 
201  int NumMyElements0 = 0;
202  if (Comm.MyPID()==0) NumMyElements0 = Map.NumGlobalElements();
203  Epetra_Map CMap0(-1, NumMyElements0, Map.IndexBase(), Comm);
204  Epetra_Import importer(CMap0, ContiguousMap);
205  Epetra_IntVector elementSizes0(CMap0);
206  Epetra_IntVector elementIDs0(CMap0);
207  Epetra_IntVector elementColorValues0(CMap0);
208  elementSizes0.Import(elementSizes, importer, Insert);
209  elementIDs0.Import(elementIDs, importer, Insert);
210  elementColorValues0.Import(elementColorValues, importer, Insert);
211 
212  Epetra_BlockMap MapOnPE0(-1,NumMyElements0, elementIDs0.Values(),
213  elementSizes0.Values(), Map.IndexBase(), Comm);
214 
215  Epetra_Import importer1(MapOnPE0, Map);
216  Epetra_MapColoring ColoringOnPE0(MapOnPE0);
217  ColoringOnPE0.Import(C2, importer1, Insert);
218 
219  for (i=0; i<MapOnPE0.NumMyElements(); i++)
220  assert(ColoringOnPE0[i]==elementColorValues0[i]);
221 
222  if (veryVerbose)
223  std::cout << "Same Map Coloring on PE 0 only" << std::endl;
224  if (veryVerbose1)
225  std::cout << ColoringOnPE0 << std::endl;
226  Epetra_MapColoring C3(Map);
227  C3.Export(ColoringOnPE0, importer1, Insert);
228  for (i=0; i<Map.NumMyElements(); i++) assert(C3[i]==C2[i]);
229  if (veryVerbose)
230  std::cout << "Same Map Coloring after Import/Export exercise" << std::endl;
231  if (veryVerbose1)
232  std::cout << ColoringOnPE0 << std::endl;
233 
234 
235  if (verbose) std::cout << "Checked OK\n\n" << std::endl;
236 
237  if (verbose1) {
238  if (verbose) std::cout << "Test ostream << operator" << std::endl << std::flush;
239  std::cout << C0 << std::endl;
240  }
241 
242 
243  delete [] elementColors;
244  for (i=0; i<maxcolor; i++) if (colorLIDs[i]!=0) delete [] colorLIDs[i];
245  delete [] colorLIDs;
246  delete [] colorCount;
247 
248  delete Map1;
249  delete Map2;
250 
251 
252 #ifdef EPETRA_MPI
253  MPI_Finalize();
254 #endif
255 
256  return returnierr;
257 }
258 
int NumGlobalElements() const
Number of elements across all processors.
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
Epetra_MapColoring: A class for coloring Epetra_Map and Epetra_BlockMap objects.
int ElementSize() const
Returns the size of elements in the map; only valid if map has constant element size.
Epetra_BlockMap * GenerateBlockMap(int Color) const
Generates an Epetra_BlockMap of the GIDs associated with the specified color.
Epetra_IntVector: A class for constructing and using dense integer vectors on a parallel computer...
int MyGlobalElements(int *MyGlobalElementList) const
Puts list of global elements on this processor into the user-provided array.
static void SetTracebackMode(int TracebackModeValue)
Set the value of the Epetra_Object error traceback report mode.
int * ElementSizeList() const
List of the element sizes corresponding to the array MyGlobalElements().
int MyPID() const
Return my process ID.
Epetra_MpiComm: The Epetra MPI Communication Class.
std::string Epetra_Version()
Epetra_Import: This class builds an import object for efficient importing of off-processor elements...
Definition: Epetra_Import.h:63
int * Values() const
Returns a pointer to an array containing the values of this vector.
int NumElementsWithColor(int Color) const
Returns number of map elements on calling processor having specified Color.
int * ColorLIDList(int Color) const
Returns pointer to array of Map LIDs associated with the specified color.
int DefaultColor() const
Returns default color.
int IndexBase() const
Index base for this map.
int NumMyElements() const
Number of elements on the calling processor.
int Export(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
Exports an Epetra_DistObject using the Epetra_Import object.
int * ElementColors() const
Returns pointer to array of the colors associated with the LIDs on the calling processor.
int GID(int LID) const
Returns global ID of local ID, return IndexBase-1 if not found on this processor. ...
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
int NumProc() const
Returns total number of processes (always returns 1 for SerialComm).
int NumColors() const
Returns number of colors on the calling processor.
Epetra_SerialComm: The Epetra Serial Communication Class.
Epetra_Map * GenerateMap(int Color) const
Generates an Epetra_Map of the GIDs associated with the specified color.
int main(int argc, char *argv[])
int Import(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
Imports an Epetra_DistObject using the Epetra_Import object.