Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/Comm/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 #include <limits.h>
43 
44 // Epetra_Comm Test routine
45 #include "../epetra_test_err.h"
46 #include "Epetra_Time.h"
47 #include "Epetra_Util.h"
48 #include "Epetra_Distributor.h"
49 #include "Epetra_SerialComm.h"
51 #include "Epetra_Version.h"
52 
53 #ifdef EPETRA_MPI
54 #include <mpi.h>
55 #include "Epetra_MpiComm.h"
56 
57 int checkMpiDataClass(bool verbose);
58 #endif
59 
60 int checkSerialDataClass(bool verbose);
61 int checkCommMethods(Epetra_Comm& petracomm,
62  bool verbose, bool verbose1,
63  int& NumProc, int& rank);
64 int checkRankAndSize(Epetra_Comm& petracomm, bool verbose, int rank, int size);
65 void checkBarrier(Epetra_Comm& petracomm, bool verbose, int rank);
66 
68  Epetra_Comm& Comm);
69 
70 int main(int argc, char* argv[]) {
71  bool verbose = false; // used to set verbose false on non-root processors
72  bool verbose1 = false; // user's command-line argument
73  // Check if we should print results to standard out
74  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose1 = true;
75 
76  int ierr = 0;
77  int returnierr = 0;
78  int size = 1;
79  int rank = 0;
80 
81  if (verbose1)
82  cout << Epetra_Version() << endl << endl;
83 
84  // Test Epetra_SerialComm
85  if(verbose1) cout << "Testing Epetra_SerialComm..." << endl;
86  Epetra_SerialComm serialcomm;
87  if (verbose1) cout << serialcomm << endl;
88  ierr = checkRankAndSize(serialcomm, verbose1, rank, size);
89  EPETRA_TEST_ERR(ierr,returnierr);
90  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
91  // method testing
92  int numProc = serialcomm.NumProc();
93  ierr = checkCommMethods(serialcomm, verbose, verbose1, numProc, rank);
94  EPETRA_TEST_ERR(ierr,returnierr);
95  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
96  // clone
97  if(verbose1) cout << "SerialComm Clone.." << endl;
98  Epetra_Comm* cloned_serialcomm = serialcomm.Clone();
99  ierr = checkCommMethods(*cloned_serialcomm, verbose, verbose1, numProc, rank);
100  delete cloned_serialcomm;
101  EPETRA_TEST_ERR(ierr,returnierr);
102  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
103  // check inner data class
104  ierr = checkSerialDataClass(verbose1);
105  EPETRA_TEST_ERR(ierr,returnierr);
106  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
107 
108  Epetra_Distributor* serialdistr = serialcomm.CreateDistributor();
109  ierr = checkDistributor(serialdistr, serialcomm);
110  delete serialdistr;
111  EPETRA_TEST_ERR(ierr, returnierr);
112 
113  // Test Epetra_MpiComm
114 #ifdef EPETRA_MPI
115  // Initialize MPI
116  if(verbose1) cout << "Testing Epetra_MpiComm..." << endl;
117  MPI_Init(&argc,&argv);
118  MPI_Comm_size(MPI_COMM_WORLD, &size);
119  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
120  Epetra_MpiComm petracomm( MPI_COMM_WORLD );
121  ierr = checkRankAndSize(petracomm, verbose1, rank, size);
122  EPETRA_TEST_ERR(ierr,returnierr);
123  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
124  MPI_Comm MPIComm1 = petracomm.Comm();
125  int size1, rank1;
126  MPI_Comm_size(MPIComm1, &size1);
127  MPI_Comm_rank(MPIComm1, &rank1);
128  if (verbose1) cout << petracomm << ". Using MPI_Comm from Petra_Comm:"
129  << " Processor " << rank1 << " of " << size1
130  << " (should be the same)." << endl;
131  EPETRA_TEST_ERR(!(rank1==rank),ierr);
132  EPETRA_TEST_ERR(!(size1==size),ierr);
133  checkBarrier(petracomm, verbose1, rank);
134 
135  // method testing
136  numProc = petracomm.NumProc();
137  ierr = checkCommMethods(petracomm, verbose, verbose1, numProc, rank);
138  EPETRA_TEST_ERR(ierr,returnierr);
139  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
140 
141  // clone
142  if(verbose1) cout << "MpiComm Clone.." << endl;
143  Epetra_Comm* cloned_mpicomm = petracomm.Clone();
144  ierr = checkCommMethods(*cloned_mpicomm, verbose, verbose1, numProc, rank);
145  delete cloned_mpicomm;
146  EPETRA_TEST_ERR(ierr,returnierr);
147  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
148 
149  // check inner data class
150  petracomm.Barrier();
151  ierr = checkMpiDataClass(verbose1);
152  EPETRA_TEST_ERR(ierr,returnierr);
153  if (verbose1 && ierr==0) cout << "Checked OK\n\n" <<endl;
154 
155  Epetra_Distributor* plldistr = petracomm.CreateDistributor();
156  ierr = checkDistributor(plldistr, petracomm);
157  delete plldistr;
158  EPETRA_TEST_ERR(ierr, returnierr);
159 
160  petracomm.Barrier();
161  MPI_Finalize();
162 #endif
163 
164  return(returnierr);
165 }
166 
167 //=============================================================================
168 void checkBarrier(Epetra_Comm& petracomm, bool verbose, int rank) {
169  // Do some timing to test barrier function
170  int MyPID = petracomm.MyPID();
171  Epetra_Time before_barrier(petracomm);
172  Epetra_Time after_barrier(petracomm);
173  // Give each processor rank+1 amount of work
174  // Time before barrier should increase roughly linearly
175  // Time after barrier should be same for all processors
176  double sum = 0.0;
177  for (int j=0; j<rank+1; j++)
178  for (int i=0; i<1000000; i++)
179  sum += ((double )rand())/((double) RAND_MAX);
180  sum /= rank+1;
181  if (verbose) cout << "Processor " << MyPID
182  << " Time to reach barrier: "
183  << before_barrier.ElapsedTime() << endl;
184  petracomm.Barrier();
185  if (verbose) cout << "Processor " << MyPID << " Sum result "
186  << sum << " Time to get beyond barrier: "
187  << after_barrier.ElapsedTime() << endl;
188 
189  petracomm.Barrier();
190 }
191 
192 //=============================================================================
193 int checkRankAndSize(Epetra_Comm& petracomm, bool verbose, int rank, int size) {
194  int ierr = 0;
195  //if(verbose) cout << "CRS Breakpoint 1" << endl;
196  int MyPID = petracomm.MyPID();
197  //if(verbose) cout << "CRS Breakpoint 2" << endl;
198  int NumProc = petracomm.NumProc();
199  EPETRA_TEST_ERR(!(MyPID==rank),ierr);
200  EPETRA_TEST_ERR(!(NumProc==size),ierr);
201  petracomm.Barrier();
202  return(ierr);
203 }
204 
205 //=============================================================================
206 int checkCommMethods(Epetra_Comm& petracomm, bool verbose, bool verbose1, int& NumProc, int& rank) {
207  int i,j;
208  int forierr = 0;
209  int ierr = 0;
210 
211  verbose = (petracomm.MyPID() == 0); // Turn verbose on;
212  // it is always false in main.
213 
214  // Some vars needed for the following tests
215  int count = 4;
216  int* iInputs = new int[count]; // General array for int type tests
217  for (i=0; i<count; i++)
218  iInputs[i] = 10*(i + rank - 2) + rank;
219  // if these values are changed, the expected maxs, mins, sums, etc must also change.
220  //NOTE: Broadcst() does not use these values. The lines that need to be changed are located
221  //in the "Values for ****** tests" sections directly below.
222 
223  double* dInputs = new double[count]; // General array for double type tests
224  for (i=0; i<count; i++)
225  dInputs[i] = pow(2.0,i-rank);
226  // if these values are changed, the expected maxs, mins, sums, etc must also change.
227  //NOTE: Broadcst() does not use these values. The lines that need to be changed are located
228  //in the "Values for ****** tests" sections directly below.
229 
230 
231  // Values for Broadcast tests
232  int* iVals = new int[count];
233  if (rank == 0) {
234  for (i=0; i<count; i++)
235  iVals[i] = i; // if these values are changed, the values in iBVals must also be changed
236  }
237 
238  int* iBVals = new int[count]; // Values to be checked against the values broadcast to the non root processors
239  for (i=0; i<count; i++)
240  iBVals[i] = i; // if these values are changed, the values in iVals must also be changed
241  double* dVals = new double[count];
242  if (rank == 0) {
243  for (i=0; i<count; i++)
244  dVals[i] = double(i); // if these values are changed, the values in dBVals must also be changed
245  }
246  double* dBVals = new double[count]; // Values to be checked against the values broadcast to the non root processors
247  for (i=0; i<count; i++)
248  dBVals[i] = i; // if these values are changed, the values in dVals must also be changed
249 
250  long long* llVals = new long long[count];
251  if (rank == 0) {
252  for (i=0; i<count; i++)
253  llVals[i] = (long long)i+INT_MAX; // if these values are changed, the values in llBVals must also be changed
254  }
255  long long* llBVals = new long long[count]; // Values to be checked against the values broadcast to the non root processors
256  for (i=0; i<count; i++)
257  llBVals[i] = (long long)i+INT_MAX; // if these values are changed, the values in dVals must also be changed
258 
259  const char *cConst = "Heidi, do you want a cookie?";
260  int cCount = strlen(cConst)+1;
261  char* cVals = new char[cCount];
262  if (rank == 0) {
263  strcpy(cVals, cConst); // if these values are changed,
264  cVals[cCount-1] = '\0'; // the values in cBVals must also be changed
265  }
266  char* cBVals = new char[cCount]; // Values to be checked against the values
267  // broadcast to the non root processors
268  strcpy(cBVals, cConst); // if these values are changed,
269  cBVals[cCount-1] = '\0'; // the values in cVals must also be changed
270 
271  // Values for MaxAll tests
272  int* iMyGlobalMaxs = new int[count];
273  for (i=0; i<count; i++)
274  iMyGlobalMaxs[i]=10 * (i + NumProc-1 -2) + NumProc-1; // if these values are changed, iInput must be changed
275  //as well as all other values dependent on iInput
276  double* dMyGlobalMaxs = new double[count];
277  for (i=0; i<count; i++)
278  dMyGlobalMaxs[i]= pow(2.0,i); //if these values are changed, dInput must be changed
279  //as well as all other values dependent on dInput
280 
281 
282  // Values for MinAll tests
283  int* iMyGlobalMins = new int[count];
284  for (i=0; i<count; i++)
285  iMyGlobalMins[i]= 10 * (i - 2); //if these values are changed, iInput must be changed
286  //as well as all other values dependent on iInput
287  double* dMyGlobalMins = new double[count];
288  for (i=0; i<count; i++)
289  dMyGlobalMins[i]= pow(2.0,i-(NumProc-1)); //if these values are changed, dInput must be changed
290  //as well as all other values dependent on dInput
291 
292 
293  // Values for SumAll tests
294  int* iMyGlobalSums = new int[count];
295  for (i=0; i<count; i++){
296  iMyGlobalSums[i]=0;
297  for (j=0; j<NumProc; j++)
298  iMyGlobalSums[i] += 10*(i+j-2) + j;// if these values are changed, iInput must be changed
299  } //as well as all other values dependent on iInput
300 
301  double* dMyGlobalSums = new double[count];
302  for (i=0; i<count; i++){
303  dMyGlobalSums[i]=0;
304  for (j=0; j<NumProc; j++)
305  dMyGlobalSums[i] += pow(2.0,i-j);// if these values are changed, dInput must be changed
306  } //as well as all other values dependent on dInput
307 
308 
309  // Values for ScanSum tests
310  int* iMyScanSums = new int[count];
311  for (i=0; i<count; i++)
312  iMyScanSums[i] = int((rank+1)*(10*(2*i+rank-4)+rank)*.5);// if these values are changed,
313  //iInput must be changed as well as
314  //all other values dependent on iInput
315  double* dMyScanSums = new double[count];
316  for (i=0; i<count; i++) {
317  dMyScanSums[i] = 0;
318  for (j=0; j<=rank; j++)
319  dMyScanSums[i] += pow(2.0,i-j); //if these values are changed, dInput must be changed
320  } //as well as all other values dependent on dInput
321 
322 
323  // Values for Gather tests
324  int totalVals = count*NumProc;
325  int* iMyOrderedVals = new int[totalVals];
326  double* dMyOrderedVals = new double[totalVals];
327  int k=0;
328  for (j=0; j<NumProc; j++) {
329  for (i=0; i<count; i++) {
330  iMyOrderedVals[k] = 10*(i + j - 2) + j;; // if these values are changed, iInput must be changed
331  //as well as all other values dependent on iInput
332  dMyOrderedVals[k] = pow(2.0,i-j); // if these values are changed, dInput must be changed
333  //as well as all other values dependent on dInput
334  k++;
335  }
336  }
337  petracomm.Barrier();
338 
339 
340  // Method testing section
341  // Test the Broadcast functions
342  EPETRA_TEST_ERR(petracomm.Broadcast(iVals,count,0),ierr);
343  if (verbose1) {
344  if (rank == 0)
345  cout << "The values on the root processor are: ";
346  else
347  cout << "The values on processor " << rank << " are: ";
348  for (i=0; i<count; i++)
349  cout << iVals[i] << " ";
350  cout << endl;
351  }
352  // ierr = 0; need to track errors the whole way through the file - this line of code seems like a bad idea
353  forierr = 0;
354  for (i=0; i<count; i++)
355  forierr += !(iVals[i] == iBVals[i]); // otherwise Broadcast didn't occur properly
356  EPETRA_TEST_ERR(forierr,ierr);
357  delete [] iVals;
358  delete [] iBVals;
359  petracomm.Barrier();
360  if (verbose) cout << endl << "Broadcast (type int) test passed!" << endl << endl;// If test gets to here the test passed,
361  //only output on one node
362  petracomm.Barrier();
363 
364  EPETRA_TEST_ERR(petracomm.Broadcast(dVals,count,0),ierr);
365  if (verbose1) {
366  if (rank == 0)
367  cout << "The values on the root processor are: ";
368  else
369  cout << "The values on processor " << rank << " are: ";
370  for (i=0; i<count; i++)
371  cout << dVals[i] << " ";
372  cout << endl;
373  }
374  forierr = 0;
375  for (i=0; i<count; i++)
376  forierr += !(dVals[i] == dBVals[i]); // otherwise Broadcast didn't occur properly
377  EPETRA_TEST_ERR(forierr,ierr);
378  delete [] dVals;
379  delete [] dBVals;
380  petracomm.Barrier();
381  if (verbose) cout << endl << "Broadcast (type double) test passed!" << endl << endl;// If test gets to here the test passed,
382  //only output on one node
383  petracomm.Barrier();
384 
385  EPETRA_TEST_ERR(petracomm.Broadcast(llVals,count,0),ierr);
386  if (verbose1) {
387  if (rank == 0)
388  cout << "The values on the root processor are: ";
389  else
390  cout << "The values on processor " << rank << " are: ";
391  for (i=0; i<count; i++)
392  cout << llVals[i] << " ";
393  cout << endl;
394  }
395  forierr = 0;
396  for (i=0; i<count; i++)
397  forierr += !(llVals[i] == llBVals[i]); // otherwise Broadcast didn't occur properly
398  EPETRA_TEST_ERR(forierr,ierr);
399  delete [] llVals;
400  delete [] llBVals;
401  petracomm.Barrier();
402  if (verbose) cout << endl << "Broadcast (type long long) test passed!" << endl << endl;// If test gets to here the test passed,
403  //only output on one node
404  petracomm.Barrier();
405 
406  EPETRA_TEST_ERR(petracomm.Broadcast(cVals,cCount,0),ierr);
407  if (verbose1) {
408  if (rank == 0)
409  cout << "The values on the root processor are: " << cVals << endl;
410  else
411  cout << "The values on processor " << rank << " are: " << cVals << endl;
412  }
413  forierr = 0;
414  for (i=0; i<cCount; i++)
415  forierr += !(cVals[i] == cBVals[i]); // otherwise Broadcast didn't work.
416  EPETRA_TEST_ERR(forierr,ierr);
417  delete [] cVals;
418  delete [] cBVals;
419  petracomm.Barrier();
420  // If test gets to here the test passed,
421  if (verbose)
422  cout << endl << "Broadcast (type char) test passed!" << endl << endl;
423  //only output on one node
424  petracomm.Barrier();
425 
426  // Test the MaxAll functions
427  int* iGlobalMaxs = new int[count];
428  if (verbose1) {
429  cout << "The values on processor " << rank << " are: ";
430  for (i=0; i<count; i++)
431  cout << iInputs[i] << " ";
432  cout << endl;
433  }
434  EPETRA_TEST_ERR(petracomm.MaxAll(iInputs,iGlobalMaxs,count),ierr);
435  petracomm.Barrier();
436 
437  if (verbose1) {
438  cout << "The max values according to processor " << rank << " are: ";
439  for (i=0; i<count; i++)
440  cout << iGlobalMaxs[i] << " ";
441  cout << endl;
442  }
443  forierr = 0;
444  for (i=0; i<count; i++)
445  forierr += !(iMyGlobalMaxs[i] == iGlobalMaxs[i]);
446  EPETRA_TEST_ERR(forierr,ierr);
447 
448  delete [] iGlobalMaxs;
449  delete [] iMyGlobalMaxs;
450  petracomm.Barrier();
451  if (verbose) cout << endl << "MaxAll (type int) test passed!" << endl << endl;// If test gets to here the test passed,
452  //only output on one node
453  petracomm.Barrier();
454 
455  double* dGlobalMaxs = new double[count];
456  if (verbose1) {
457  cout << "The values on processor " << rank << " are: ";
458  for (i=0; i<count; i++)
459  cout << dInputs[i] << " ";
460  cout << endl;
461  }
462  EPETRA_TEST_ERR(petracomm.MaxAll(dInputs,dGlobalMaxs,count),ierr);
463  petracomm.Barrier();
464 
465  if (verbose1) {
466  cout << "The max values according to processor " << rank << " are: ";
467  for (i=0; i<count; i++)
468  cout << dGlobalMaxs[i] << " ";
469  cout << endl;
470  }
471  forierr = 0;
472  for (i=0; i<count; i++)
473  forierr += !(Epetra_Util::Chop(dMyGlobalMaxs[i] - dGlobalMaxs[i]) == 0);
474  EPETRA_TEST_ERR(forierr,ierr);
475  delete [] dGlobalMaxs;
476  delete [] dMyGlobalMaxs;
477  petracomm.Barrier();
478  if (verbose) cout << endl << "MaxAll (type double) test passed!" << endl << endl;// If test gets to here the test passed,
479  //only output on one node
480  petracomm.Barrier();
481 
482 
483  // Test the MinAll functions
484  int* iGlobalMins = new int[count];
485  if (verbose1) {
486  cout << "The values on processor " << rank << " are: ";
487  for (i=0; i<count; i++)
488  cout << iInputs[i] << " ";
489  cout << endl;
490  }
491  EPETRA_TEST_ERR(petracomm.MinAll(iInputs,iGlobalMins,count),ierr);
492  petracomm.Barrier();
493 
494  if (verbose1) {
495  cout << "The min values according to processor " << rank << " are: ";
496  for (i=0; i<count; i++)
497  cout << iGlobalMins[i] << " ";
498  cout << endl;
499  }
500  forierr = 0;
501  for (i=0; i<count; i++)
502  forierr += !(iMyGlobalMins[i] == iGlobalMins[i]); // otherwise calculated min is wrong
503  EPETRA_TEST_ERR(forierr,ierr);
504  delete [] iGlobalMins;
505  delete [] iMyGlobalMins;
506  petracomm.Barrier();
507  if (verbose) cout << endl << "MinAll (type int) test passed!" << endl << endl;// If test gets to here the test passed,
508  //only output on one node
509  petracomm.Barrier();
510 
511  double* dGlobalMins = new double[count];
512  if (verbose1) {
513  cout << "The values on processor " << rank << " are: ";
514  for (i=0; i<count; i++)
515  cout << dInputs[i] << " ";
516  cout << endl;
517  }
518  EPETRA_TEST_ERR(petracomm.MinAll(dInputs,dGlobalMins,count),ierr);
519  petracomm.Barrier();
520 
521  if (verbose1) {
522  cout << "The min values according to processor " << rank << " are: ";
523  for (i=0; i<count; i++)
524  cout << dGlobalMins[i] << " ";
525  cout << endl;
526  }
527  forierr = 0;
528  for (i=0; i<count; i++)
529  forierr += !(Epetra_Util::Chop(dMyGlobalMins[i] - dGlobalMins[i]) == 0); // otherwise calculated min is wrong
530  EPETRA_TEST_ERR(forierr,ierr);
531  delete [] dGlobalMins;
532  delete [] dMyGlobalMins;
533  petracomm.Barrier();
534  if (verbose) cout << endl << "MinAll (type double) test passed!" << endl << endl;// If test gets to here the test passed,
535  //only output on one node
536  petracomm.Barrier();
537 
538 
539  // Test the SumAll functions
540  int* iGlobalSums = new int[count];
541  if (verbose1) {
542  cout << "The values on processor " << rank << " are: ";
543  for (i=0; i<count; i++)
544  cout << iInputs[i] << " ";
545  cout << endl;
546  }
547  EPETRA_TEST_ERR(petracomm.SumAll(iInputs,iGlobalSums,count),ierr);
548  petracomm.Barrier();
549 
550  if (verbose1) {
551  cout << "The sums of all values according to processor " << rank << " are: ";
552  for (i=0; i<count; i++)
553  cout << iGlobalSums[i] << " ";
554  cout << endl;
555  }
556  forierr = 0;
557  for (i=0; i<count; i++)
558  forierr += !(iMyGlobalSums[i] == iGlobalSums[i]); // otherwise calculated sum is wrong
559  EPETRA_TEST_ERR(forierr,ierr);
560  delete [] iGlobalSums;
561  delete [] iMyGlobalSums;
562  petracomm.Barrier();
563  if (verbose) cout << endl << "SumAll (type int) test passed!" << endl << endl;// If test gets to here the test passed,
564  //only output on one node
565  petracomm.Barrier();
566 
567  double* dGlobalSums = new double[count];
568  if (verbose1) {
569  cout << "The values on processor " << rank << " are: ";
570  for (i=0; i<count; i++)
571  cout << dInputs[i] << " ";
572  cout << endl;
573  }
574  EPETRA_TEST_ERR(petracomm.SumAll(dInputs,dGlobalSums,count),ierr);
575  petracomm.Barrier();
576 
577  if (verbose1) {
578  cout << "The sums of all values according to processor " << rank << " are: ";
579  for (i=0; i<count; i++)
580  cout << dGlobalSums[i] << " ";
581  cout << endl;
582  }
583  forierr = 0;
584  for (i=0; i<count; i++)
585  forierr += !(Epetra_Util::Chop(dMyGlobalSums[i] - dGlobalSums[i]) == 0); // otherwise calculated sum is wrong
586  EPETRA_TEST_ERR(forierr,ierr);
587 
588  delete [] dGlobalSums;
589  delete [] dMyGlobalSums;
590  petracomm.Barrier();
591  if (verbose) cout << endl << "SumAll (type double) test passed!" << endl << endl;// If test gets to here the test passed,
592  //only output on one node
593  petracomm.Barrier();
594 
595 
596  // Test the ScanSum functions
597  int* iScanSums = new int[count];
598  if (verbose1) {
599  cout << "The values on processor " << rank << " are: ";
600  for (i=0; i<count; i++)
601  cout << iInputs[i] << " ";
602  cout << endl;
603  }
604 
605  EPETRA_TEST_ERR(petracomm.ScanSum(iInputs,iScanSums,count),ierr);
606  petracomm.Barrier();
607 
608  if (verbose1) {
609  cout << "The sums of all values on processors 0 - " << rank << " are: ";
610  for (i=0; i<count; i++) {
611  cout << iScanSums[i] << " ";
612  }
613  cout << endl;
614  }
615  forierr = 0;
616  for (i=0; i<count; i++)
617  forierr += !(iMyScanSums[i] == iScanSums[i]);
618  EPETRA_TEST_ERR(forierr,ierr);
619  delete [] iScanSums;
620  delete [] iMyScanSums;
621  petracomm.Barrier();
622  if (verbose) cout << endl << "ScanSum (type int) test passed!" << endl << endl;// If test gets to here the test passed,
623  //only output on one node
624  petracomm.Barrier();
625 
626  double* dScanSums = new double[count];
627  if (verbose1) {
628  cout << "The values on processor " << rank << " are: ";
629  for (i=0; i<count; i++)
630  cout << dInputs[i] << " ";
631  cout << endl;
632  }
633 
634  EPETRA_TEST_ERR(petracomm.ScanSum(dInputs,dScanSums,count),ierr);
635  petracomm.Barrier();
636 
637  if (verbose1) {
638  cout << "The sums of all values on processors 0 - " << rank << " are: ";
639  for (i=0; i<count; i++) {
640  cout << dScanSums[i] << " ";
641  }
642  cout << endl;
643  }
644  forierr = 0;
645  for (i=0; i<count; i++)
646  forierr += !(Epetra_Util::Chop(dMyScanSums[i] - dScanSums[i])== 0);
647  EPETRA_TEST_ERR(forierr,ierr);
648  delete [] dScanSums;
649  delete [] dMyScanSums;
650  petracomm.Barrier();
651  if (verbose) cout << endl << "ScanSum (type double) test passed!" << endl << endl;// If test gets to here the test passed,
652  //only output on one node
653  petracomm.Barrier();
654 
655 
656  // Test the Gather functions
657  int* iOrderedVals = new int[totalVals];
658  if (verbose1) {
659  cout << "The values on processor " << rank << " are: ";
660  for (i=0; i<count; i++)
661  cout << iInputs[i] << " ";
662  cout << endl;
663  }
664 
665  EPETRA_TEST_ERR(petracomm.GatherAll(iInputs,iOrderedVals,count),ierr);
666  petracomm.Barrier();
667 
668  if (verbose1) {
669  cout << "The combined list of all values from all processors according to processor " << rank << " is: ";
670  for (i=0; i<totalVals; i++) {
671  cout << iOrderedVals[i] << " ";
672  }
673  cout << endl;
674  }
675  forierr = 0;
676  for (i=0; i<totalVals; i++)
677  forierr += !(iMyOrderedVals[i] == iOrderedVals[i]);
678  EPETRA_TEST_ERR(forierr,ierr);
679  delete [] iOrderedVals;
680  delete [] iMyOrderedVals;
681  petracomm.Barrier();
682  if (verbose) cout << endl << "GatherAll (type int) test passed!" << endl << endl;// If test gets to here the test passed,
683  //only output on one node
684  petracomm.Barrier();
685 
686  double* dOrderedVals = new double[totalVals];
687  if (verbose1) {
688  cout << "The values on processor " << rank << " are: ";
689  for (i=0; i<count; i++)
690  cout << dInputs[i] << " ";
691  cout << endl;
692  }
693 
694  EPETRA_TEST_ERR(petracomm.GatherAll(dInputs,dOrderedVals,count),ierr);
695  petracomm.Barrier();
696 
697  if (verbose1) {
698  cout << "The combined list of all values from all processors according to processor " << rank << " is: ";
699  for (i=0; i<totalVals; i++) {
700  cout << dOrderedVals[i] << " ";
701  }
702  cout << endl;
703  }
704  forierr = 0;
705  for (i=0; i<totalVals; i++)
706  forierr += !(Epetra_Util::Chop(dMyOrderedVals[i] - dOrderedVals[i]) == 0);
707  EPETRA_TEST_ERR(forierr,ierr);
708  delete [] dOrderedVals;
709  delete [] dMyOrderedVals;
710  petracomm.Barrier();
711  if (verbose) cout << endl << "GatherAll (type double) test passed!" << endl << endl;// If test gets to here the test passed,
712  //only output on one node
713  petracomm.Barrier();
714 
715  delete[] dInputs;
716  delete[] iInputs;
717 
718  return(ierr);
719 }
720 
721 //=============================================================================
722 int checkSerialDataClass(bool verbose) {
723  int ierr = 0;
724  if(verbose) cout << "Testing Reference Counting... ";
726  int c1count = c1.ReferenceCount();
727  const Epetra_SerialCommData* c1addr = c1.DataPtr();
728  EPETRA_TEST_ERR(!(c1count==1),ierr); // count should be 1
729  if(verbose) cout << "Default constructor. \nc1= " << c1count << " " << c1addr << endl;
730 
731  Epetra_SerialComm* c2 = new Epetra_SerialComm(c1);
732  int c2count = c2->ReferenceCount();
733  const Epetra_SerialCommData* c2addr = c2->DataPtr();
734  int c1countold = c1count;
735  c1count = c1.ReferenceCount();
736  EPETRA_TEST_ERR(!(c2count==c1count && c1count==(c1countold+1)),ierr); // both counts should be 2
737  EPETRA_TEST_ERR(!(c1addr==c2addr),ierr); // addresses should be same
738  if(verbose) cout << "Copy constructor(heap). \nc1= " << c1count << " " << c1addr
739  << "\nc2= " << c2count << " " << c2addr << endl;
740  delete c2;
741  c1countold = c1count;
742  c1count = c1.ReferenceCount();
743  EPETRA_TEST_ERR(!(c1count==c1countold-1),ierr); // count should have decremented (to 1)
744  EPETRA_TEST_ERR(!(c1addr==c1.DataPtr()),ierr); // c1addr should be unchanged
745  if(verbose) cout << "c2 Deleted. \nc1= " << c1count << " " << c1addr << endl;
746  { // inside own set of brackets so that c2a will be automatically at end of brackets
747  // so that we can test to make sure objects on the stack deallocate correctly
748  Epetra_SerialComm c2a(c1);
749  c2count = c2a.ReferenceCount();
750  c2addr = c2a.DataPtr();
751  c1countold = c1count;
752  c1count = c1.ReferenceCount();
753  EPETRA_TEST_ERR(!(c2count==c1count && c1count==c1countold+1),ierr); // both counts should be 2
754  EPETRA_TEST_ERR(!(c1addr==c2addr),ierr); // addresses should be same
755  if(verbose) cout << "Copy constructor(stack). \nc1= " << c1count << " " << c1addr
756  << "\nc2a= " << c2count << " " << c2addr << endl;
757  }
758  c1countold = c1count;
759  c1count = c1.ReferenceCount();
760  EPETRA_TEST_ERR(!(c1count==c1countold-1),ierr); // count should have decremented (to 1)
761  EPETRA_TEST_ERR(!(c1addr==c1.DataPtr()),ierr); // c1addr should be unchanged
762  if(verbose) cout << "c2a Destroyed. \nc1= " << c1count << " " << c1addr << endl;
763  if(verbose) cout << "Assignment operator, post construction" << endl;
765  int c3count = c3.ReferenceCount();
766  const Epetra_SerialCommData* c3addr = c3.DataPtr();
767  EPETRA_TEST_ERR(!(c3count==1),ierr); // c3count should be 1 initially
768  EPETRA_TEST_ERR(!(c1addr!=c3addr),ierr); // c1 and c3 should have different ptr addresses
769  if(verbose)cout << "Prior to assignment: \nc1=" << c1count << " " << c1addr
770  << "\nc3=" << c3count << " " << c3addr << endl;
771  c3 = c1;
772  c3count = c3.ReferenceCount();
773  c3addr = c3.DataPtr();
774  c1countold = c1count;
775  c1count = c1.ReferenceCount();
776  EPETRA_TEST_ERR(!(c3count==c1count && c1count==c1countold+1),ierr); // both counts should be 2
777  EPETRA_TEST_ERR(!(c1addr==c3addr),ierr); // addresses should be same
778  if(verbose)cout << "After assignment: \nc1=" << c1count << " " << c1addr
779  << "\nc3=" << c3count << " " << c3addr << endl;
780  return(ierr);
781 }
782 
783 //=============================================================================
784 #ifdef EPETRA_MPI
785 int checkMpiDataClass(bool verbose) {
786  int ierr = 0;
787  if(verbose) cout << "Testing Reference Counting... ";
788  Epetra_MpiComm c1( MPI_COMM_WORLD );
789  int c1count = c1.ReferenceCount();
790  const Epetra_MpiCommData* c1addr = c1.DataPtr();
791  EPETRA_TEST_ERR(!(c1count==1),ierr); // count should be 1
792  if(verbose) cout << "Default constructor. \nc1= " << c1count << " " << c1addr << endl;
793 
794  Epetra_MpiComm* c2 = new Epetra_MpiComm(c1);
795  int c2count = c2->ReferenceCount();
796  const Epetra_MpiCommData* c2addr = c2->DataPtr();
797  int c1countold = c1count;
798  c1count = c1.ReferenceCount();
799  EPETRA_TEST_ERR(!(c2count==c1count && c1count==(c1countold+1)),ierr); // both counts should be 2
800  EPETRA_TEST_ERR(!(c1addr==c2addr),ierr); // addresses should be same
801  if(verbose) cout << "Copy constructor(heap). \nc1= " << c1count << " " << c1addr
802  << "\nc2= " << c2count << " " << c2addr << endl;
803  delete c2;
804  c1countold = c1count;
805  c1count = c1.ReferenceCount();
806  EPETRA_TEST_ERR(!(c1count==c1countold-1),ierr); // count should have decremented (to 1)
807  EPETRA_TEST_ERR(!(c1addr==c1.DataPtr()),ierr); // c1addr should be unchanged
808  if(verbose) cout << "c2 Deleted. \nc1= " << c1count << " " << c1addr << endl;
809  { // inside own set of brackets so that c2a will be automatically at end of brackets
810  // so that we can test to make sure objects on the stack deallocate correctly
811  Epetra_MpiComm c2a(c1);
812  c2count = c2a.ReferenceCount();
813  c2addr = c2a.DataPtr();
814  c1countold = c1count;
815  c1count = c1.ReferenceCount();
816  EPETRA_TEST_ERR(!(c2count==c1count && c1count==c1countold+1),ierr); // both counts should be 2
817  EPETRA_TEST_ERR(!(c1addr==c2addr),ierr); // addresses should be same
818  if(verbose) cout << "Copy constructor(stack). \nc1= " << c1count << " " << c1addr
819  << "\nc2a= " << c2count << " " << c2addr << endl;
820  }
821  c1countold = c1count;
822  c1count = c1.ReferenceCount();
823  EPETRA_TEST_ERR(!(c1count==c1countold-1),ierr); // count should have decremented (to 1)
824  EPETRA_TEST_ERR(!(c1addr==c1.DataPtr()),ierr); // c1addr should be unchanged
825  if(verbose) cout << "c2a Destroyed. \nc1= " << c1count << " " << c1addr << endl;
826  if(verbose) cout << "Assignment operator, post construction" << endl;
827  Epetra_MpiComm c3( MPI_COMM_WORLD );
828  int c3count = c3.ReferenceCount();
829  const Epetra_MpiCommData* c3addr = c3.DataPtr();
830  EPETRA_TEST_ERR(!(c3count==1),ierr); // c3count should be 1 initially
831  EPETRA_TEST_ERR(!(c1addr!=c3addr),ierr); // c1 and c3 should have different ptr addresses
832  if(verbose)cout << "Prior to assignment: \nc1=" << c1count << " " << c1addr
833  << "\nc3=" << c3count << " " << c3addr << endl;
834  c3 = c1;
835  c3count = c3.ReferenceCount();
836  c3addr = c3.DataPtr();
837  c1countold = c1count;
838  c1count = c1.ReferenceCount();
839  EPETRA_TEST_ERR(!(c3count==c1count && c1count==c1countold+1),ierr); // both counts should be 2
840  EPETRA_TEST_ERR(!(c1addr==c3addr),ierr); // addresses should be same
841  if(verbose)cout << "After assignment: \nc1=" << c1count << " " << c1addr
842  << "\nc3=" << c3count << " " << c3addr << endl;
843  return(ierr);
844 }
845 #endif
846 
848  Epetra_Comm& Comm)
849 {
850  int numprocs = Comm.NumProc();
851 
852  int numExportIDs = numprocs;
853  int* exportPIDs = new int[numExportIDs];
854  for(int p=0; p<numExportIDs; ++p) {
855  exportPIDs[p] = p;
856  }
857 
858  bool deterministic = true;
859  int numRemoteIDs = 0;
860 
861  int err = distr->CreateFromSends(numExportIDs, exportPIDs,
862  deterministic, numRemoteIDs);
863 
864  //numRemoteIDs should equal numExportIDs.
865 
866  int returnValue = numRemoteIDs == numExportIDs ? 0 : -99;
867 
868  delete [] exportPIDs;
869 
870  if (returnValue + err != 0) {
871  return(returnValue + err);
872  }
873 
874  int* exportIDs = new int[numExportIDs];
875  for(int i=0; i<numExportIDs; ++i) {
876  exportIDs[i] = i+1;
877  }
878 
879  int len_imports = 0;
880  char* imports = NULL;
881 
882  err = distr->Do((char*)exportIDs, sizeof(int),
883  len_imports, imports);
884 
885  delete [] exportIDs;
886 
887  if (len_imports > 0) {
888  delete [] imports;
889  }
890 
891  return(err);
892 }
893 
894 /*
895  end of file cxx_main.cpp
896 */
const Epetra_MpiCommData * DataPtr() const
Returns a pointer to the MpiCommData instance this MpiComm uses.
Epetra_Distributor: The Epetra Gather/Scatter Setup Base Class.
virtual int GatherAll(double *MyVals, double *AllVals, int Count) const =0
Epetra_Comm All Gather function.
#define EPETRA_TEST_ERR(a, b)
Epetra_Comm * Clone() const
Clone method.
double ElapsedTime(void) const
Epetra_Time elapsed time function.
virtual int CreateFromSends(const int &NumExportIDs, const int *ExportPIDs, bool Deterministic, int &NumRemoteIDs)=0
Create Distributor object using list of process IDs to which we export.
int NumProc() const
Returns total number of processes.
void Barrier() const
Epetra_MpiComm Barrier function.
int ReferenceCount() const
Returns the reference count of MpiCommData.
virtual int MinAll(double *PartialMins, double *GlobalMins, int Count) const =0
Epetra_Comm Global Min function.
Epetra_MpiComm: The Epetra MPI Communication Class.
std::string Epetra_Version()
virtual void Barrier() const =0
Epetra_Comm Barrier function.
Epetra_Comm * Clone() const
Clone method.
virtual int MyPID() const =0
Return my process ID.
const Epetra_SerialCommData * DataPtr() const
Returns a pointer to the SerialCommData instance this SerialComm uses.
Epetra_Time: The Epetra Timing Class.
Definition: Epetra_Time.h:75
virtual int MaxAll(double *PartialMaxs, double *GlobalMaxs, int Count) const =0
Epetra_Comm Global Max function.
virtual int SumAll(double *PartialSums, double *GlobalSums, int Count) const =0
Epetra_Comm Global Sum function.
Epetra_Distributor * CreateDistributor() const
Create a distributor object.
void checkBarrier(Epetra_Comm &petracomm, bool verbose, int rank)
Epetra_SerialCommData: The Epetra Serial Communication Data Class.
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:73
MPI_Comm Comm() const
Extract MPI Communicator from a Epetra_MpiComm object.
int NumProc() const
Returns total number of processes (always returns 1 for SerialComm).
virtual int Broadcast(double *MyVals, int Count, int Root) const =0
Epetra_Comm Broadcast function.
virtual int Do(char *export_objs, int obj_size, int &len_import_objs, char *&import_objs)=0
Execute plan on buffer of export objects in a single step.
Epetra_SerialComm: The Epetra Serial Communication Class.
Epetra_MpiCommData: The Epetra Mpi Communication Data Class.
virtual int NumProc() const =0
Returns total number of processes.
int checkMpiDataClass(bool verbose)
int ReferenceCount() const
Returns the reference count of SerialCommData.
int checkRankAndSize(Epetra_Comm &petracomm, bool verbose, int rank, int size)
static double Chop(const double &Value)
Epetra_Util Chop method. Return zero if input Value is less than ChopValue.
Definition: Epetra_Util.cpp:65
int main(int argc, char *argv[])
int checkDistributor(Epetra_Distributor *distr, Epetra_Comm &Comm)
Epetra_Distributor * CreateDistributor() const
Create a distributor object.
virtual int ScanSum(double *MyVals, double *ScanSums, int Count) const =0
Epetra_Comm Scan Sum function.
int checkCommMethods(Epetra_Comm &petracomm, bool verbose, bool verbose1, int &NumProc, int &rank)
int checkSerialDataClass(bool verbose)