Epetra Package Browser (Single Doxygen Collection)  Development
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  int defaultColor = C0.DefaultColor();
127  for (i=0; i<Map.NumMyElements(); i++) {
128  assert(C0[i]==defaultColor);
129  assert(C0(Map.GID64(i))==defaultColor);
130  if (i%2==0) C0[i] = i%6+5+i%14; // cycle through 5...23 on even elements
131  else C0(Map.GID64(i)) = i%5+1; // cycle through 1...5 on odd elements
132  elementColors[i] = C0[i]; // Record color of ith element for use below
133  colorCount[C0[i]]++; // Count how many of each color for checking below
134  }
135 
136  if (veryVerbose)
137  std::cout << "Original Map Coloring using element-by-element definitions" << std::endl;
138  if (veryVerbose1)
139  std::cout << C0 << std::endl;
140 
141  int numColors = 0;
142  for (i=0; i<maxcolor; i++)
143  if (colorCount[i]>0) {
144  numColors++;
145  colorLIDs[i] = new int[colorCount[i]];
146  }
147  for (i=0; i<maxcolor; i++) colorCount[i] = 0;
148  for (i=0; i<Map.NumMyElements(); i++) colorLIDs[C0[i]][colorCount[C0[i]]++] = i;
149 
150 
151 
152  int newDefaultColor = -1;
153  Epetra_MapColoring C1(Map, elementColors, newDefaultColor);
154  if (veryVerbose)
155  std::cout << "Same Map Coloring using one-time construction" << std::endl;
156  if (veryVerbose1)
157  std::cout << C1 << std::endl;
158  assert(C1.DefaultColor()==newDefaultColor);
159  for (i=0; i<Map.NumMyElements(); i++) assert(C1[i]==C0[i]);
160 
161  Epetra_MapColoring C2(C1);
162  if (veryVerbose)
163  std::cout << "Same Map Coloring using copy constructor" << std::endl;
164  if (veryVerbose1)
165  std::cout << C1 << std::endl;
166  for (i=0; i<Map.NumMyElements(); i++) assert(C2[i]==C0[i]);
167  assert(C2.DefaultColor()==newDefaultColor);
168 
169  assert(numColors==C2.NumColors());
170 
171  for (i=0; i<maxcolor; i++) {
172  int curNumElementsWithColor = C2.NumElementsWithColor(i);
173  assert(colorCount[i]==curNumElementsWithColor);
174  int * curColorLIDList = C2.ColorLIDList(i);
175  if (curNumElementsWithColor==0) {
176  assert(curColorLIDList==0);
177  }
178  else
179  for (int j=0; j<curNumElementsWithColor; j++) assert(curColorLIDList[j]==colorLIDs[i][j]);
180  }
181  int curColor = 1;
182  Epetra_Map * Map1 = C2.GenerateMap(curColor);
183  Epetra_BlockMap * Map2 = C2.GenerateBlockMap(curColor);
184 
185  assert(Map1->NumMyElements()==colorCount[curColor]);
186  assert(Map2->NumMyElements()==colorCount[curColor]);
187 
188  for (i=0; i<Map1->NumMyElements(); i++) {
189  assert(Map1->GID64(i)==Map.GID64(colorLIDs[curColor][i]));
190  assert(Map2->GID64(i)==Map.GID64(colorLIDs[curColor][i]));
191  assert(Map2->ElementSize(i)==Map.ElementSize(colorLIDs[curColor][i]));
192  }
193 
194  // Now test data redistribution capabilities
195 
196 
197  Epetra_Map ContiguousMap(-1LL, Map.NumMyElements(), Map.IndexBase64(), Comm);
198  // This vector contains the element sizes for the original map.
199  Epetra_IntVector elementSizes(Copy, ContiguousMap, Map.ElementSizeList());
200  Epetra_LongLongVector elementIDs(Copy, ContiguousMap, Map.MyGlobalElements64());
201  Epetra_IntVector elementColorValues(Copy, ContiguousMap, C2.ElementColors());
202 
203 
204  long long NumMyElements0 = 0;
205  if (Comm.MyPID()==0) NumMyElements0 = Map.NumGlobalElements64();
206  Epetra_Map CMap0(-1LL, NumMyElements0, Map.IndexBase64(), Comm);
207  Epetra_Import importer(CMap0, ContiguousMap);
208  Epetra_IntVector elementSizes0(CMap0);
209  Epetra_LongLongVector elementIDs0(CMap0);
210  Epetra_IntVector elementColorValues0(CMap0);
211  elementSizes0.Import(elementSizes, importer, Insert);
212  elementIDs0.Import(elementIDs, importer, Insert);
213  elementColorValues0.Import(elementColorValues, importer, Insert);
214 
215  Epetra_BlockMap MapOnPE0(
216  -1LL,NumMyElements0, elementIDs0.Values(),
217  elementSizes0.Values(), Map.IndexBase64(), Comm
218  );
219 
220  Epetra_Import importer1(MapOnPE0, Map);
221  Epetra_MapColoring ColoringOnPE0(MapOnPE0);
222  ColoringOnPE0.Import(C2, importer1, Insert);
223 
224  for (i=0; i<MapOnPE0.NumMyElements(); i++)
225  assert(ColoringOnPE0[i]==elementColorValues0[i]);
226 
227  if (veryVerbose)
228  std::cout << "Same Map Coloring on PE 0 only" << std::endl;
229  if (veryVerbose1)
230  std::cout << ColoringOnPE0 << std::endl;
231  Epetra_MapColoring C3(Map);
232  C3.Export(ColoringOnPE0, importer1, Insert);
233  for (i=0; i<Map.NumMyElements(); i++) assert(C3[i]==C2[i]);
234  if (veryVerbose)
235  std::cout << "Same Map Coloring after Import/Export exercise" << std::endl;
236  if (veryVerbose1)
237  std::cout << ColoringOnPE0 << std::endl;
238 
239 
240  if (verbose) std::cout << "Checked OK\n\n" << std::endl;
241 
242  if (verbose1) {
243  if (verbose) std::cout << "Test ostream << operator" << std::endl << std::flush;
244  std::cout << C0 << std::endl;
245  }
246 
247 
248  delete [] elementColors;
249  for (i=0; i<maxcolor; i++) if (colorLIDs[i]!=0) delete [] colorLIDs[i];
250  delete [] colorLIDs;
251  delete [] colorCount;
252 
253  delete Map1;
254  delete Map2;
255 
256 
257 #ifdef EPETRA_MPI
258  MPI_Finalize();
259 #endif
260 
261  return returnierr;
262 }
263 
Epetra_BlockMap::NumMyElements
int NumMyElements() const
Number of elements on the calling processor.
Definition: Epetra_BlockMap.h:555
Epetra_Version.h
Epetra_BlockMap::ElementSize
int ElementSize() const
Returns the size of elements in the map; only valid if map has constant element size.
Definition: Epetra_BlockMap.h:573
Epetra_MapColoring::GenerateMap
Epetra_Map * GenerateMap(int Color) const
Generates an Epetra_Map of the GIDs associated with the specified color.
Definition: Epetra_MapColoring.cpp:245
Epetra_BlockMap::IndexBase64
long long IndexBase64() const
Definition: Epetra_BlockMap.h:592
Epetra_MapColoring.h
Epetra_IntVector
Epetra_IntVector: A class for constructing and using dense integer vectors on a parallel computer.
Definition: Epetra_IntVector.h:124
Copy
Definition: Epetra_DataAccess.h:55
Epetra_SerialComm::MyPID
int MyPID() const
Return my process ID.
Definition: Epetra_SerialComm.h:432
Epetra_Version
std::string Epetra_Version()
Definition: Epetra_Version.h:47
Epetra_SerialComm::NumProc
int NumProc() const
Returns total number of processes (always returns 1 for SerialComm).
Definition: Epetra_SerialComm.h:435
main
int main(int argc, char *argv[])
Definition: test/MapColoring_LL/cxx_main.cpp:60
Insert
Definition: Epetra_CombineMode.h:68
Epetra_SerialComm.h
Epetra_MapColoring::GenerateBlockMap
Epetra_BlockMap * GenerateBlockMap(int Color) const
Generates an Epetra_BlockMap of the GIDs associated with the specified color.
Definition: Epetra_MapColoring.cpp:295
Epetra_MpiComm.h
Epetra_DistObject::Export
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.
Definition: Epetra_DistObject.cpp:198
Epetra_LongLongVector::Values
long long * Values() const
Returns a pointer to an array containing the values of this vector.
Definition: Epetra_LongLongVector.h:249
Epetra_MapColoring::ElementColors
int * ElementColors() const
Returns pointer to array of the colors associated with the LIDs on the calling processor.
Definition: Epetra_MapColoring.h:222
Epetra_MpiComm
Epetra_MpiComm: The Epetra MPI Communication Class.
Definition: Epetra_MpiComm.h:64
Epetra_SerialComm
Epetra_SerialComm: The Epetra Serial Communication Class.
Definition: Epetra_SerialComm.h:61
Epetra_IntVector.h
Epetra_BlockMap.h
Epetra_BlockMap
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
Definition: Epetra_BlockMap.h:194
Epetra_Import.h
Epetra_LongLongVector
Epetra_LongLongVector: A class for constructing and using dense integer vectors on a parallel compute...
Definition: Epetra_LongLongVector.h:126
Epetra_MapColoring::NumElementsWithColor
int NumElementsWithColor(int Color) const
Returns number of map elements on calling processor having specified Color.
Definition: Epetra_MapColoring.cpp:205
Epetra_IntVector::Values
int * Values() const
Returns a pointer to an array containing the values of this vector.
Definition: Epetra_IntVector.h:247
Epetra_BlockMap::NumGlobalElements64
long long NumGlobalElements64() const
Definition: Epetra_BlockMap.h:552
Epetra_MapColoring
Epetra_MapColoring: A class for coloring Epetra_Map and Epetra_BlockMap objects.
Definition: Epetra_MapColoring.h:105
Epetra_BlockMap::ElementSizeList
int * ElementSizeList() const
List of the element sizes corresponding to the array MyGlobalElements().
Definition: Epetra_BlockMap.cpp:1009
Epetra_MapColoring::NumColors
int NumColors() const
Returns number of colors on the calling processor.
Definition: Epetra_MapColoring.h:193
Epetra_MapColoring::ColorLIDList
int * ColorLIDList(int Color) const
Returns pointer to array of Map LIDs associated with the specified color.
Definition: Epetra_MapColoring.cpp:214
Epetra_Map.h
Epetra_DistObject::Import
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.
Definition: Epetra_DistObject.cpp:113
Epetra_BlockMap::MyGlobalElements64
long long * MyGlobalElements64() const
Definition: Epetra_BlockMap.cpp:908
Epetra_LongLongVector.h
Epetra_Object::SetTracebackMode
static void SetTracebackMode(int TracebackModeValue)
Set the value of the Epetra_Object error traceback report mode.
Definition: Epetra_Object.cpp:77
Epetra_BlockMap::GID64
long long GID64(int LID) const
Definition: Epetra_BlockMap.cpp:1252
Epetra_Map
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
Epetra_MapColoring::DefaultColor
int DefaultColor() const
Returns default color.
Definition: Epetra_MapColoring.h:206
Epetra_Import
Epetra_Import: This class builds an import object for efficient importing of off-processor elements.
Definition: Epetra_Import.h:63
Epetra_Time.h