Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/BlockMap/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_Time.h"
46 #include "Epetra_BlockMap.h"
47 #ifdef EPETRA_MPI
48 #include "Epetra_MpiComm.h"
49 #include <mpi.h>
50 #endif
51 #include "Epetra_SerialComm.h"
52 #include "checkmap.h"
53 #include "../epetra_test_err.h"
54 #include "Epetra_Version.h"
55 
56 int main(int argc, char *argv[]) {
57  bool verbose = false;
58  // Check if we should print results to standard out
59  if (argc > 1)
60  if ((argv[1][0] == '-') && (argv[1][1] == 'v'))
61  verbose = true;
62 
63  int i;
64  int ierr = 0;
65  int returnierr = 0;
66 
67 #ifdef EPETRA_MPI
68 
69  // Initialize MPI
70  MPI_Init(&argc,&argv);
71  Epetra_MpiComm Comm(MPI_COMM_WORLD);
72 #else
73  Epetra_SerialComm Comm;
74 #endif
75 
76  if (!verbose) {
77  Comm.SetTracebackMode(0); // This should shut down any error traceback reporting
78  }
79  int MyPID = Comm.MyPID();
80  int NumProc = Comm.NumProc();
81 
82  int verbose_int = verbose ? 1 : 0;
83  Comm.Broadcast(&verbose_int, 1, 0);
84  verbose = verbose_int==1 ? true : false;
85 
86  if (verbose && MyPID==0)
87  cout << Epetra_Version() << endl << endl;
88 
89  if (verbose) cout << Comm << endl << flush;
90  Comm.Barrier();
91  bool verbose1 = verbose;
92  if (verbose) verbose = (MyPID==0);
93 
94  int NumMyElements = 10000;
95  int NumGlobalElements = NumMyElements*NumProc+EPETRA_MIN(NumProc,3);
96  if (MyPID < 3) NumMyElements++;
97  int IndexBase = 0;
98  int ElementSize = 7;
99  bool DistributedGlobal = (NumGlobalElements>NumMyElements);
100  bool IsOneToOne = true;
101  Epetra_BlockMap * Map;
102 
103  // Test exceptions
104 
105  if (verbose)
106  cout << "*******************************************************************************************" << endl
107  << " Testing Exceptions (Expect error messages if EPETRA_NO_ERROR_REPORTS is not defined" << endl
108  << "*******************************************************************************************" << endl
109  << endl << endl;
110 
111  try {
112  if (verbose) cout << "Checking Epetra_BlockMap(-2, ElementSize, IndexBase, Comm)" << endl;
113  Epetra_BlockMap TestMap(-2, ElementSize, IndexBase, Comm);
114  }
115  catch (int Error) {
116  if (Error != -1) {
117  if (Error != 0) {
118  EPETRA_TEST_ERR(Error,returnierr);
119  if (verbose) cout << "Error code should be -1" << endl;
120  }
121  else { // Error == 0
122  cout << "Error code = " << Error << "Should be -1" << endl;
123  returnierr += 1;
124  }
125  }
126  else if (verbose) cout << "Checked OK\n\n" << endl;
127  }
128 
129  try {
130  if (verbose) cout << "Checking Epetra_BlockMap(2, 3, ElementSize, IndexBase, Comm)" << endl;
131  Epetra_BlockMap TestMap(2, 3, ElementSize, IndexBase, Comm);
132  }
133  catch (int Error) {
134  if (Error != -4) {
135  if (Error != 0) {
136  EPETRA_TEST_ERR(Error,returnierr);
137  if (verbose) cout << "Error code should be -4" << endl;
138  }
139  else { // Error == 0
140  cout << "Error code = " << Error << "Should be -4" << endl;
141  returnierr += 1;
142  }
143  }
144  else if (verbose) cout << "Checked OK\n\n" << endl;
145  }
146 
147  if (verbose) cerr << flush;
148  if (verbose) cout << flush;
149  Comm.Barrier();
150  if (verbose)
151  cout << endl << endl
152  << "*******************************************************************************************" << endl
153  << " Testing valid constructor now......................................................" << endl
154  << "*******************************************************************************************" << endl
155  << endl << endl;
156  // Test Epetra-defined uniform linear distribution constructor
157  Map = new Epetra_BlockMap(NumGlobalElements, ElementSize, IndexBase, Comm);
158  if (verbose) cout << "Checking Epetra_BlockMap(NumGlobalElements, ElementSize, IndexBase, Comm)" << endl;
159  ierr = checkmap(*Map, NumGlobalElements, NumMyElements, 0, ElementSize, 0,
160  NumGlobalElements*ElementSize, NumMyElements*ElementSize,
161  IndexBase, Comm, DistributedGlobal,IsOneToOne);
162 
163  EPETRA_TEST_ERR(ierr,returnierr);
164  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
165 
166  delete Map;
167 
168  // Test User-defined linear distribution constructor
169  Map = new Epetra_BlockMap(NumGlobalElements, NumMyElements, ElementSize, IndexBase, Comm);
170 
171  if (verbose) cout << "Checking Epetra_BlockMap(NumGlobalElements, NumMyElements, ElementSize, IndexBase, Comm)" << endl;
172  ierr = checkmap(*Map, NumGlobalElements, NumMyElements, 0, ElementSize, 0,
173  NumGlobalElements*ElementSize, NumMyElements*ElementSize,
174  IndexBase, Comm, DistributedGlobal,IsOneToOne);
175 
176  EPETRA_TEST_ERR(ierr,returnierr);
177  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
178 
179  delete Map;
180 
181  // Test User-defined arbitrary distribution constructor and fill MyGlobalElements
182  // such that the map is not one-to-one.
183  int NumMyElems = 5;
184  int NumGlobalElems = (Comm.NumProc()+1)*NumMyElems;
185  int myFirstElem = Comm.MyPID()*NumMyElems;
186  if (Comm.MyPID() == 0) NumMyElems *= 2;
187 
188  int* myElems = new int[NumMyElems];
189  for(int ii=0; ii<NumMyElems; ++ii) {
190  myElems[ii] = myFirstElem + ii;
191  }
192 
193  Map = new Epetra_BlockMap(NumGlobalElems, NumMyElems, myElems, 1, 0, Comm);
194 
195  if (verbose) cout << "Checking non-oneToOne Epetra_BlockMap(...)"<<endl;
196  ierr = Map->IsOneToOne() == false ? 0 : -1;
197 
198  //this Map is 1-to-1 if we're running on 1 processor, otherwise it
199  //should not be 1-to-1.
200  if (Comm.NumProc() > 1) {
201  EPETRA_TEST_ERR(ierr,returnierr);
202  }
203  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
204 
205  delete [] myElems;
206  delete Map;
207 
208  // Test User-defined arbitrary distribution constructor
209  // Generate Global Element List. Do in reverse for fun!
210 
211  int * MyGlobalElements = new int[NumMyElements];
212  int MaxMyGID = (Comm.MyPID()+1)*NumMyElements-1+IndexBase;
213  if (Comm.MyPID()>2)
214  MaxMyGID+=3;
215  for (i = 0; i<NumMyElements; i++)
216  MyGlobalElements[i] = MaxMyGID-i;
217 
218  Map = new Epetra_BlockMap(NumGlobalElements, NumMyElements, MyGlobalElements, ElementSize,
219  IndexBase, Comm);
220 
221  if (verbose) cout << "Checking Epetra_BlockMap(NumGlobalElements, NumMyElements, MyGlobalElements, ElementSize, IndexBase, Comm)" << endl;
222  ierr = checkmap(*Map, NumGlobalElements, NumMyElements, MyGlobalElements, ElementSize, 0,
223  NumGlobalElements*ElementSize, NumMyElements*ElementSize,
224  IndexBase, Comm, DistributedGlobal,IsOneToOne);
225 
226  EPETRA_TEST_ERR(ierr,returnierr);
227  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
228 
229  Epetra_BlockMap * Map3 = new Epetra_BlockMap(*Map);// A map to test the SameAs method later
230 
231  delete Map;
232 
233  int * ElementSizeList = new int[NumMyElements];
234  int NumMyEquations = 0;
235  int NumGlobalEquations = 0;
236  for (i = 0; i<NumMyElements; i++) {
237  ElementSizeList[i] = i%6 + 2; // elementsizes go from 2 to 7
238  NumMyEquations += ElementSizeList[i];
239  }
240  ElementSize = 7; // Set to maximum for use in checkmap
241  NumGlobalEquations = Comm.NumProc()*NumMyEquations;
242 
243  // Adjust NumGlobalEquations based on processor ID
244  if (Comm.NumProc() > 3) {
245  if (Comm.MyPID()>2)
246  NumGlobalEquations += 3*((NumMyElements)%6+2);
247  else
248  NumGlobalEquations -= (Comm.NumProc()-3)*((NumMyElements-1)%6+2);
249  }
250  Map = new Epetra_BlockMap(NumGlobalElements, NumMyElements, MyGlobalElements, ElementSizeList,
251  IndexBase, Comm);
252  if (verbose) cout << "Checking Epetra_BlockMap(NumGlobalElements, NumMyElements, MyGlobalElements, ElementSizeList, IndexBase, Comm)" << endl;
253  ierr = checkmap(*Map, NumGlobalElements, NumMyElements, MyGlobalElements, ElementSize, ElementSizeList,
254  NumGlobalEquations, NumMyEquations,
255  IndexBase, Comm, DistributedGlobal,IsOneToOne);
256 
257  EPETRA_TEST_ERR(ierr,returnierr);
258  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
259 
260  // Test Copy constructor
261  Epetra_BlockMap * Map1 = new Epetra_BlockMap(*Map);
262 
263  // Test SameAs() method
264  bool same = Map1->SameAs(*Map);
265  EPETRA_TEST_ERR(!(same==true),returnierr);// should return true since Map1 is a copy of Map
266 
267  Epetra_BlockMap * Map2 = new Epetra_BlockMap(NumGlobalElements,NumMyElements,MyGlobalElements,ElementSizeList,IndexBase,Comm);
268  same = Map2->SameAs(*Map);
269  EPETRA_TEST_ERR(!(same==true),returnierr); // Map and Map2 were created with the same sets of parameters
270  delete Map2;
271 
272  // now test SameAs() on some maps that are different
273 
274  Map2 = new Epetra_BlockMap(NumGlobalElements,NumMyElements,MyGlobalElements,ElementSizeList,IndexBase-1,Comm);
275  same = Map2->SameAs(*Map);
276  EPETRA_TEST_ERR(!(same==false),returnierr); // IndexBases are different
277  delete Map2;
278 
279  int *ElementSizeList1 = new int[NumMyElements];
280  for (i=0; i<NumMyElements; i++)
281  ElementSizeList1[i] = i%5 + 2; // element sizes go from 2 to 6
282  Map2 = new Epetra_BlockMap(NumGlobalElements,NumMyElements,MyGlobalElements,ElementSizeList1,IndexBase,Comm);
283  same = Map2->SameAs(*Map);
284  EPETRA_TEST_ERR(!(same==false),returnierr); // ElementSizes are different
285  delete [] ElementSizeList1;
286  delete Map2;
287 
288  same = Map3->SameAs(*Map);
289  EPETRA_TEST_ERR(!(same==false),returnierr); // Map3 saved from an earlier test
290  delete Map3;
291 
292  // Back to testing copy constructor
293  if (verbose) cout << "Checking Epetra_BlockMap(*Map)" << endl;
294  ierr = checkmap(*Map1, NumGlobalElements, NumMyElements, MyGlobalElements, ElementSize, ElementSizeList,
295  NumGlobalEquations, NumMyEquations,
296  IndexBase, Comm, DistributedGlobal,IsOneToOne);
297 
298  EPETRA_TEST_ERR(ierr,returnierr);
299  if (verbose && ierr==0) cout << "Checked OK\n\n" <<endl;
300 
301  if (verbose1) {
302  if (verbose) cout << "Test ostream << operator" << endl << flush;
303  }
304  // Build a small map for test cout. Use 10 elements from current map
305  int * MyEls = Map->MyGlobalElements();
306  int * MySz = Map->ElementSizeList();
307  int IndBase = Map->IndexBase();
308  int MyLen = EPETRA_MIN(10+Comm.MyPID(),Map->NumMyElements());
309  Epetra_BlockMap * SmallMap = new Epetra_BlockMap(-1, MyLen, MyEls, MySz, IndBase, Comm);
310  if (verbose1) {
311  cout << *SmallMap;
312  }
313  delete SmallMap;
314 
315  delete Map;
316  delete Map1;
317 
318 
319  //create a map where proc 1 has no local elements, then check to make sure that
320  //if NumMyElements == 0, then MaxMyGID < MinMyGID.
321 
322  if (MyPID == 1) {
323  Map1 = new Epetra_BlockMap(-1, 0, (int*)0, (int*)0, IndexBase, Comm);
324  }
325  else {
326  Map1 = new Epetra_BlockMap(-1, NumMyElements, MyGlobalElements, ElementSizeList, IndexBase, Comm);
327  }
328 
329  int numMyElems = Map1->NumMyElements();
330  if (MyPID == 1) {
331  EPETRA_TEST_ERR(!(numMyElems == 0), returnierr);
332  int maxgid = Map1->MaxMyGID();
333  int mingid = Map1->MinMyGID();
334  EPETRA_TEST_ERR( !(maxgid<mingid), returnierr);
335  }
336 
337  delete[] ElementSizeList;
338  delete[] MyGlobalElements;
339  delete Map1;
340 
341  // test reference counting
342  ierr = 0;
343 
344  if (verbose)
345  cout << endl << endl
346  << "*******************************************************************************************" << endl
347  << " Testing reference counting now....................................................." << endl
348  << "*******************************************************************************************" << endl << endl;
349 
350  Epetra_BlockMap b1(NumGlobalElements, NumMyElements, ElementSize, IndexBase, Comm);
351  int b1count = b1.ReferenceCount();
352  const Epetra_BlockMapData* b1addr = b1.DataPtr();
353  EPETRA_TEST_ERR(!(b1count==1),ierr); // count should be 1
354  if(verbose) cout << "Default constructor. \nb1= " << b1count << " " << b1addr << endl;
355 
356  Epetra_BlockMap* b2 = new Epetra_BlockMap(b1);
357  int b2count = b2->ReferenceCount();
358  const Epetra_BlockMapData* b2addr = b2->DataPtr();
359  int b1countold = b1count;
360  b1count = b1.ReferenceCount();
361  EPETRA_TEST_ERR(!(b2count==b1count && b1count==(b1countold+1)),ierr); // both counts should be 2
362  EPETRA_TEST_ERR(!(b1addr==b2addr),ierr); // addresses should be same
363  if(verbose) cout << "Copy constructor. \nb1= " << b1count << " " << b1addr << "\nb2= " << b2count << " " << b2addr << endl;
364 
365  delete b2;
366  b1countold = b1count;
367  b1count = b1.ReferenceCount();
368  EPETRA_TEST_ERR(!(b1count==b1countold-1), ierr); // count should have decremented (to 1)
369  EPETRA_TEST_ERR(!(b1addr==b1.DataPtr()), ierr); // b1addr should be unchanged
370  if(verbose) cout << "b2 destroyed. \nb1= " << b1count << " " << b1addr << endl;
371 
372  { // inside of braces to test stack deallocation.
373  if(verbose) cout << "Assignment operator, post construction" << endl;
374  Epetra_BlockMap b3(NumGlobalElements, NumMyElements, ElementSize, IndexBase-1, Comm);
375  int b3count = b3.ReferenceCount();
376  const Epetra_BlockMapData* b3addr = b3.DataPtr();
377  EPETRA_TEST_ERR(!(b3count==1),ierr); // b3count should be 1 initially
378  EPETRA_TEST_ERR(!(b1addr!=b3addr),ierr); // b1 and b3 should have different ptr addresses
379  if(verbose) cout << "Prior to assignment: \nb1= " << b1count << " " << b1addr << "\nb3= " << b3count << " " << b3addr << endl;
380  b3 = b1;
381  b3count = b3.ReferenceCount();
382  b3addr = b3.DataPtr();
383  b1countold = b1count;
384  b1count = b1.ReferenceCount();
385  EPETRA_TEST_ERR(!(b3count==b1count && b1count==b1countold+1),ierr); // both counts should be 2
386  EPETRA_TEST_ERR(!(b1addr==b3addr),ierr); // addresses should be same
387  if(verbose) cout << "After assignment: \nb1= " << b1count << " " << b1addr << "\nb3= " << b3count << " " << b3addr << endl;
388  }
389  b1countold = b1count;
390  b1count = b1.ReferenceCount();
391  EPETRA_TEST_ERR(!(b1count==b1countold-1), ierr); // count should have decremented (to 1)
392  EPETRA_TEST_ERR(!(b1addr==b1.DataPtr()), ierr); // b1addr should be unchanged
393  if (verbose) cout << "b3 destroyed. \nb1= " << b1count << " " << b1addr << endl;
394 
395  EPETRA_TEST_ERR(ierr,returnierr);
396  if (verbose && (ierr == 0)) cout << "Checked OK\n\n" <<endl;
397  // done with reference counting testing
398 
399  // test subcommunicators
400  ierr=0;
401  if (verbose)
402  cout << endl << endl
403  << "*******************************************************************************************" << endl
404  << " Testing subcommunicators now......................................................." << endl
405  << "*******************************************************************************************" << endl << endl;
406 
407  // Create a map where everything is on proc 0
408  if (MyPID != 0) {
409  Map1 = new Epetra_BlockMap(-1, 0, 1, IndexBase, Comm);
410  }
411  else {
412  Map1 = new Epetra_BlockMap(-1, NumMyElements, 1, IndexBase, Comm);
413  }
414 
415  // Remove empty processes ...
416  Map2=0; Map2 = Map1->RemoveEmptyProcesses();
417  if(MyPID==0) {EPETRA_TEST_ERR(Map2 == 0,ierr);}
418  else {EPETRA_TEST_ERR(Map2 != 0,ierr);}
419 
420  // Replace comm
421  const Epetra_Comm * TempComm = Map2 ? &Map2->Comm() : 0;
422  Map3=0; Map3 = Map1->ReplaceCommWithSubset(TempComm);
423  if(MyPID==0) {EPETRA_TEST_ERR(Map3 == 0,ierr);}
424  else {EPETRA_TEST_ERR(Map3 != 0,ierr);}
425 
426  delete Map1; delete Map2; delete Map3;
427 
428  EPETRA_TEST_ERR(ierr,returnierr);
429  if (verbose && (ierr == 0)) cout << "Checked OK\n\n" <<endl;
430  // done with testing subcommunicators
431 
432 #ifdef EPETRA_MPI
433  MPI_Finalize();
434 #endif
435 
436  return returnierr;
437 }
Epetra_BlockMapData: The Epetra BlockMap Data Class.
bool SameAs(const Epetra_BlockMap &Map) const
Returns true if this and Map are identical maps.
#define EPETRA_TEST_ERR(a, b)
int MyGlobalElements(int *MyGlobalElementList) const
Puts list of global elements on this processor into the user-provided array.
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.
int * ElementSizeList() const
List of the element sizes corresponding to the array MyGlobalElements().
Epetra_BlockMap * RemoveEmptyProcesses() const
Return a new BlockMap with processes with zero elements removed.
#define EPETRA_MIN(x, y)
int MyPID() const
Return my process ID.
Epetra_MpiComm: The Epetra MPI Communication Class.
std::string Epetra_Version()
bool IsOneToOne() const
int IndexBase() const
Index base for this map.
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)
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_BlockMap * ReplaceCommWithSubset(const Epetra_Comm *Comm) const
Replace this BlockMap&#39;s communicator with a subset communicator.
int MaxMyGID() const
Returns the maximum global ID owned by this processor.
int MinMyGID() const
Returns the minimum global ID owned by this processor.
const Epetra_Comm & Comm() const
Access function for Epetra_Comm communicator.
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[])
int Broadcast(double *MyVals, int Count, int Root) const
Epetra_SerialComm Broadcast function.