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