33#include <Teuchos_DefaultComm.hpp>
34#include <Teuchos_XMLObject.hpp>
35#include <Teuchos_FileInputSource.hpp>
43using Teuchos::ParameterList;
46using Teuchos::ArrayRCP;
47using Teuchos::XMLObject;
54using std::ostringstream;
57#define ERRMSG(msg) if (rank == 0){ std::cerr << "FAIL: " << msg << std::endl; }
58#define EXC_ERRMSG(msg, e) \
59if (rank==0){ std::cerr << "FAIL: " << msg << std::endl << e.what() << std::endl;}
62 Teuchos::ParameterList & plist)
65 Teuchos::XMLParameterListReader reader;
66 plist = reader.toParameterList(xml);
74 Teuchos::ParameterList zoltan2Parameters;
77 if (plist.isSublist(
"Zoltan2Parameters")) {
79 ParameterList &sub = plist.sublist(
"Zoltan2Parameters");
80 zoltan2Parameters.setParameters(sub);
85 queue<ParameterList> &problems,
86 queue<ParameterList> &comparisons,
87 const RCP<
const Teuchos::Comm<int> > & comm)
89 int rank = comm->getRank();
92 Teuchos::FileInputSource inputSource(inputFileName);
94 std::cout <<
"Input file source: " << inputFileName << std::endl;
101 xmlInput = inputSource.getObject();
110 for(
int i = 0; i < xmlInput.numChildren(); i++)
115 if(plist.name() ==
"Comparison") {
116 comparisons.emplace(plist);
119 problems.emplace(plist);
128 std::ostringstream & msg,
129 const ParameterList &problem_parameters) {
130 #define ANALYZE_METRICS(adapterClass, metricAnalyzerClass) \
131 RCP<EvaluateBaseClass<adapterClass>> pCast = \
132 rcp_dynamic_cast<EvaluateBaseClass<adapterClass>>( \
133 evaluateFactory->getEvaluateClass()); \
134 if(pCast == Teuchos::null) throw std::logic_error( \
135 "Bad evaluate class cast in analyzeMetrics!" ); \
136 metricAnalyzerClass analyzer(pCast); \
137 return analyzer.analyzeMetrics( \
138 problem_parameters.sublist("Metrics"), msg);
140 #define ANALYZE_METRICS_PARTITIONING(adapterClass) \
141 ANALYZE_METRICS(adapterClass, \
142 MetricAnalyzerEvaluatePartition<adapterClass>)
144 #define ANALYZE_METRICS_ORDERING(adapterClass) \
145 ANALYZE_METRICS(adapterClass, \
146 MetricAnalyzerEvaluateOrdering<adapterClass>)
148 if(evaluateFactory->getProblemName() ==
"partitioning") {
151 else if(evaluateFactory->getProblemName() ==
"ordering") {
155 throw std::logic_error(
156 "analyzeMetrics not implemented for this problem type!" );
162 RCP<ProblemFactory> problemFactory) {
163 #define GET_LOCAL_ORDERING(adapterClass) \
164 return (rcp_dynamic_cast<OrderingProblem<adapterClass>>( \
165 problemFactory->getProblem()))->getLocalOrderingSolution();
171 #define GET_PROBLEM_PARTS(adapterClass) \
172 return (rcp_dynamic_cast<PartitioningProblem<adapterClass>>( \
173 problemFactory->getProblem()))->getSolution().getPartListView();
179 #define GET_IDS_VIEW(adapterClass) \
180 return dynamic_cast<adapterClass*>( \
181 adapterFactory->getMainAdapter())->getIDsView(Ids);
183 throw std::logic_error(
"getIDsView() failed to match adapter name" );
187 const ParameterList &problem_parameters,
188 bool bHasComparisons,
189 RCP<ComparisonHelper> & comparison_helper,
190 const RCP<
const Teuchos::Comm<int> > & comm)
199 int rank = comm->getRank();
200 if(!problem_parameters.isParameter(
"kind"))
203 std::cout <<
"Problem kind not provided" << std::endl;
207 if(!problem_parameters.isParameter(
"InputAdapterParameters"))
210 std::cout <<
"Input adapter parameters not provided" << std::endl;
214 if(!problem_parameters.isParameter(
"Zoltan2Parameters"))
217 std::cout <<
"Zoltan2 problem parameters not provided" << std::endl;
223 std::cout <<
"\n\nRunning test: " << problem_parameters.name() << std::endl;
231 comparison_helper->AddSource(problem_parameters.name(), comparison_source);
232 comparison_source->addTimer(
"adapter construction time");
233 comparison_source->addTimer(
"problem construction time");
234 comparison_source->addTimer(
"solve time");
239 const ParameterList &adapterPlist =
240 problem_parameters.sublist(
"InputAdapterParameters");
241 comparison_source->timers[
"adapter construction time"]->start();
247 comparison_source->timers[
"adapter construction time"]->stop();
249 comparison_source->adapterFactory = adapterFactory;
255 string adapter_name = adapterPlist.get<
string>(
"input adapter");
257 ParameterList zoltan2_parameters =
258 const_cast<ParameterList &
>(problem_parameters.sublist(
"Zoltan2Parameters"));
260 std::cout << std::endl;
263 comparison_source->timers[
"problem construction time"]->start();
264 std::string problem_kind = problem_parameters.get<std::string>(
"kind");
266 std::cout <<
"Creating a new " << problem_kind <<
" problem." << std::endl;
273 #ifdef HAVE_ZOLTAN2_MPI
279 std::cout <<
"Using input adapter type: " + adapter_name << std::endl;
282 comparison_source->problemFactory = problemFactory;
287 comparison_source->timers[
"solve time"]->start();
289 problemFactory->getProblem()->solve();
291 comparison_source->timers[
"solve time"]->stop();
293 std::cout << problem_kind +
" problem solved." << std::endl;
298 if(problem_kind ==
"partitioning") {
302 i < adapterFactory->getMainAdapter()->getLocalNumIDs(); i++) {
303 std::cout << rank <<
" LID " << i
304 <<
" GID " << kddIDs[i]
310 if (adapter_name ==
"XpetraCrsGraph") {
315 dynamic_cast<const xCG_xCG_t *
>(adapterFactory->getMainAdapter());
319 const gno_t *vertexIds;
324 for (
int edim = 0; edim < ewgtDim; edim++) {
328 for (
lno_t i=0; i < localNumObj; i++)
329 for (
offset_t j=offsets[i]; j < offsets[i+1]; j++)
330 std::cout << edim <<
" " << vertexIds[i] <<
" "
331 << adjIds[j] <<
" " <<
weights[stride*j] << std::endl;
339 bool bSuccess =
true;
340 if(problem_parameters.isSublist(
"Metrics") || bHasComparisons) {
348 std::cout <<
"Create evaluate class for: " + problem_kind << std::endl;
352 comparison_source->evaluateFactory = evaluateFactory;
354 std::ostringstream msgSummary;
356 evaluateFactory->getEvaluateClass()->printMetrics(msgSummary);
358 std::cout << msgSummary.str();
361 std::ostringstream msgResults;
362 if (!
analyzeMetrics(evaluateFactory, msgResults, problem_parameters)) {
364 std::cout <<
"MetricAnalyzer::analyzeMetrics() "
365 <<
"returned false and the test is FAILED." << std::endl;
368 std::cout << msgResults.str();
373 if (problem_kind ==
"ordering") {
374 std::cout <<
"\nLet's examine the solution..." << std::endl;
378 std::cout <<
"Number of column blocks: "
381 if (localOrderingSolution->
havePerm()) {
382 std::cout <<
"permutation: {";
384 std::cout <<
" " << x;
385 std::cout <<
"}" << std::endl;
389 std::cout <<
"inverse permutation: {";
391 std::cout <<
" " << x;
392 std::cout <<
"}" << std::endl;
396 std::cout <<
"separator range: {";
398 std::cout <<
" " << x;
399 std::cout <<
"}" << std::endl;
403 std::cout <<
"separator tree: {";
405 std::cout <<
" " << x;
406 std::cout <<
"}" << std::endl;
413 comparison_source->printTimers();
425bool mainExecute(
int narg,
char *arg[], RCP<
const Comm<int> > &comm)
430 int rank = comm->getRank();
436 string inputFileName(
"");
438 inputFileName = arg[1];
441 std::cout <<
"\nFAILED to specify xml input file!" << std::endl;
442 std::ostringstream msg;
443 msg <<
"\nStandard use of test_driver.cpp:\n";
444 msg <<
"mpiexec -n <procs> ./Zoltan2_test_driver.exe <input_file.xml>\n";
445 std::cout << msg.str() << std::endl;
453 queue<ParameterList> problems, comparisons;
463 const ParameterList inputParameters = problems.front();
464 if(inputParameters.name() !=
"InputParameters")
467 std::cout <<
"InputParameters not defined. Testing FAILED." << std::endl;
488 while (!problems.empty()) {
492 if (!
run(uinput, problems.front(), !comparisons.empty(),
493 comparison_manager, comm)) {
494 std::cout <<
"Problem run returned false" << std::endl;
505 while (!comparisons.empty()) {
506 if (!comparison_manager->Compare(comparisons.front(),comm)) {
508 std::cout <<
"Comparison manager returned false so the "
509 <<
"test should fail." << std::endl;
518 std::cout <<
"\nFAILED to load input data source. Skipping "
519 "all tests." << std::endl;
529 Tpetra::ScopeGuard tscope(&narg, &arg);
530 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
533 int rank = comm->getRank();
538 catch(std::logic_error &e) {
543 std::cout <<
"Test driver for rank " << rank
544 <<
" caught the following exception: " << e.what() << std::endl;
548 catch(std::runtime_error &e) {
549 std::cout <<
"Test driver for rank " << rank
550 <<
" caught the following exception: " << e.what() << std::endl;
553 catch(std::exception &e) {
554 std::cout <<
"Test driver for rank " << rank
555 <<
" caught the following exception: " << e.what() << std::endl;
565 reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
566 &result, &resultReduced);
571 std::cout <<
"Test Driver with " << comm->getSize()
572 <<
" processes has reduced result code " << resultReduced
573 <<
" and is exiting in the "
574 << ((resultReduced == 0 ) ?
"PASSED" :
"FAILED") <<
" state."
Generate Adapter for testing purposes.
Defines the BasicIdentifierAdapter class.
Store and compare solution sets from different problems.
Returns a pointer to new test classes. Is not responsible for memory management!
Defines Parameter related enumerators, declares functions.
Tpetra::Map ::global_ordinal_type zgno_t
keep typedefs that commonly appear in many places localized
#define Z2_TEST_UPCAST(adptr, TEMPLATE_ACTION)
Defines XpetraCrsGraphAdapter class.
Defines the XpetraCrsMatrixAdapter class.
Defines the XpetraMultiVectorAdapter.
A class for comparing solutions, metrics, and timing data of Zoltan2 problems.
A class used to save problem solutions and timers.
typename InputTraits< User >::lno_t lno_t
typename InputTraits< User >::scalar_t scalar_t
typename InputTraits< userTypes_t >::gno_t gno_t
ArrayRCP< index_t > & getPermutationRCPConst(bool inverse=false) const
Get (local) permuted GIDs by const RCP.
bool haveSeparators() const
Do we have the separators?
index_t getNumSeparatorBlocks() const
Get number of separator column blocks.
bool haveSeparatorTree() const
Do we have the separator tree?
bool haveSeparatorRange() const
Do we have the separator range?
bool haveInverse() const
Do we have the inverse permutation?
ArrayRCP< index_t > & getSeparatorRangeRCPConst() const
Get separator range by const RCP.
bool havePerm() const
Do we have the direct permutation?
ArrayRCP< index_t > & getSeparatorTreeRCPConst() const
Get separator tree by const RCP.
size_t getLocalNumVertices() const override
Returns the number of vertices on this process.
void getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const override
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const override
int getNumWeightsPerEdge() const override
Returns the number (0 or greater) of edge weights.
void getVertexIDsView(const gno_t *&ids) const override
ProblemFactory class contains 1 static factory method.
ProblemFactory class contains 1 static factory method.
Zoltan2::XpetraCrsGraphAdapter< xcrsGraph_t, tMVector_t > xCG_xCG_t
void createAllParameters(Teuchos::ParameterList &pList)
Create a list of all Zoltan2 parameters and validators.
#define GET_PROBLEM_PARTS(adapterClass)
bool analyzeMetrics(RCP< EvaluateFactory > evaluateFactory, std::ostringstream &msg, const ParameterList &problem_parameters)
const zpart_t * getPartListView(RCP< ProblemFactory > problemFactory)
#define GET_IDS_VIEW(adapterClass)
#define ANALYZE_METRICS_ORDERING(adapterClass)
#define EXC_ERRMSG(msg, e)
LocalOrderingSolution< zlno_t > * getLocalOrderingSolution(RCP< ProblemFactory > problemFactory)
bool getParameterLists(const string &inputFileName, queue< ParameterList > &problems, queue< ParameterList > &comparisons, const RCP< const Teuchos::Comm< int > > &comm)
bool run(const UserInputForTests &uinput, const ParameterList &problem_parameters, bool bHasComparisons, RCP< ComparisonHelper > &comparison_helper, const RCP< const Teuchos::Comm< int > > &comm)
void getIDsView(RCP< AdapterFactory > adapterFactory, const zgno_t *&Ids)
bool mainExecute(int narg, char *arg[], RCP< const Comm< int > > &comm)
void xmlToModelPList(const Teuchos::XMLObject &xml, Teuchos::ParameterList &plist)
#define ANALYZE_METRICS_PARTITIONING(adapterClass)
#define GET_LOCAL_ORDERING(adapterClass)