Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_ComparisonHelper.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Zoltan2: A package of combinatorial algorithms for scientific computing
4 //
5 // Copyright 2012 NTESS and the Zoltan2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
14 #pragma once
15 
16 #include "Zoltan2_TestHelpers.hpp"
18 #include <Zoltan2_Typedefs.hpp>
21 #include <AdapterForTests.hpp>
22 #include <Teuchos_DefaultComm.hpp>
23 #include <Teuchos_Time.hpp>
24 
25 #include <sstream>
26 #include <string>
27 #include <map>
28 #include <iostream>
29 
30 using Teuchos::Comm;
31 using Teuchos::RCP;
32 using Teuchos::ParameterList;
33 using Teuchos::Time;
34 using std::string;
35 using std::map;
36 using std::pair;
37 using Teuchos::reduceAll;
38 using namespace Zoltan2_TestingFramework;
39 
44 {
45 public:
46  /* \brief Add a timer by name to the comparison sources timers map.
47  * \param name is the name of the timer to be defined
48  */
49  void addTimer(const std::string &name)
50  {
51  timers.insert(std::pair<const std::string &, RCP<Time> >(name,rcp(new Time(name))));
52  timers[name]->enable();
53  }
54 
55  void printTimers()
56  {
57  for(auto it = timers.begin(); it != timers.end(); ++it) {
58  std::cout << it->first << " " << it->second->totalElapsedTime()
59  << std::endl;
60  }
61  }
62 
63  // TODO: Add method to print a timer summary: max/min/avg over all procs
64 
65  std::map<const std::string, RCP<Time> > timers;
66 
67  RCP<EvaluateFactory> evaluateFactory;
68  RCP<ProblemFactory> problemFactory;
69  RCP<AdapterFactory> adapterFactory;
70 };
71 
75 {
76 public:
77  /* \brief Compare the solutions, metrics or timers of two Zoltan2 solutions.
78  * \param pList is a parameter list defining the comparison
79  * \param comm is the process communicator
80  */
81  bool Compare(const ParameterList &pList, const RCP<const Comm<int> > &comm);
82 
83  /* \brief Add a new source by name to the comparison source map.
84  * \param name is the name of the new source
85  * \param source a problem source that to be used for comparison to another source
86  */
87  void AddSource(const string &name, RCP<ComparisonSource> source);
88 
89  /* \brief Return the total number of saved sources.
90  */
91  size_t getNumberOfSources() const
92  {
93  return this->sources.size();
94  }
95 
96 private:
97  map<const string,RCP<const ComparisonSource> > sources;
98 
99  /* \brief Method called to compare two solutions
100  * \param p1 is the name of problem 1
101  * \param p2 is the name of problem 2
102  * \param comm is the process communicator
103  */
104  bool CompareSolutions(const string &p1,
105  const string &p2,
106  const RCP<const Comm<int> > &comm);
107 
108  /* \brief Safely get parts list by adapter type
109  * \param problemFactory is the ProblemFactory
110  */
111  const zpart_t * getPartListView(RCP<ProblemFactory> problemFactory) const;
112 
113  /* \brief Method called to compare two paritioning solutions
114  * \param sourceA is a ptr to problem A's comparison source
115  * \param sourceB is a ptr to problem B's comparison source
116  * \param comm is the process communicator
117  */
118  bool ComparePartitionSolutions(const ComparisonSource * sourceA,
119  const ComparisonSource * sourceB,
120  const RCP<const Comm<int> > &comm);
121 
122  /* \brief Method called to compare two coloring solutions
123  * \param sourceA is a ptr to problem A's comparison source
124  * \param sourceB is a ptr to problem B's comparison source
125  * \param comm is the process communicator
126  */
127  bool CompareColoringSolutions(const ComparisonSource * sourceA,
128  const ComparisonSource * sourceB,
129  const RCP<const Comm<int> > &comm);
130 
131  /* \brief Method called to compare two ordering solutions
132  * \param sourceA is a ptr to problem A's comparison source
133  * \param sourceB is a ptr to problem B's comparison source
134  * \param comm is the process communicator
135  */
136  bool CompareOrderingSolutions(const ComparisonSource * sourceA,
137  const ComparisonSource * sourceB,
138  const RCP<const Comm<int> > &comm);
139 
140  /* \brief Safely get metric info by adapter type.
141  * \param problemFactory is the ProblemFactory
142  * \param metricInfo will be filled
143  * \param metricsPlist are the parameters to read
144  */
145  void loadMetricInfo(RCP<EvaluateFactory> problemFactory,
146  std::vector<MetricAnalyzerInfo> & metricInfo,
147  const ParameterList &metricsPlist);
148 
149  /* \brief Method called to compare the metrics/timers of two problems.
150  * \param metricsPlist is a parameter list defining the comparison
151  * \param comm is the process communicator
152  */
153  bool CompareMetrics(const ParameterList &metricsPlist,
154  const RCP<const Comm<int> > &comm);
155 
156  /* \brief Method that compares two metrics and returns a pass/fail message.
157  * \param[in] comm is the process communicator
158  * \param[in] metric is the metric to be compared to a reference metric
159  * \param[in] ref_metric is the reference metric for comparison
160  * \param[in] metricPlist is the parameter list defining the metric tolerances
161  * \param[out] msg is a returned pass/fail message
162  *
163  * \return boolean value indicated pass/fail status
164  */
165  static bool
166  metricComparisonTest(const RCP<const Comm<int> > &comm,
167  const MetricAnalyzerInfo & metric,
168  const MetricAnalyzerInfo &ref_metric,
169  const Teuchos::ParameterList & metricPlist,
170  std::ostringstream &msg);
171 
172  /* \brief Method that compares two timers and returns a pass/fail message.
173  * \param[in] comm is the process communicator
174  * \param[in] time is the timer data to be compared to a reference metric
175  * \param[in] ref_time is the reference timer for comparison
176  * \param[in] metricPlist is the parameter list defining the timer tolerances
177  * \param[out] msg is a returned pass/fail message
178  *
179  * \return boolean value indicated pass/fail status
180  */
181  static bool
182  timerComparisonTest(const RCP<const Comm<int> > &comm,
183  const double time,
184  const double ref_time,
185  const Teuchos::ParameterList & metricPlist,
186  std::ostringstream &msg);
187 
188  /* \brief Method for inserting data from all timers to a map of clocked times
189  * param[in] timers a map of timers
190  *
191  * \return a map with clocked times from timers
192  */
193  static std::map<const string, const double>
194  timerDataToMap(const map<const std::string, RCP<Time> > &timers);
195 
196 
197  /* \brief Method for extracting all methods to compare from a parameter list
198  * param[in] plist a parameter list defining 1 or more metric/timer comparisons
199  *
200  * \return a queue of metric comparison definitions
201  */
202  static std::queue<ParameterList>
203  getMetricsToCompare(const ParameterList & pList);
204 
205  static void
206  reduceWithMessage(const RCP<const Comm<int> > &comm,
207  const std::string &msg_in,
208  int &local_status, std::ostringstream &msg);
209 
210 };
211 
212 
213 void ComparisonHelper::AddSource(const string &name,
214  RCP<ComparisonSource> source)
215 {
216  typedef std::pair<const string &, RCP<const ComparisonSource> > pair_t;
217  this->sources.insert(pair_t(name, source));
218 }
219 
220 bool ComparisonHelper::Compare(const ParameterList &pList,
221  const RCP<const Comm<int> > &comm)
222 {
223  if(pList.isParameter("A") && pList.isParameter("B")) {
224  // comparing solutions
225  string pA = pList.get<string>("A");
226  if(this->sources.find(pA) == this->sources.end())
227  {
228  std::cout << "\nProblem: " + pA + ", was not saved for comparison.";
229  std::cout << "\nThis typically indicates that an error ";
230  std::cout << "occurred while running the problem.";
231  std::cout << "\nSolution comparison FAILED." << std::endl;
232  return false;
233  }
234 
235  string pB = pList.get<string>("B");
236  if(this->sources.find(pB) == this->sources.end()) {
237  std::cout << "\nProblem: " + pB + ", was not saved for comparison.";
238  std::cout << "\nThis typically indicates that an error ";
239  std::cout << "occurred while running the problem.";
240  std::cout << "\nSolution comparison FAILED." << std::endl;
241  return false;
242  }
243 
244  bool bResult = this->CompareSolutions(pA, pB, comm);
245  return bResult;
246  }
247  else if (pList.isParameter("Problem") && pList.isParameter("Reference")) {
248  // comparing metrics/timers
249  string prb = pList.get<string>("Problem");
250  if(this->sources.find(prb) == this->sources.end()) {
251  std::cout << "\nProblem: " + prb + ", was not saved for comparison.";
252  std::cout << "\nThis typically indicates that an error ";
253  std::cout << "occurred while running the problem.";
254  std::cout << "\nMetric comparison FAILED." << std::endl;
255  return false;
256  }
257 
258  string ref = pList.get<string>("Reference");
259  if(this->sources.find(ref) == this->sources.end()) {
260  std::cout << "\nReference: " + ref + ", was not saved for comparison.";
261  std::cout << "\nThis typically indicates that an error ";
262  std::cout << "occurred while running the problem.";
263  std::cout << "\nMetric comparison FAILED." << std::endl;
264  return false;
265  }
266 
267  bool bResult = this->CompareMetrics(pList, comm);
268  return bResult;
269  }
270  else if (pList.isParameter("A") || pList.isParameter("B"))
271  {
272  if(comm->getRank() == 0)
273  {
274  std::cout << "Problem A or Problem B is not specified -- check input.";
275  std::cout <<"\nSolution comparison FAILED." << std::endl;
276  }
277  }
278  else if (pList.isParameter("Problem") || pList.isParameter("Reference")) {
279  if(comm->getRank() == 0) {
280  std::cout << "Problem or reference is not specified -- check input.";
281  std::cout <<"\nMetric comparison FAILED." << std::endl;
282  }
283  }
284  else {
285  if (comm->getRank() == 0) {
286  std::cout << "ComparisonHelper did not understand how to read the xml. ";
287  std::cout << "Test FAILED." << std::endl;
288  }
289  }
290  return false;
291 }
292 
293 bool ComparisonHelper::CompareSolutions(const string &p1,
294  const string &p2,
295  const RCP<const Comm<int> > &comm)
296 {
297  if(comm->getRank() == 0) printf( "\nComparing: %s and %s\n",
298  p1.c_str(), p2.c_str());
299  auto A = this->sources[p1];
300  auto B = this->sources[p2];
301  if(A->problemFactory->getProblemName() != B->problemFactory->getProblemName()) {
302  std::cout << "Problem A and B are of a different kind and cannot be compared.";
303  std::cout <<"\nSolution comparison FAILED." << std::endl;
304  }
305  else {
306  if(A->problemFactory->getProblemName() == "partitioning") {
307  return this->ComparePartitionSolutions(A.getRawPtr(), B.getRawPtr(), comm);
308  }
309  else if(A->problemFactory->getProblemName() == "coloring") {
310  return this->CompareColoringSolutions(A.getRawPtr(), B.getRawPtr(), comm);
311  }
312  else if(A->problemFactory->getProblemName() == "ordering"){
313  return this->CompareOrderingSolutions(A.getRawPtr(), B.getRawPtr(), comm);
314  }
315  else {
316  std::cout << "Problem kind: " << A->problemFactory->getProblemName() <<
317  " not recognized. Check spelling.";
318  std::cout <<"\nSolution comparison FAILED." << std::endl;
319  }
320  }
321  return false;
322 }
323 
324 void ComparisonHelper::reduceWithMessage(const RCP<const Comm<int> > &comm,
325  const std::string &msg_in,
326  int &local_status,
327  std::ostringstream &msg) {
328  comm->barrier();
329  int global_buff;
330  Teuchos::Ptr<int> global(&global_buff);
331  reduceAll<int,int>(*comm.get(), Teuchos::EReductionType::REDUCE_MAX,
332  local_status , global);
333  local_status = *global;
334  if (local_status == 1) {
335  msg << msg_in;
336  }
337 }
338 
339 const zpart_t * ComparisonHelper::getPartListView(
340  RCP<ProblemFactory> problemFactory) const {
341  #define GET_PROBLEM_PARTS(adapterClass) \
342  return (rcp_dynamic_cast<PartitioningProblem<adapterClass>>( \
343  problemFactory->getProblem()))->getSolution().getPartListView();
344  Z2_TEST_UPCAST(problemFactory->getAdapterType(), GET_PROBLEM_PARTS)
345 }
346 
347 bool ComparisonHelper::ComparePartitionSolutions(const ComparisonSource * sourceA,
348  const ComparisonSource * sourceB,
349  const RCP<const Comm<int> > &comm)
350 {
351  int rank = comm->getRank();
352  std::ostringstream status;
353  int failed = 0;
354 
355  if(sourceA->adapterFactory->getMainAdapter()->getLocalNumIDs()
356  != sourceB->adapterFactory->getMainAdapter()->getLocalNumIDs()) {
357  failed = 1;
358  }
359 
360  ComparisonHelper::reduceWithMessage(comm,
361  "Number of parts in Solution A != Solution B. \
362  Partitioning solution comparison FAILED.",
363  failed, status);
364 
365  if (!failed) {
366  for(size_t i = 0;
367  i < sourceA->adapterFactory->getMainAdapter()->getLocalNumIDs(); i++) {
368  if(!failed && getPartListView(sourceA->problemFactory)[i] !=
369  getPartListView(sourceB->problemFactory)[i]) {
370  failed = 1;
371  ComparisonHelper::reduceWithMessage(comm,
372  "Solution sets A and B have different values for getPartListView(). "
373  "Solution comparison FAILED.", failed, status);
374  }
375  }
376  }
377 
378  if(!failed) {
379  status << "Solution sets A and B are the same. ";
380  status << "Solution set comparison PASSED.";
381  }
382 
383  if(rank == 0) {
384  std::cout << status.str() << std::endl;
385  }
386  return (failed == 0);
387 }
388 
389 
390 bool ComparisonHelper::CompareColoringSolutions(const ComparisonSource * sourceA,
391  const ComparisonSource * sourceB,
392  const RCP<const Comm<int> > &comm)
393 {
394  int rank = comm->getRank();
395  std::ostringstream status;
396  int failed = 0;
397 
398  // TO DO - implement coloring comparison
399  /*
400  if(sourceA->problemFactory->getNumColors()
401  != sourceB->problemFactory->getNumColors()) {
402  failed = 1;
403  }
404 
405  ComparisonHelper::reduceWithMessage(comm,
406  "Number of colors for Solution A != Solution B. \
407  Coloring solution comparison FAILED.",
408  failed, status);
409 
410  if (!failed) {
411  if(sourceA->problemFactory->getColorsSize()
412  != sourceB->problemFactory->getColorsSize()) {
413  failed = 1;
414  }
415  ComparisonHelper::reduceWithMessage(comm,
416  "Size of colors array for Solution A != Solution B. \
417  Coloring solution comparison FAILED.",
418  failed,
419  status);
420  }
421 
422  if (!failed) {
423  for(size_t i = 0; i < sourceA->problemFactory->getColorsSize(); i++) {
424  if (sourceA->problemFactory->getColors()[i] !=
425  sourceB->problemFactory->getColors()[i]) {
426  failed = 1; // fail
427  }
428  }
429  ComparisonHelper::reduceWithMessage(comm,
430  "Coloring solution comparison FAILED.",
431  failed,
432  status);
433  }
434  */
435 
436  if (!failed) {
437  status << "Solution sets A and B are the same. ";
438  status << "Solution set comparison PASSED.";
439  }
440 
441  if (rank == 0) {
442  std::cout << status.str() << std::endl;
443  }
444  return (failed == 0);
445 }
446 
447 bool ComparisonHelper::CompareOrderingSolutions(const ComparisonSource * sourceA,
448  const ComparisonSource * sourceB,
449  const RCP<const Comm<int> > &comm)
450 {
451  int rank = comm->getRank();
452  std::ostringstream status;
453  int failed = 0;
454 
455  // TO DO - implement ordering comparison
456 
457  if(!failed) {
458  status << "Solution sets A and B are the same. ";
459  status << "Solution set comparison PASSED.";
460  }
461 
462  if(rank == 0) {
463  std::cout << status.str() << std::endl;
464  }
465  return (failed == 0);
466 }
467 
468 // Utility function for safe type conversion of adapter
469 void ComparisonHelper::loadMetricInfo(RCP<EvaluateFactory> evaluateFactory,
470  std::vector<MetricAnalyzerInfo> & metricInfo,
471  const ParameterList &metricsPlist) {
472 
473  #define LOAD_METRIC_INFO(adapterClass, metricAnalyzerClass) \
474  RCP<EvaluateBaseClass<adapterClass>> pCast = \
475  rcp_dynamic_cast<EvaluateBaseClass<adapterClass>>(evaluateFactory->getEvaluateClass()); \
476  if(pCast == Teuchos::null) throw std::logic_error( \
477  "Bad evaluate class cast in loadMetricInfo!" ); \
478  metricAnalyzerClass analyzer(pCast); \
479  analyzer.LoadMetricInfo(metricInfo, metricsPlist.sublist("Metrics"));
480 
481  #define LOAD_METRIC_INFO_PARTITIONING(adapterClass) \
482  LOAD_METRIC_INFO(adapterClass, MetricAnalyzerEvaluatePartition<adapterClass>)
483 
484  #define LOAD_METRIC_INFO_ORDERING(adapterClass) \
485  LOAD_METRIC_INFO(adapterClass, MetricAnalyzerEvaluateOrdering<adapterClass>)
486 
487  if(evaluateFactory->getProblemName() == "partitioning") {
488  Z2_TEST_UPCAST(evaluateFactory->getAdapterType(), LOAD_METRIC_INFO_PARTITIONING)
489  }
490  else if(evaluateFactory->getProblemName() == "ordering") {
491  Z2_TEST_UPCAST(evaluateFactory->getAdapterType(), LOAD_METRIC_INFO_ORDERING)
492  }
493  else {
494  throw std::logic_error(
495  "loadMetricInfo not implemented for this problem type!" );
496  }
497 }
498 
499 // compare metrics
500 bool ComparisonHelper::CompareMetrics(const ParameterList &metricsPlist, const RCP<const Comm<int> > &comm)
501 {
502  int rank = comm->getRank();
503 
504  //get sources for problema nd reference
505  const string prb_name = metricsPlist.get<string>("Problem");
506  const string ref_name = metricsPlist.get<string>("Reference");
507  if(rank == 0) {
508  std::cout << "\nMetric/Timer comparison of: " << prb_name << " and ";
509  std::cout << ref_name <<" (reference source)\n";
510  }
511 
512  // get sources
513  RCP<const ComparisonSource> sourcePrb = this->sources[prb_name];
514  RCP<const ComparisonSource> sourceRef = this->sources[ref_name];
515 
516  // get timing data
517  std::map< const string, const double> prb_timers = this->timerDataToMap(sourcePrb->timers);
518  std::map< const string, const double> ref_timers = this->timerDataToMap(sourceRef->timers);
519 
520  // get all of the metrics to be tested
521  std::queue<ParameterList> metrics = ComparisonHelper::getMetricsToCompare(metricsPlist);
522 
523  // run comparison
524  int all_tests_pass = 1;
525  string metric_name;
526 
527  while(!metrics.empty()) {
528  // print their names...
529  std::ostringstream msg;
530  metric_name = metrics.front().name();
531 
532  if (metric_name == "Metrics") { // special key word means compare the metrics list
533  std::vector<MetricAnalyzerInfo> metricInfoSetPrb;
534  std::vector<MetricAnalyzerInfo> metricInfoSetRef;
535 
536  loadMetricInfo(sourcePrb.get()->evaluateFactory, metricInfoSetPrb, metricsPlist);
537  loadMetricInfo(sourceRef.get()->evaluateFactory, metricInfoSetRef, metricsPlist);
538 
539  // there is some redundancy here because the metric info holds both the questions and the results
540  // this happened because I wanted to reuse the MetricAnalyzer code for loading metric checks or comparisons
541  // we can iterate over either to get the questions
542  for (size_t n = 0; n < metricInfoSetPrb.size(); ++n) {
543  if(!ComparisonHelper::metricComparisonTest(comm, metricInfoSetPrb[n], metricInfoSetRef[n], metrics.front(), msg)) {
544  all_tests_pass = 0;
545  }
546  if(rank == 0) {
547  std::cout << msg.str() << std::endl;
548  }
549  }
550  }
551  else if(prb_timers.find(metric_name) != prb_timers.end() && ref_timers.find(metric_name) != ref_timers.end()) {
552  if(rank == 0) std::cout << "\ncomparing timer: " << metric_name << std::endl;
553  if(!ComparisonHelper::timerComparisonTest(comm,
554  prb_timers.at(metric_name),
555  ref_timers.at(metric_name),
556  metrics.front(), msg)) {
557  all_tests_pass = 0;
558  if (rank == 0) {
559  std::cout << "timer comparison test caused a FAILED event." << std::endl;
560  }
561  }
562 
563  if(rank == 0) {
564  std::cout << msg.str() << std::endl;
565  }
566  }
567 
568  metrics.pop();
569  }
570 
571  if(rank == 0) {
572  if(all_tests_pass == 1) {
573  std::cout << "\nAll metric/timer comparisons PASSED." << std::endl;
574  }
575  else {
576  std::cout << "\nMetric/timer metric comparisons FAILED." << std::endl;
577  }
578  }
579 
580  return (all_tests_pass == 1);
581 }
582 
583 std::map<const string, const double> ComparisonHelper::timerDataToMap(const map<const std::string, RCP<Time> > &timers)
584 {
585  typedef std::pair<const string,const double> pair_t;
586  std::map<const string, const double> time_data;
587  for (auto &i : timers) {
588  time_data.insert(pair_t(i.first, i.second->totalElapsedTime()));
589  }
590  return time_data;
591 }
592 
593 bool ComparisonHelper::metricComparisonTest(const RCP<const Comm<int> > &comm,
594  const MetricAnalyzerInfo & metric,
595  const MetricAnalyzerInfo & ref_metric,
596  const Teuchos::ParameterList & metricPlist,
597  std::ostringstream &msg)
598 {
599  // run a comparison of min and max against a given metric
600  // return an error message on failure
601  bool pass = true;
602  string test_name = metricPlist.name() + " test";
603  double ref_value = ref_metric.theValue;
604  double value = metric.theValue;
605 
606  if (ref_value == 0) {
607  throw std::logic_error( "The parameter list had a 0 value for the reference value so a percentage cannot be calculated." );
608  }
609  double percentRatio = value / ref_value;
610 
611  // want to reduce value to max value for all procs
612 
613  if (ref_metric.bFoundLowerBound) {
614  double min = ref_metric.lowerValue;
615  if (percentRatio < min) {
616  msg << test_name << " FAILED: " << ref_metric.parameterDescription << ": " << value << " is " << percentRatio << " percent of the reference value " << ref_value << ", which less than specified allowable minimum percent, " << min << ".\n";
617  pass = false;
618  }
619  else {
620  msg << test_name << " PASSED: " << ref_metric.parameterDescription << ": " << value << " is " << percentRatio << " percent of the reference value " << ref_value << ", which is greater than specified allowable minimum percent, " << min << ".\n";
621  }
622  }
623 
624  if (ref_metric.bFoundUpperBound) {
625  double max = ref_metric.upperValue;
626  if (percentRatio > max) {
627  msg << test_name << " FAILED: " << ref_metric.parameterDescription << ": " << value << " is " << percentRatio << " percent of the reference value " << ref_value << ", which is greater than specified allowable maximum percent, " << max << ".\n";
628  pass = false;
629  }
630  else {
631  msg << test_name << " PASSED: " << ref_metric.parameterDescription << ": " << value << " is " << percentRatio << " percent of the reference value " << ref_value << ", which is less than specified allowable maximum percent, " << max << ".\n";
632  }
633  }
634 
635  return pass;
636 }
637 // BDD, to do: print metrics even for pass
638 // reduce max metric to process 0
639 // print only on process 0 --- duh.
640 bool ComparisonHelper::timerComparisonTest(const RCP<const Comm<int> > &comm,
641  const double time,
642  const double ref_time,
643  const Teuchos::ParameterList & metricPlist,
644  std::ostringstream &msg)
645 {
646  // Reduce time from test
647  double global_time;
648  Teuchos::Ptr<double> global(&global_time);
649  comm->barrier();
650  reduceAll<int, double>(*comm.get(),Teuchos::EReductionType::REDUCE_MAX,time,global);
651 
652  // Reduce time from reference
653  double global_ref_time;
654  Teuchos::Ptr<double> globalRef(&global_ref_time);
655  comm->barrier();
656  reduceAll<int, double>(*comm.get(),Teuchos::EReductionType::REDUCE_MAX,ref_time,globalRef);
657 
658  // run a comparison of min and max against a given metric
659  // return an error message on failure
660  bool pass = true;
661  string test_name = metricPlist.name() + " test";
662 
663  if (metricPlist.isParameter("lower")) {
664  double min = metricPlist.get<double>("lower")*global_ref_time;
665 
666  if (global_time < min) {
667  msg << test_name << " FAILED: Minimum time, "
668  << time <<
669  "[s], less than specified allowable minimum time, " << min <<"[s]"<< ".\n";
670  pass = false;
671  }
672  else {
673  msg << test_name << " PASSED: Minimum time, "
674  << time <<
675  "[s], greater than specified allowable minimum time, " << min <<"[s]"<< ".\n";
676  }
677  }
678 
679  if (metricPlist.isParameter("upper" ) && pass != false) {
680  double max = metricPlist.get<double>("upper") * global_ref_time;
681  if (global_time > max) {
682  msg << test_name << " FAILED: Maximum time, "
683  << global_time <<
684  "[s], greater than specified allowable maximum time, " << max <<"[s]"<< ".\n";
685  pass = false;
686  }
687  else {
688  msg << test_name << " PASSED: Maximum time, "
689  << global_time <<
690  "[s], less than specified allowable maximum time, " << max <<"[s]"<< ".\n";
691  }
692  }
693 
694  return pass;
695 }
696 
697 std::queue<ParameterList> ComparisonHelper::getMetricsToCompare(const ParameterList &pList)
698 {
699  // extract all of the metrics to be tested
700  std::queue<ParameterList> metrics;
701  for(auto it = pList.begin(); it != pList.end(); ++it) {
702  if (pList.isSublist(it->first)) {
703  metrics.push(pList.sublist(it->first));
704  }
705  }
706  return metrics;
707 }
708 
keep typedefs that commonly appear in many places localized
const zpart_t * getPartListView(RCP< ProblemFactory > problemFactory)
Returns a pointer to new test classes. Is not responsible for memory management!
void AddSource(const string &name, RCP< ComparisonSource > source)
RCP< ProblemFactory > problemFactory
common code used by tests
#define LOAD_METRIC_INFO_ORDERING(adapterClass)
bool Compare(const ParameterList &pList, const RCP< const Comm< int > > &comm)
int zpart_t
#define Z2_TEST_UPCAST(adptr, TEMPLATE_ACTION)
#define LOAD_METRIC_INFO_PARTITIONING(adapterClass)
size_t getNumberOfSources() const
RCP< EvaluateFactory > evaluateFactory
A class for comparing solutions, metrics, and timing data of Zoltan2 problems.
std::map< const std::string, RCP< Time > > timers
static const std::string pass
A class used to save problem solutions and timers.
void addTimer(const std::string &name)
Generate Adapter for testing purposes.
RCP< AdapterFactory > adapterFactory
#define GET_PROBLEM_PARTS(adapterClass)