Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/Map_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_Map (and Epetra_LocalMap) Test routine
44 #include "Epetra_ConfigDefs.h"
45 #include "Epetra_Time.h"
46 #include "Epetra_Map.h"
47 #include "Epetra_LocalMap.h"
48 #ifdef EPETRA_MPI
49 #include "Epetra_MpiComm.h"
50 #include <mpi.h>
51 #else
52 #include "Epetra_SerialComm.h"
53 #endif
54 #include "checkmap.h"
55 #include "../epetra_test_err.h"
56 #include "Epetra_Version.h"
57 
58 int checkMapDataClass(Epetra_Comm& Comm, int verbose);
59 int checkLocalMapDataClass(Epetra_Comm& Comm, int verbose);
60 
61 int main(int argc, char *argv[]) {
62 
63  int ierr=0, returnierr=0;
64 
65 #ifdef EPETRA_MPI
66  MPI_Init(&argc,&argv);
67  Epetra_MpiComm Comm(MPI_COMM_WORLD);
68 #else
69  Epetra_SerialComm Comm;
70 #endif
71 
72  bool verbose = false;
73 
74  // Check if we should print results to standard out
75  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;
76 
77 
78  if (!verbose) {
79  Comm.SetTracebackMode(0); // This should shut down any error traceback reporting
80  }
81  int MyPID = Comm.MyPID();
82  int NumProc = Comm.NumProc();
83 
84  if (verbose && MyPID==0)
85  cout << Epetra_Version() << endl << endl;
86 
87  if (verbose) cout << Comm << endl;
88 
89  bool verbose1 = verbose;
90  if (verbose) verbose = (MyPID==0);
91 
92  int NumMyElements = 10000;
93  int NumMyElements1 = NumMyElements; // Used for local map
94  long long NumGlobalElements = ((long long)NumMyElements)*NumProc+EPETRA_MIN(NumProc,3);
95  if (MyPID < 3) NumMyElements++;
96  long long IndexBase = 0;
97  bool DistributedGlobal = (NumGlobalElements>NumMyElements);
98 
99  Epetra_Map* Map;
100 
101  // Test exceptions
102 
103  if (verbose)
104  cout << "*******************************************************************************************" << endl
105  << " Testing Exceptions (Expect error messages if EPETRA_NO_ERROR_REPORTS is not defined" << endl
106  << "*******************************************************************************************" << endl
107  << endl << endl;
108 
109  try {
110  if (verbose) cout << "Checking Epetra_Map(-2, IndexBase, Comm)" << endl;
111  Map = new Epetra_Map((long long)-2, IndexBase, Comm);
112  }
113  catch (int Error) {
114  if (Error!=-1) {
115  if (Error!=0) {
116  EPETRA_TEST_ERR(Error,returnierr);
117  if (verbose) cout << "Error code should be -1" << endl;
118  }
119  else {
120  cout << "Error code = " << Error << "Should be -1" << endl;
121  returnierr+=1;
122  }
123  }
124  else if (verbose) cout << "Checked OK\n\n" << endl;
125  }
126 
127  try {
128  if (verbose) cout << "Checking Epetra_Map(2, 3, IndexBase, Comm)" << endl;
129  Map = new Epetra_Map((long long)2, 3, IndexBase, Comm);
130  }
131  catch (int Error) {
132  if (Error!=-4) {
133  if (Error!=0) {
134  EPETRA_TEST_ERR(Error,returnierr);
135  if (verbose) cout << "Error code should be -4" << endl;
136  }
137  else {
138  cout << "Error code = " << Error << "Should be -4" << endl;
139  returnierr+=1;
140  }
141  }
142  else if (verbose) cout << "Checked OK\n\n" << endl;
143  }
144 
145  if (verbose) cerr << flush;
146  if (verbose) cout << flush;
147  Comm.Barrier();
148  if (verbose)
149  cout << endl << endl
150  << "*******************************************************************************************" << endl
151  << " Testing valid constructor now......................................................" << endl
152  << "*******************************************************************************************" << endl
153  << endl << endl;
154 
155  // Test Epetra-defined uniform linear distribution constructor
156  Map = new Epetra_Map(NumGlobalElements, IndexBase, Comm);
157  if (verbose) cout << "Checking Epetra_Map(NumGlobalElements, IndexBase, Comm)" << endl;
158  ierr = checkmap(*Map, NumGlobalElements, NumMyElements, 0,
159  IndexBase, Comm, DistributedGlobal);
160 
161  EPETRA_TEST_ERR(ierr,returnierr);
162  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
163 
164  delete Map;
165 
166  // Test User-defined linear distribution constructor
167  Map = new Epetra_Map(NumGlobalElements, NumMyElements, IndexBase, Comm);
168 
169  if (verbose) cout << "Checking Epetra_Map(NumGlobalElements, NumMyElements, IndexBase, Comm)" << endl;
170  ierr = checkmap(*Map, NumGlobalElements, NumMyElements, 0,
171  IndexBase, Comm, DistributedGlobal);
172 
173  EPETRA_TEST_ERR(ierr,returnierr);
174  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
175  delete Map;
176 
177  // Test User-defined arbitrary distribution constructor
178  // Generate Global Element List. Do in reverse for fun!
179 
180  long long * MyGlobalElements = new long long[NumMyElements];
181  long long MaxMyGID = (Comm.MyPID()+1)*NumMyElements-1+IndexBase;
182  if (Comm.MyPID()>2) MaxMyGID+=3;
183  for (int i = 0; i<NumMyElements; i++) MyGlobalElements[i] = MaxMyGID-i;
184 
185  Map = new Epetra_Map(NumGlobalElements, NumMyElements, MyGlobalElements,
186  IndexBase, Comm);
187  if (verbose) cout << "Checking Epetra_Map(NumGlobalElements, NumMyElements, MyGlobalElements, IndexBase, Comm)" << endl;
188  ierr = checkmap(*Map, NumGlobalElements, NumMyElements, MyGlobalElements,
189  IndexBase, Comm, DistributedGlobal);
190 
191  EPETRA_TEST_ERR(ierr,returnierr);
192  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
193  // Test Copy constructor
194  Epetra_Map* Map1 = new Epetra_Map(*Map);
195 
196  // Test SameAs() method
197  bool same = Map1->SameAs(*Map);
198  EPETRA_TEST_ERR(!(same==true),ierr);// should return true since Map1 is a copy of Map
199 
200  Epetra_BlockMap* Map2 = new Epetra_Map(NumGlobalElements, NumMyElements, MyGlobalElements, IndexBase, Comm);
201  same = Map2->SameAs(*Map);
202  EPETRA_TEST_ERR(!(same==true),ierr); // Map and Map2 were created with the same sets of parameters
203  delete Map2;
204 
205  // now test SameAs() on a map that is different
206 
207  Map2 = new Epetra_Map(NumGlobalElements, NumMyElements, MyGlobalElements, IndexBase-1, Comm);
208  same = Map2->SameAs(*Map);
209  EPETRA_TEST_ERR(!(same==false),ierr); // IndexBases are different
210  delete Map2;
211 
212  // Back to testing copy constructor
213  if (verbose) cout << "Checking Epetra_Map(*Map)" << endl;
214  ierr = checkmap(*Map1, NumGlobalElements, NumMyElements, MyGlobalElements,
215  IndexBase, Comm, DistributedGlobal);
216 
217  EPETRA_TEST_ERR(ierr,returnierr);
218  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
219  Epetra_Map* SmallMap = 0;
220  if (verbose1) {
221  // Build a small map for test cout. Use 10 elements from current map
222  long long* MyEls = Map->MyGlobalElements64();
223  long long IndBase = Map->IndexBase64();
224  int MyLen = EPETRA_MIN(10+Comm.MyPID(),Map->NumMyElements());
225  SmallMap = new Epetra_Map((long long)-1, MyLen, MyEls, IndBase, Comm);
226  }
227 
228  delete [] MyGlobalElements;
229  delete Map;
230  delete Map1;
231 
232  // Test reference-counting in Epetra_Map
233  if (verbose) cout << "Checking Epetra_Map reference counting" << endl;
234  ierr = checkMapDataClass(Comm, verbose);
235  EPETRA_TEST_ERR(ierr,returnierr);
236  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
237 
238  // Test LocalMap constructor
239  Epetra_LocalMap* LocalMap = new Epetra_LocalMap((long long)NumMyElements1, IndexBase, Comm);
240  if (verbose) cout << "Checking Epetra_LocalMap(NumMyElements1, IndexBase, Comm)" << endl;
241  ierr = checkmap(*LocalMap, NumMyElements1, NumMyElements1, 0, IndexBase, Comm, false);
242 
243  EPETRA_TEST_ERR(ierr,returnierr);
244  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
245  // Test Copy constructor
246  Epetra_LocalMap* LocalMap1 = new Epetra_LocalMap(*LocalMap);
247  if (verbose) cout << "Checking Epetra_LocalMap(*LocalMap)" << endl;
248  ierr = checkmap(*LocalMap1, NumMyElements1, NumMyElements1, 0, IndexBase, Comm, false);
249 
250  EPETRA_TEST_ERR(ierr,returnierr);
251  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
252  delete LocalMap1;
253  delete LocalMap;
254 
255  // Test reference-counting in Epetra_LocalMap
256  if (verbose) cout << "Checking Epetra_LocalMap reference counting" << endl;
257  ierr = checkLocalMapDataClass(Comm, verbose);
258  EPETRA_TEST_ERR(ierr,returnierr);
259  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
260 
261  // Test output
262  if (verbose1) {
263  if (verbose) cout << "Test ostream << operator" << endl << flush;
264  cout << *SmallMap;
265  delete SmallMap;
266  }
267 
268 #ifdef EPETRA_MPI
269  MPI_Finalize();
270 #endif
271 
272  return returnierr;
273 }
274 
275 int checkMapDataClass(Epetra_Comm& Comm, int verbose) {
276  int returnierr = 0;
277  long long NumGlobalElements = 1000;
278  long long IndexBase = 0;
279 
280  Epetra_Map m1(NumGlobalElements, IndexBase, Comm);
281  int m1count = m1.ReferenceCount();
282  const Epetra_BlockMapData* m1addr = m1.DataPtr();
283  EPETRA_TEST_ERR(!(m1count==1),returnierr); // count should be 1
284  if(verbose) cout << "Default constructor. \nm1= " << m1count << " " << m1addr << endl;
285 
286  Epetra_Map* m2 = new Epetra_Map(m1);
287  int m2count = m2->ReferenceCount();
288  const Epetra_BlockMapData* m2addr = m2->DataPtr();
289  int m1countold = m1count;
290  m1count = m1.ReferenceCount();
291  EPETRA_TEST_ERR(!(m2count==m1count && m1count==(m1countold+1)),returnierr); // both counts should be 2
292  EPETRA_TEST_ERR(!(m1addr==m2addr),returnierr); // addresses should be same
293  if(verbose) cout << "Copy constructor. \nm1= " << m1count << " " << m1addr
294  << "\nm2= " << m2count << " " << m2addr << endl;
295 
296  delete m2;
297  m1countold = m1count;
298  m1count = m1.ReferenceCount();
299  EPETRA_TEST_ERR(!(m1count == m1countold-1), returnierr); // count should have decremented (to 1)
300  EPETRA_TEST_ERR(!(m1addr == m1.DataPtr()),returnierr); // m1addr should be unchanged
301  if(verbose) cout << "m2 destroyed. \nm1= " << m1count << " " << m1addr << endl;
302 
303  { // inside of braces to test stack deallocation.
304  if(verbose) cout << "Assignment operator, post construction" << endl;
305  Epetra_Map m3(NumGlobalElements, IndexBase+1, Comm);
306  int m3count = m3.ReferenceCount();
307  const Epetra_BlockMapData* m3addr = m3.DataPtr();
308  EPETRA_TEST_ERR(!(m3count==1),returnierr); // m3count should be 1 initially
309  EPETRA_TEST_ERR(!(m1addr!=m3addr),returnierr); // m1 and m3 should have different ptr addresses
310  if(verbose) cout << "Prior to assignment: \nm1= " << m1count << " " << m1addr
311  << "\nm3= " << m3count << " " << m3addr << endl;
312  m3 = m1;
313  m3count = m3.ReferenceCount();
314  m3addr = m3.DataPtr();
315  m1countold = m1count;
316  m1count = m1.ReferenceCount();
317  EPETRA_TEST_ERR(!(m3count==m1count && m1count==m1countold+1),returnierr); // both counts should be 2
318  EPETRA_TEST_ERR(!(m1addr==m3addr),returnierr); // addresses should be same
319  if(verbose) cout << "After assignment: \nm1= " << m1count << " " << m1addr
320  << "\nm3= " << m3count << " " << m3addr << endl;
321  }
322  m1countold = m1count;
323  m1count = m1.ReferenceCount();
324  EPETRA_TEST_ERR(!(m1count==m1countold-1), returnierr); // count should have decremented (to 1)
325  EPETRA_TEST_ERR(!(m1addr== m1.DataPtr()),returnierr); // m1addr should be unchanged
326  if(verbose) cout << "m3 destroyed. \nm1= " << m1count << " " << m1addr << endl;
327 
328  return(returnierr);
329 }
330 
331 int checkLocalMapDataClass(Epetra_Comm& Comm, int verbose) {
332  int returnierr = 0;
333  long long NumMyElements = 100;
334  long long IndexBase = 0;
335 
336  Epetra_LocalMap m1(NumMyElements, IndexBase, Comm);
337  int m1count = m1.ReferenceCount();
338  const Epetra_BlockMapData* m1addr = m1.DataPtr();
339  EPETRA_TEST_ERR(!(m1count==1),returnierr); // count should be 1
340  if(verbose) cout << "Default constructor. \nm1= " << m1count << " " << m1addr << endl;
341 
342  Epetra_LocalMap* m2 = new Epetra_LocalMap(m1);
343  int m2count = m2->ReferenceCount();
344  const Epetra_BlockMapData* m2addr = m2->DataPtr();
345  int m1countold = m1count;
346  m1count = m1.ReferenceCount();
347  EPETRA_TEST_ERR(!(m2count==m1count && m1count==(m1countold+1)),returnierr); // both counts should be 2
348  EPETRA_TEST_ERR(!(m1addr==m2addr),returnierr); // addresses should be same
349  if(verbose) cout << "Copy constructor. \nm1= " << m1count << " " << m1addr
350  << "\nm2= " << m2count << " " << m2addr << endl;
351 
352  delete m2;
353  m1countold = m1count;
354  m1count = m1.ReferenceCount();
355  EPETRA_TEST_ERR(!(m1count == m1countold-1), returnierr); // count should have decremented (to 1)
356  EPETRA_TEST_ERR(!(m1addr == m1.DataPtr()),returnierr); // m1addr should be unchanged
357  if(verbose) cout << "m2 destroyed. \nm1= " << m1count << " " << m1addr << endl;
358 
359  { // inside of braces to test stack deallocation.
360  if(verbose) cout << "Assignment operator, post construction" << endl;
361  Epetra_LocalMap m3(NumMyElements, IndexBase+1, Comm);
362  int m3count = m3.ReferenceCount();
363  const Epetra_BlockMapData* m3addr = m3.DataPtr();
364  EPETRA_TEST_ERR(!(m3count==1),returnierr); // m3count should be 1 initially
365  EPETRA_TEST_ERR(!(m1addr!=m3addr),returnierr); // m1 and m3 should have different ptr addresses
366  if(verbose) cout << "Prior to assignment: \nm1= " << m1count << " " << m1addr
367  << "\nm3= " << m3count << " " << m3addr << endl;
368  m3 = m1;
369  m3count = m3.ReferenceCount();
370  m3addr = m3.DataPtr(); // cast int* to int
371  m1countold = m1count;
372  m1count = m1.ReferenceCount();
373  EPETRA_TEST_ERR(!(m3count==m1count && m1count==m1countold+1),returnierr); // both counts should be 2
374  EPETRA_TEST_ERR(!(m1addr==m3addr),returnierr); // addresses should be same
375  if(verbose) cout << "After assignment: \nm1= " << m1count << " " << m1addr
376  << "\nm3= " << m3count << " " << m3addr << endl;
377  }
378  m1countold = m1count;
379  m1count = m1.ReferenceCount();
380  EPETRA_TEST_ERR(!(m1count==m1countold-1), returnierr); // count should have decremented (to 1)
381  EPETRA_TEST_ERR(!(m1addr== m1.DataPtr()),returnierr); // m1addr should be unchanged
382  if(verbose) cout << "m3 destroyed. \nm1= " << m1count << " " << m1addr << endl;
383 
384  return(returnierr);
385 }
Epetra_BlockMapData: The Epetra BlockMap Data Class.
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
bool SameAs(const Epetra_BlockMap &Map) const
Returns true if this and Map are identical maps.
#define EPETRA_TEST_ERR(a, b)
long long IndexBase64() const
int ReferenceCount() const
Returns the reference count of BlockMapData.
static void SetTracebackMode(int TracebackModeValue)
Set the value of the Epetra_Object error traceback report mode.
#define EPETRA_MIN(x, y)
int MyPID() const
Return my process ID.
int checkMapDataClass(Epetra_Comm &Comm, int verbose)
Epetra_MpiComm: The Epetra MPI Communication Class.
int checkLocalMapDataClass(Epetra_Comm &Comm, int verbose)
std::string Epetra_Version()
int NumMyElements() const
Number of elements on the calling processor.
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:73
int checkmap(Epetra_BlockMap &Map, int NumGlobalElements, int NumMyElements, int *MyGlobalElements, int ElementSize, int *ElementSizeList, int NumGlobalPoints, int NumMyPoints, int IndexBase, Epetra_Comm &Comm, bool DistributedGlobal, bool IsOneToOne)
int ReferenceCount() const
Get reference count.
Definition: Epetra_Data.cpp:71
Epetra_BlockMap: A class for partitioning block element vectors and matrices.
int NumProc() const
Returns total number of processes (always returns 1 for SerialComm).
Epetra_SerialComm: The Epetra Serial Communication Class.
void Barrier() const
Epetra_SerialComm Barrier function.
const Epetra_BlockMapData * DataPtr() const
Returns a pointer to the BlockMapData instance this BlockMap uses.
int main(int argc, char *argv[])
long long * MyGlobalElements64() const
Epetra_LocalMap: A class for replicating vectors and matrices across multiple processors.