Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/MapColoring_LL/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_LongLongVector.h"
50 #include "Epetra_Import.h"
51 #ifdef EPETRA_MPI
52 #include "Epetra_MpiComm.h"
53 #include <mpi.h>
54 #endif
55 
56 #include "Epetra_SerialComm.h"
57 #include "Epetra_Time.h"
58 #include "Epetra_Version.h"
59 
60 int main(int argc, char *argv[]) {
61 
62  int i, returnierr=0;
63 
64 #ifdef EPETRA_MPI
65  // Initialize MPI
66  MPI_Init(&argc,&argv);
67  Epetra_MpiComm Comm(MPI_COMM_WORLD);
68 #else
69  Epetra_SerialComm Comm;
70 #endif
71 
72  // Uncomment to debug in parallel int tmp; if (Comm.MyPID()==0) cin >> tmp; Comm.Barrier();
73 
74  bool verbose = false;
75  bool veryVerbose = false;
76 
77  // Check if we should print results to standard out
78  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;
79 
80  // Check if we should print lots of results to standard out
81  if (argc>2) if (argv[2][0]=='-' && argv[2][1]=='v') veryVerbose = true;
82 
83  if (verbose && Comm.MyPID()==0)
84  std::cout << Epetra_Version() << std::endl << std::endl;
85 
86  if (!verbose) Comm.SetTracebackMode(0); // This should shut down any error traceback reporting
87 
88  if (verbose) std::cout << Comm << std::endl << std::flush;
89 
90  bool verbose1 = verbose;
91  if (verbose) verbose = (Comm.MyPID()==0);
92 
93  bool veryVerbose1 = veryVerbose;
94  if (veryVerbose) veryVerbose = (Comm.MyPID()==0);
95 
96  int NumMyElements = 100;
97  if (veryVerbose1) NumMyElements = 10;
98  NumMyElements += Comm.MyPID();
99  int MaxNumMyElements = NumMyElements+Comm.NumProc()-1;
100  int * ElementSizeList = new int[NumMyElements];
101  long long * MyGlobalElements = new long long[NumMyElements];
102 
103  for (i = 0; i<NumMyElements; i++) {
104  MyGlobalElements[i] = (long long)(Comm.MyPID()*MaxNumMyElements+i)*2;
105  ElementSizeList[i] = i%6 + 2; // elementsizes go from 2 to 7
106  }
107 
108  Epetra_BlockMap Map(
109  -1LL, NumMyElements, MyGlobalElements, ElementSizeList,
110  0, Comm
111  );
112 
113  delete [] ElementSizeList;
114  delete [] MyGlobalElements;
115 
116  Epetra_MapColoring C0(Map);
117 
118  int * elementColors = new int[NumMyElements];
119 
120  int maxcolor = 24;
121  int * colorCount = new int[maxcolor];
122  int ** colorLIDs = new int*[maxcolor];
123  for (i=0; i<maxcolor; i++) colorCount[i] = 0;
124  for (i=0; i<maxcolor; i++) colorLIDs[i] = 0;
125 
126  for (i=0; i<Map.NumMyElements(); i++) {
127  assert(C0[i]==C0.DefaultColor());
128  assert(C0(Map.GID64(i))==C0.DefaultColor());
129  if (i%2==0) C0[i] = i%6+5+i%14; // cycle through 5...23 on even elements
130  else C0(Map.GID64(i)) = i%5+1; // cycle through 1...5 on odd elements
131  elementColors[i] = C0[i]; // Record color of ith element for use below
132  colorCount[C0[i]]++; // Count how many of each color for checking below
133  }
134 
135  if (veryVerbose)
136  std::cout << "Original Map Coloring using element-by-element definitions" << std::endl;
137  if (veryVerbose1)
138  std::cout << C0 << std::endl;
139 
140  int numColors = 0;
141  for (i=0; i<maxcolor; i++)
142  if (colorCount[i]>0) {
143  numColors++;
144  colorLIDs[i] = new int[colorCount[i]];
145  }
146  for (i=0; i<maxcolor; i++) colorCount[i] = 0;
147  for (i=0; i<Map.NumMyElements(); i++) colorLIDs[C0[i]][colorCount[C0[i]]++] = i;
148 
149 
150 
151  int newDefaultColor = -1;
152  Epetra_MapColoring C1(Map, elementColors, newDefaultColor);
153  if (veryVerbose)
154  std::cout << "Same Map Coloring using one-time construction" << std::endl;
155  if (veryVerbose1)
156  std::cout << C1 << std::endl;
157  assert(C1.DefaultColor()==newDefaultColor);
158  for (i=0; i<Map.NumMyElements(); i++) assert(C1[i]==C0[i]);
159 
160  Epetra_MapColoring C2(C1);
161  if (veryVerbose)
162  std::cout << "Same Map Coloring using copy constructor" << std::endl;
163  if (veryVerbose1)
164  std::cout << C1 << std::endl;
165  for (i=0; i<Map.NumMyElements(); i++) assert(C2[i]==C0[i]);
166  assert(C2.DefaultColor()==newDefaultColor);
167 
168  assert(numColors==C2.NumColors());
169 
170  for (i=0; i<maxcolor; i++) {
171  int curNumElementsWithColor = C2.NumElementsWithColor(i);
172  assert(colorCount[i]==curNumElementsWithColor);
173  if (curNumElementsWithColor==0) {
174  assert(C2.ColorLIDList(i)==0);
175  }
176  else
177  for (int j=0; j<curNumElementsWithColor; j++) assert(C2.ColorLIDList(i)[j]==colorLIDs[i][j]);
178  }
179  int curColor = 1;
180  Epetra_Map * Map1 = C2.GenerateMap(curColor);
181  Epetra_BlockMap * Map2 = C2.GenerateBlockMap(curColor);
182 
183  assert(Map1->NumMyElements()==colorCount[curColor]);
184  assert(Map2->NumMyElements()==colorCount[curColor]);
185 
186  for (i=0; i<Map1->NumMyElements(); i++) {
187  assert(Map1->GID64(i)==Map.GID64(colorLIDs[curColor][i]));
188  assert(Map2->GID64(i)==Map.GID64(colorLIDs[curColor][i]));
189  assert(Map2->ElementSize(i)==Map.ElementSize(colorLIDs[curColor][i]));
190  }
191 
192  // Now test data redistribution capabilities
193 
194 
195  Epetra_Map ContiguousMap(-1LL, Map.NumMyElements(), Map.IndexBase64(), Comm);
196  // This vector contains the element sizes for the original map.
197  Epetra_IntVector elementSizes(Copy, ContiguousMap, Map.ElementSizeList());
198  Epetra_LongLongVector elementIDs(Copy, ContiguousMap, Map.MyGlobalElements64());
199  Epetra_IntVector elementColorValues(Copy, ContiguousMap, C2.ElementColors());
200 
201 
202  long long NumMyElements0 = 0;
203  if (Comm.MyPID()==0) NumMyElements0 = Map.NumGlobalElements64();
204  Epetra_Map CMap0(-1LL, NumMyElements0, Map.IndexBase64(), Comm);
205  Epetra_Import importer(CMap0, ContiguousMap);
206  Epetra_IntVector elementSizes0(CMap0);
207  Epetra_LongLongVector elementIDs0(CMap0);
208  Epetra_IntVector elementColorValues0(CMap0);
209  elementSizes0.Import(elementSizes, importer, Insert);
210  elementIDs0.Import(elementIDs, importer, Insert);
211  elementColorValues0.Import(elementColorValues, importer, Insert);
212 
213  Epetra_BlockMap MapOnPE0(
214  -1LL,NumMyElements0, elementIDs0.Values(),
215  elementSizes0.Values(), Map.IndexBase64(), Comm
216  );
217 
218  Epetra_Import importer1(MapOnPE0, Map);
219  Epetra_MapColoring ColoringOnPE0(MapOnPE0);
220  ColoringOnPE0.Import(C2, importer1, Insert);
221 
222  for (i=0; i<MapOnPE0.NumMyElements(); i++)
223  assert(ColoringOnPE0[i]==elementColorValues0[i]);
224 
225  if (veryVerbose)
226  std::cout << "Same Map Coloring on PE 0 only" << std::endl;
227  if (veryVerbose1)
228  std::cout << ColoringOnPE0 << std::endl;
229  Epetra_MapColoring C3(Map);
230  C3.Export(ColoringOnPE0, importer1, Insert);
231  for (i=0; i<Map.NumMyElements(); i++) assert(C3[i]==C2[i]);
232  if (veryVerbose)
233  std::cout << "Same Map Coloring after Import/Export exercise" << std::endl;
234  if (veryVerbose1)
235  std::cout << ColoringOnPE0 << std::endl;
236 
237 
238  if (verbose) std::cout << "Checked OK\n\n" << std::endl;
239 
240  if (verbose1) {
241  if (verbose) std::cout << "Test ostream << operator" << std::endl << std::flush;
242  std::cout << C0 << std::endl;
243  }
244 
245 
246  delete [] elementColors;
247  for (i=0; i<maxcolor; i++) if (colorLIDs[i]!=0) delete [] colorLIDs[i];
248  delete [] colorLIDs;
249  delete [] colorCount;
250 
251  delete Map1;
252  delete Map2;
253 
254 
255 #ifdef EPETRA_MPI
256  MPI_Finalize();
257 #endif
258 
259  return returnierr;
260 }
261 
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_LongLongVector: A class for constructing and using dense integer vectors on a parallel compute...
Epetra_BlockMap * GenerateBlockMap(int Color) const
Generates an Epetra_BlockMap of the GIDs associated with the specified color.
long long IndexBase64() const
Epetra_IntVector: A class for constructing and using dense integer vectors on a parallel computer...
long long NumGlobalElements64() const
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 NumMyElements() const
Number of elements on the calling processor.
long long GID64(int LID) const
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.
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[])
long long * MyGlobalElements64() const
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.
long long * Values() const
Returns a pointer to an array containing the values of this vector.