Zoltan2
AdapterForTests.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 
51 #ifndef ADAPTERFORTESTS
52 #define ADAPTERFORTESTS
53 
54 #include <Zoltan2_Parameters.hpp>
55 #include <UserInputForTests.hpp>
56 
59 
65 
66 #ifdef HAVE_ZOLTAN2_PAMGEN
68 #endif
69 
70 #include <Teuchos_DefaultComm.hpp>
71 #include <Teuchos_XMLObject.hpp>
72 #include <Teuchos_FileInputSource.hpp>
73 
74 #include <Tpetra_MultiVector.hpp>
75 #include <Tpetra_CrsMatrix.hpp>
76 
77 #include <string>
78 #include <iostream>
79 #include <vector>
80 
81 using namespace std;
82 using Teuchos::RCP;
83 using Teuchos::ArrayRCP;
84 using Teuchos::ArrayView;
85 using Teuchos::Array;
86 using Teuchos::Comm;
87 using Teuchos::rcp;
88 using Teuchos::arcp;
89 using Teuchos::rcp_const_cast;
90 using Teuchos::ParameterList;
91 using std::string;
92 using namespace Zoltan2_TestingFramework;
93 
94 // helper struct to store both an adapter and the coordinate adapter
96 {
97  Zoltan2::BaseAdapterRoot * adapter = nullptr; // generic base class
98  EAdapterType adapterType; // convert back to proper adapter type
99 };
100 
102 {
103  AdapterWithTemplateName main; // the main adapter - never null
105 };
106 
107 /* \brief A class for constructing Zoltan2 input adapters */
109 public:
120  UserInputForTests *uinput, const ParameterList &pList,
121  const RCP<const Comm<int> > &comm);
122 
123  ~AdapterFactory(); // handles deleting BaseAdapterRoot * data for adapter
124 
126  return adaptersSet.main.adapter;
127  }
128 
130  return adaptersSet.main.adapterType;
131  }
132 
134  return adaptersSet.coordinate.adapter;
135  }
136 
138  return adaptersSet.coordinate.adapterType;
139  }
140 
141 private:
143 
153  getBasicIdentiferAdapterForInput(UserInputForTests *uinput,
154  const ParameterList &pList, const RCP<const Comm<int> > &comm);
155 
165  getXpetraMVAdapterForInput(UserInputForTests *uinput,
166  const ParameterList &pList, const RCP<const Comm<int> > &comm);
167 
177  getXpetraCrsGraphAdapterForInput(UserInputForTests *uinput,
178  const ParameterList &pList, const RCP<const Comm<int> > &comm);
179 
189  getXpetraCrsMatrixAdapterForInput(UserInputForTests *uinput,
190  const ParameterList &pList, const RCP<const Comm<int> > &comm);
191 
201  getBasicVectorAdapterForInput(UserInputForTests *uinput,
202  const ParameterList &pList, const RCP<const Comm<int> > &comm);
203 
213  getPamgenMeshAdapterForInput(UserInputForTests *uinput,
214  const ParameterList &pList, const RCP<const Comm<int> > &comm);
215 
224  template <typename T>
225  void InitializeVectorData(const RCP<T> &data,
226  vector<const zscalar_t *> &coords,
227  vector<int> & strides,
228  int stride);
229 
230 #ifdef HAVE_EPETRA_DATA_TYPES
231 
239  template <typename T>
240  void InitializeEpetraVectorData(const RCP<T> &data,
241  vector<const zscalar_t *> &coords,
242  vector<int> & strides,
243  int stride);
244 #endif
245 };
246 
247 
249  UserInputForTests *uinput,
250  const ParameterList &pList,
251  const RCP<const Comm<int> > &comm)
252 {
253  if(!pList.isParameter("input adapter"))
254  {
255  std::cerr << "Input adapter unspecified" << std::endl;
256  return;
257  }
258 
259  // pick method for chosen adapter
260  std::string input_adapter_name = pList.get<string>("input adapter");
261 
262  if(input_adapter_name == "BasicIdentifier")
263  adaptersSet.main = getBasicIdentiferAdapterForInput(uinput, pList, comm);
264  else if(input_adapter_name == "XpetraMultiVector")
265  adaptersSet.main = getXpetraMVAdapterForInput(uinput, pList, comm);
266  else if(input_adapter_name == "XpetraCrsGraph")
267  adaptersSet = getXpetraCrsGraphAdapterForInput(uinput,pList, comm);
268  else if(input_adapter_name == "XpetraCrsMatrix")
269  adaptersSet = getXpetraCrsMatrixAdapterForInput(uinput,pList, comm);
270  else if(input_adapter_name == "BasicVector")
271  adaptersSet.main = getBasicVectorAdapterForInput(uinput,pList, comm);
272  else if(input_adapter_name == "PamgenMesh")
273  adaptersSet.main = getPamgenMeshAdapterForInput(uinput,pList, comm);
274 
275  if(adaptersSet.main.adapter == nullptr) {
276  throw std::logic_error("AdapterFactory failed to create adapter!");
277  }
278 }
279 
281  if( adaptersSet.main.adapter ) {
282  delete adaptersSet.main.adapter;
283  }
284 
285  if( adaptersSet.coordinate.adapter ) {
286  delete adaptersSet.coordinate.adapter;
287  }
288 }
289 
290 
292  AdapterFactory::getBasicIdentiferAdapterForInput(UserInputForTests *uinput,
293  const ParameterList &pList,
294  const RCP<const Comm<int> > &comm)
295 {
297 
298  if(!pList.isParameter("data type"))
299  {
300  std::cerr << "Input data type unspecified" << std::endl;
301  return result;
302  }
303 
304  string input_type = pList.get<string>("data type"); // get the input type
305 
306  if (!uinput->hasInputDataType(input_type))
307  {
308  std::cerr << "Input type: " + input_type + " unavailable or misspelled."
309  << std::endl; // bad type
310  return result;
311  }
312 
313  vector<const zscalar_t *> weights;
314  std::vector<int> weightStrides;
315  const zgno_t *globalIds = NULL;
316  size_t localCount = 0;
317 
318  // get weights if any
319  // get weights if any
320  if(uinput->hasUIWeights())
321  {
322  RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
323  // copy to weight
324  size_t cols = vtx_weights->getNumVectors();
325  for (size_t i = 0; i< cols; i++) {
326  weights.push_back(vtx_weights->getData(i).getRawPtr());
327  weightStrides.push_back((int)vtx_weights->getStride());
328  }
329  }
330 
331  if(input_type == "coordinates")
332  {
333  RCP<tMVector_t> data = uinput->getUICoordinates();
334  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
335  localCount = data->getLocalLength();
336  }
337  else if(input_type == "tpetra_vector")
338  {
339  RCP<tVector_t> data = uinput->getUITpetraVector();
340  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
341  localCount = data->getLocalLength();
342  }
343  else if(input_type == "tpetra_multivector")
344  {
345  int nvec = pList.get<int>("vector_dimension");
346  RCP<tMVector_t> data = uinput->getUITpetraMultiVector(nvec);
347  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
348  localCount = data->getLocalLength();
349  }
350  else if(input_type == "tpetra_crs_graph")
351  {
352  RCP<tcrsGraph_t> data = uinput->getUITpetraCrsGraph();
353  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
354  localCount = data->getNodeNumCols();
355  }
356  else if(input_type == "tpetra_crs_matrix")
357  {
358  RCP<tcrsMatrix_t> data = uinput->getUITpetraCrsMatrix();
359  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
360  localCount = data->getNodeNumCols();
361  }
362  else if(input_type == "xpetra_vector")
363  {
364  RCP<xVector_t> data = uinput->getUIXpetraVector();
365  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
366  localCount = data->getLocalLength();
367  }
368  else if(input_type == "xpetra_multivector")
369  {
370  int nvec = pList.get<int>("vector_dimension");
371  RCP<xMVector_t> data = uinput->getUIXpetraMultiVector(nvec);
372  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
373  localCount = data->getLocalLength();
374  }
375  else if(input_type == "xpetra_crs_graph")
376  {
377  RCP<xcrsGraph_t> data = uinput->getUIXpetraCrsGraph();
378  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
379  localCount = data->getNodeNumCols();
380  }
381  else if(input_type == "xpetra_crs_matrix")
382  {
383  RCP<xcrsMatrix_t> data = uinput->getUIXpetraCrsMatrix();
384  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
385  localCount = data->getNodeNumCols();
386  }
387 #ifdef HAVE_EPETRA_DATA_TYPES
388  else if(input_type == "epetra_vector")
389  {
390  RCP<Epetra_Vector> data = uinput->getUIEpetraVector();
391  globalIds = (zgno_t *)data->Map().MyGlobalElements();
392  localCount = data->MyLength();
393  }
394  else if(input_type == "epetra_multivector")
395  {
396  int nvec = pList.get<int>("vector_dimension");
397  RCP<Epetra_MultiVector> data = uinput->getUIEpetraMultiVector(nvec);
398  globalIds = (zgno_t *)data->Map().MyGlobalElements();
399  localCount = data->MyLength();
400  }
401  else if(input_type == "epetra_crs_graph")
402  {
403  RCP<Epetra_CrsGraph> data = uinput->getUIEpetraCrsGraph();
404  globalIds = (zgno_t *)data->Map().MyGlobalElements();
405  localCount = data->NumMyCols();
406  }
407  else if(input_type == "epetra_crs_matrix")
408  {
409  RCP<Epetra_CrsMatrix> data = uinput->getUIEpetraCrsMatrix();
410  globalIds = (zgno_t *)data->Map().MyGlobalElements();
411  localCount = data->NumMyCols();
412  }
413 #endif
414 
415  result.adapterType = AT_basic_id_t;
416  result.adapter = new Zoltan2_TestingFramework::basic_id_t(zlno_t(localCount),
417  globalIds,
418  weights,weightStrides);
419  return result;
420 }
421 
422 
423 AdapterWithTemplateName AdapterFactory::getXpetraMVAdapterForInput(
424  UserInputForTests *uinput,
425  const ParameterList &pList,
426  const RCP<const Comm<int> > &comm)
427 {
429 
430  if(!pList.isParameter("data type"))
431  {
432  std::cerr << "Input data type unspecified" << std::endl;
433  return result;
434  }
435 
436  string input_type = pList.get<string>("data type");
437  if (!uinput->hasInputDataType(input_type))
438  {
439  std::cerr << "Input type:" + input_type + ", unavailable or misspelled."
440  << std::endl; // bad type
441  return result;
442  }
443 
444  vector<const zscalar_t *> weights;
445  std::vector<int> weightStrides;
446 
447  // get weights if any
448  if(uinput->hasUIWeights())
449  {
450  RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
451  // copy to weight
452  size_t weightsPerRow = vtx_weights->getNumVectors();
453  for (size_t i = 0; i< weightsPerRow; i++) {
454  weights.push_back(vtx_weights->getData(i).getRawPtr());
455  weightStrides.push_back(1);
456  }
457  }
458 
459  // set adapter
460  if(input_type == "coordinates")
461  {
462  RCP<tMVector_t> data = uinput->getUICoordinates();
463  RCP<const tMVector_t> const_data = rcp_const_cast<const tMVector_t>(data);
464  if(weights.empty())
465  result.adapter = new xMV_tMV_t(const_data);
466  else {
467  result.adapter = new xMV_tMV_t(const_data,weights,weightStrides);
468  }
469  result.adapterType = AT_xMV_tMV_t;
470  }
471  else if(input_type == "tpetra_multivector")
472  {
473  int nvec = pList.get<int>("vector_dimension");
474  RCP<tMVector_t> data = uinput->getUITpetraMultiVector(nvec);
475  RCP<const tMVector_t> const_data = rcp_const_cast<const tMVector_t>(data);
476  if(weights.empty())
477  result.adapter = new xMV_tMV_t(const_data);
478  else
479  result.adapter = new xMV_tMV_t(const_data,weights,weightStrides);
480  result.adapterType = AT_xMV_tMV_t;
481  }
482  else if(input_type == "xpetra_multivector")
483  {
484  int nvec = pList.get<int>("vector_dimension");
485  RCP<xMVector_t> data = uinput->getUIXpetraMultiVector(nvec);
486  RCP<const xMVector_t> const_data = rcp_const_cast<const xMVector_t>(data);
487  if(weights.empty())
488  result.adapter = new xMV_xMV_t(const_data);
489  else{
490  result.adapter = new xMV_xMV_t(const_data,weights,weightStrides);
491  }
492  result.adapterType = AT_xMV_xMV_t;
493  }
494 #ifdef HAVE_EPETRA_DATA_TYPES
495  else if(input_type == "epetra_multivector")
496  {
497  int nvec = pList.get<int>("vector_dimension");
498  RCP<Epetra_MultiVector> data = uinput->getUIEpetraMultiVector(nvec);
499  RCP<const Epetra_MultiVector> const_data = rcp_const_cast<const Epetra_MultiVector>(data);
500 
501  if(weights.empty())
502  result.adapter = new xMV_eMV_t(const_data);
503  else
504  result.adapter = new xMV_eMV_t(const_data,weights,weightStrides);
505  result.adapterType = AT_xMV_eMV_t;
506  }
507 #endif
508 
509  if(result.adapter == nullptr)
510  std::cerr << "Input data chosen not compatible with xpetra multi-vector adapter." << std::endl;
511 
512  return result;
513 }
514 
515 
516 AdapterWithOptionalCoordinateAdapter AdapterFactory::getXpetraCrsGraphAdapterForInput(
517  UserInputForTests *uinput,
518  const ParameterList &pList,
519  const RCP<const Comm<int> > &comm)
520 {
521 
523 
524  if(!pList.isParameter("data type"))
525  {
526  std::cerr << "Input data type unspecified" << std::endl;
527  return adapters;
528  }
529 
530  string input_type = pList.get<string>("data type");
531  if (!uinput->hasInputDataType(input_type))
532  {
533  std::cerr << "Input type: " + input_type + ", unavailable or misspelled."
534  << std::endl; // bad type
535  return adapters;
536  }
537 
538  vector<const zscalar_t *> vtx_weights;
539  vector<const zscalar_t *> edge_weights;
540  vector<int> vtx_weightStride;
541  vector<int> edge_weightStride;
542 
543  // get vtx weights if any
544  if(uinput->hasUIWeights())
545  {
546  RCP<tMVector_t> vtx_weights_tmp = uinput->getUIWeights();
547  // copy to weight
548  size_t weightsPerRow = vtx_weights_tmp->getNumVectors();
549  for (size_t i = 0; i< weightsPerRow; i++) {
550  vtx_weights.push_back(vtx_weights_tmp->getData(i).getRawPtr());
551  vtx_weightStride.push_back(1);
552  }
553  }
554 
555  // get edge weights if any
556  if(uinput->hasUIEdgeWeights())
557  {
558  RCP<tMVector_t> edge_weights_tmp = uinput->getUIEdgeWeights();
559  // copy to weight
560  size_t weightsPerRow = edge_weights_tmp->getNumVectors();
561  for (size_t i = 0; i< weightsPerRow; i++) {
562  edge_weights.push_back(edge_weights_tmp->getData(i).getRawPtr());
563  edge_weightStride.push_back(1);
564  }
565  }
566 
567  // make the coordinate adapter
568  // get an adapter for the coordinates
569  // need to make a copy of the plist and change the vector type
570  Teuchos::ParameterList pCopy(pList);
571  pCopy = pCopy.set<std::string>("data type","coordinates");
572 
573  // for coordinate adapter
574  #define SET_COORDS_INPUT_1(adapterClass) \
575  auto * ca = dynamic_cast<adapterClass*>(adapters.coordinate.adapter); \
576  if(!ca) {throw std::logic_error( "Coordinate adapter case failed!" );} \
577  ia->setCoordinateInput(ca);
578 
579  if(input_type == "tpetra_crs_graph")
580  {
581  RCP<tcrsGraph_t> data = uinput->getUITpetraCrsGraph();
582  RCP<const tcrsGraph_t> const_data = rcp_const_cast<const tcrsGraph_t>(data);
583 
584  xCG_tCG_t * ia = new xCG_tCG_t(const_data,(int)vtx_weights.size(),(int)edge_weights.size());
585  adapters.main.adapterType = AT_xCG_tCG_t;
586  adapters.main.adapter = ia;
587 
588  if(!vtx_weights.empty()) {
589  for(int i = 0; i < (int)vtx_weights.size(); i++)
590  ia->setVertexWeights(vtx_weights[i],vtx_weightStride[i],i);
591  }
592 
593  if(!edge_weights.empty()) {
594  for(int i = 0; i < (int)edge_weights.size(); i++)
595  ia->setEdgeWeights(edge_weights[i],edge_weightStride[i],i);
596  }
597 
598  if (uinput->hasUICoordinates()) {
599  adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
601  }
602  }
603  else if(input_type == "xpetra_crs_graph")
604  {
605  RCP<xcrsGraph_t> data = uinput->getUIXpetraCrsGraph();
606  RCP<const xcrsGraph_t> const_data = rcp_const_cast<const xcrsGraph_t>(data);
607 
608  xCG_xCG_t * ia = new xCG_xCG_t(const_data, (int)vtx_weights.size(), (int)edge_weights.size());
609  adapters.main.adapterType = AT_xCG_xCG_t;
610  adapters.main.adapter = ia;
611  if(!vtx_weights.empty())
612  {
613  for(int i = 0; i < (int)vtx_weights.size(); i++)
614  ia->setVertexWeights(vtx_weights[i],vtx_weightStride[i],i);
615  }
616 
617  if(!edge_weights.empty())
618  {
619  for(int i = 0; i < (int)edge_weights.size(); i++)
620  ia->setEdgeWeights(edge_weights[i],edge_weightStride[i],i);
621  }
622 
623  if (uinput->hasUICoordinates()) {
624  adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
626  }
627  }
628 #ifdef HAVE_EPETRA_DATA_TYPES
629 
630  else if(input_type == "epetra_crs_graph")
631  {
632  RCP<Epetra_CrsGraph> data = uinput->getUIEpetraCrsGraph();
633  RCP<const Epetra_CrsGraph> const_data = rcp_const_cast<const Epetra_CrsGraph>(data);
634  xCG_eCG_t * ia = new xCG_eCG_t(const_data,(int)vtx_weights.size(),(int)edge_weights.size());
635  adapters.main.adapterType = AT_xCG_eCG_t;
636  adapters.main.adapter = ia;
637  if(!vtx_weights.empty())
638  {
639  for(int i = 0; i < (int)vtx_weights.size(); i++)
640  ia->setVertexWeights(vtx_weights[i],vtx_weightStride[i],i);
641  }
642 
643  if(!edge_weights.empty())
644  {
645  for(int i = 0; i < (int)edge_weights.size(); i++)
646  ia->setEdgeWeights(edge_weights[i],edge_weightStride[i],i);
647  }
648 
649  if (uinput->hasUICoordinates()) {
650  adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
652  }
653  }
654 #endif
655 
656  if(adapters.main.adapter == nullptr) {
657  std::cerr << "Input data chosen not compatible with "
658  << "XpetraCrsGraph adapter." << std::endl;
659  return adapters;
660  }
661 
662  return adapters;
663 }
664 
665 
666 AdapterWithOptionalCoordinateAdapter AdapterFactory::getXpetraCrsMatrixAdapterForInput(
667  UserInputForTests *uinput,
668  const ParameterList &pList,
669  const RCP<const Comm<int> > &comm)
670 {
672 
673  if(!pList.isParameter("data type"))
674  {
675  std::cerr << "Input data type unspecified" << std::endl;
676  return adapters;
677  }
678 
679  string input_type = pList.get<string>("data type");
680  if (!uinput->hasInputDataType(input_type))
681  {
682  std::cerr << "Input type:" + input_type + ", unavailable or misspelled."
683  << std::endl; // bad type
684  return adapters;
685  }
686 
687  vector<const zscalar_t *> weights;
688  vector<int> strides;
689 
690  // get weights if any
691  if(uinput->hasUIWeights())
692  {
693  if(comm->getRank() == 0) cout << "Have weights...." << endl;
694  RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
695 
696  // copy to weight
697  int weightsPerRow = (int)vtx_weights->getNumVectors();
698  for (int i = 0; i< weightsPerRow; i++)
699  {
700  weights.push_back(vtx_weights->getData(i).getRawPtr());
701  strides.push_back(1);
702  }
703 
704  }
705 
706  // make the coordinate adapter
707  // get an adapter for the coordinates
708  // need to make a copy of the plist and change the vector type
709  Teuchos::ParameterList pCopy(pList);
710  pCopy = pCopy.set<std::string>("data type","coordinates");
711 
712  // for coordinate adapter
713  #define SET_COORDS_INPUT_2(adapterClass) \
714  auto * ca = dynamic_cast<adapterClass*>(adapters.coordinate.adapter); \
715  if(!ca) {throw std::logic_error( "Coordinate adapter case failed!" );} \
716  ia->setCoordinateInput(ca);
717 
718  // set adapter
719  if(input_type == "tpetra_crs_matrix")
720  {
721  if(comm->getRank() == 0) cout << "Make tpetra crs matrix adapter...." << endl;
722 
723  // get pointer to data
724  RCP<tcrsMatrix_t> data = uinput->getUITpetraCrsMatrix();
725  RCP<const tcrsMatrix_t> const_data = rcp_const_cast<const tcrsMatrix_t>(data); // const cast data
726 
727  // new adapter
728  xCM_tCM_t *ia = new xCM_tCM_t(const_data, (int)weights.size());
729  adapters.main.adapterType = AT_xCM_tCM_t;
730  adapters.main.adapter = ia;
731 
732  // if we have weights set them
733  if(!weights.empty())
734  {
735  for(int i = 0; i < (int)weights.size(); i++)
736  ia->setWeights(weights[i],strides[i],i);
737  }
738 
739  if (uinput->hasUICoordinates()) {
740  adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
742  }
743  }
744  else if(input_type == "xpetra_crs_matrix")
745  {
746  RCP<xcrsMatrix_t> data = uinput->getUIXpetraCrsMatrix();
747  RCP<const xcrsMatrix_t> const_data = rcp_const_cast<const xcrsMatrix_t>(data);
748 
749  // new adapter
750  xCM_xCM_t *ia = new xCM_xCM_t(const_data, (int)weights.size());
751  adapters.main.adapterType = AT_xCM_xCM_t;
752  adapters.main.adapter = ia;
753 
754  // if we have weights set them
755  if(!weights.empty())
756  {
757  for(int i = 0; i < (int)weights.size(); i++)
758  ia->setWeights(weights[i],strides[i],i);
759  }
760 
761  if (uinput->hasUICoordinates()) {
762  adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
764  }
765  }
766 #ifdef HAVE_EPETRA_DATA_TYPES
767  else if(input_type == "epetra_crs_matrix")
768  {
769  RCP<Epetra_CrsMatrix> data = uinput->getUIEpetraCrsMatrix();
770  RCP<const Epetra_CrsMatrix> const_data = rcp_const_cast<const Epetra_CrsMatrix>(data);
771 
772  // new adapter
773  xCM_eCM_t *ia = new xCM_eCM_t(const_data, (int)weights.size());
774  adapters.main.adapterType = AT_xCM_eCM_t;
775  adapters.main.adapter = ia;
776 
777  // if we have weights set them
778  if(!weights.empty())
779  {
780  for(int i = 0; i < (int)weights.size(); i++)
781  ia->setWeights(weights[i],strides[i],i);
782  }
783 
784  if (uinput->hasUICoordinates()) {
785  adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
787  }
788  }
789 #endif
790 
791  if(adapters.main.adapter == nullptr)
792  {
793  std::cerr << "Input data chosen not compatible with "
794  << "XpetraCrsMatrix adapter." << std::endl;
795  return adapters;
796  }
797 
798  return adapters;
799 }
800 
801 AdapterWithTemplateName AdapterFactory::getBasicVectorAdapterForInput(
802  UserInputForTests *uinput,
803  const ParameterList &pList,
804  const RCP<const Comm<int> > &comm)
805 {
806 
808 
809  if(!pList.isParameter("data type"))
810  {
811  std::cerr << "Input data type unspecified" << std::endl;
812  return result;
813  }
814 
815  string input_type = pList.get<string>("data type");
816  if (!uinput->hasInputDataType(input_type))
817  {
818  std::cerr << "Input type:" + input_type + ", unavailable or misspelled."
819  << std::endl; // bad type
820  return result;
821  }
822 
823  std::vector<const zscalar_t *> weights;
824  std::vector<int> weightStrides;
825  const zgno_t * globalIds;
826  zlno_t localCount = 0;
827 
828  // get weights if any
829  // get weights if any
830  if(uinput->hasUIWeights())
831  {
832  RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
833  // copy to weight
834  size_t cols = vtx_weights->getNumVectors();
835  for (size_t i = 0; i< cols; i++) {
836  weights.push_back(vtx_weights->getData(i).getRawPtr());
837  weightStrides.push_back(1);
838  }
839  }
840 
841  // get vector stride
842  int stride = 1;
843  if(pList.isParameter("stride"))
844  stride = pList.get<int>("stride");
845 
847 
848  if(input_type == "coordinates")
849  {
850  RCP<tMVector_t> data = uinput->getUICoordinates();
851  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
852  localCount = static_cast<zlno_t>(data->getLocalLength());
853 
854  // get strided data
855  std::vector<const zscalar_t *> coords;
856  std::vector<int> entry_strides;
857  InitializeVectorData(data,coords,entry_strides,stride);
858 
859 
860 
861  if (weights.empty()) {
862  size_t dim = coords.size(); //BDD add NULL for constructor call
863  size_t push_null = 3-dim;
864  for (size_t i = 0; i < push_null; i ++) coords.push_back(NULL);
866  zlno_t(localCount),
867  globalIds,
868  coords[0],
869  coords[1],coords[2],
870  stride, stride, stride);
871  } else if (weights.size() == 1) {
872  size_t dim = coords.size(); //BDD add NULL for constructor call
873  size_t push_null = 3-dim;
874  for (size_t i = 0; i < push_null; i ++) coords.push_back(NULL);
876  zlno_t(localCount),
877  globalIds,
878  coords[0],
879  coords[1],coords[2],
880  stride, stride, stride,
881  true,
882  weights[0],
883  weightStrides[0]);
884  } else { // More than one weight per ID
886  zlno_t(localCount),
887  globalIds,
888  coords, entry_strides,
889  weights, weightStrides);
890  }
891  }
892  else if(input_type == "tpetra_vector")
893  {
894  RCP<tVector_t> data = uinput->getUITpetraVector();
895  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
896  localCount = static_cast<zlno_t>(data->getLocalLength());
897 
898  // get strided data
899  vector<const zscalar_t *> coords;
900  vector<int> entry_strides;
901  InitializeVectorData(data,coords,entry_strides,stride);
902 
903  if(weights.empty())
904  {
905  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
906  coords[0], entry_strides[0]);
907  }else{
908  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
909  coords[0], entry_strides[0],
910  true,
911  weights[0],
912  weightStrides[0]);
913 
914  }
915 
916  }
917  else if(input_type == "tpetra_multivector")
918  {
919  int nvec = pList.get<int>("vector_dimension");
920 
921  RCP<tMVector_t> data = uinput->getUITpetraMultiVector(nvec);
922  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
923  localCount = static_cast<zlno_t>(data->getLocalLength());
924 
925  // get strided data
926  vector<const zscalar_t *> coords;
927  vector<int> entry_strides;
928  InitializeVectorData(data,coords,entry_strides,stride);
929 
930  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
931  coords, entry_strides,
932  weights,weightStrides);
933 
934  }
935  else if(input_type == "xpetra_vector")
936  {
937  RCP<xVector_t> data = uinput->getUIXpetraVector();
938  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
939  localCount = static_cast<zlno_t>(data->getLocalLength());
940 
941  // get strided data
942  vector<const zscalar_t *> coords;
943  vector<int> entry_strides;
944  InitializeVectorData(data,coords,entry_strides,stride);
945 
946  if(weights.empty())
947  {
948  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
949  coords[0], entry_strides[0]);
950  }else{
951  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
952  coords[0], entry_strides[0],
953  true,
954  weights[0],
955  weightStrides[0]);
956 
957  }
958  }
959  else if(input_type == "xpetra_multivector")
960  {
961  int nvec = pList.get<int>("vector_dimension");
962  RCP<xMVector_t> data = uinput->getUIXpetraMultiVector(nvec);
963  globalIds = (zgno_t *)data->getMap()->getNodeElementList().getRawPtr();
964  localCount = static_cast<zlno_t>(data->getLocalLength());
965 
966  // get strided data
967  vector<const zscalar_t *> coords;
968  vector<int> entry_strides;
969  InitializeVectorData(data,coords,entry_strides,stride);
970  if(comm->getRank() == 0) cout << "size of entry strides: " << entry_strides.size() << endl;
971  if(comm->getRank() == 0) cout << "size of coords: " << coords.size() << endl;
972 
973  // make vector!
974  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
975  coords, entry_strides,
976  weights,weightStrides);
977  }
978 
979 #ifdef HAVE_EPETRA_DATA_TYPES
980  else if(input_type == "epetra_vector")
981  {
982  RCP<Epetra_Vector> data = uinput->getUIEpetraVector();
983  globalIds = (zgno_t *)data->Map().MyGlobalElements();
984  localCount = static_cast<zlno_t>(data->MyLength());
985 
986  // get strided data
987  vector<const zscalar_t *> coords;
988  vector<int> entry_strides;
989  InitializeEpetraVectorData(data,coords,entry_strides,stride);
990  if(weights.empty())
991  {
992  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
993  coords[0], entry_strides[0]);
994  }else{
995  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
996  coords[0], entry_strides[0],
997  true,
998  weights[0],
999  weightStrides[0]);
1000 
1001  }
1002 
1003  // delete [] epetravectors;
1004  }
1005  else if(input_type == "epetra_multivector")
1006  {
1007  int nvec = pList.get<int>("vector_dimension");
1008  RCP<Epetra_MultiVector> data = uinput->getUIEpetraMultiVector(nvec);
1009  globalIds = (zgno_t *)data->Map().MyGlobalElements();
1010  localCount = data->MyLength();
1011 
1012  vector<const zscalar_t *> coords;
1013  vector<int> entry_strides;
1014  InitializeEpetraVectorData(data,coords,entry_strides,stride);
1015 
1016  // make vector!
1017  result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
1018  coords, entry_strides,
1019  weights,weightStrides);
1020  }
1021 
1022 #endif
1023 
1024  return result;
1025 }
1026 
1027 template <typename T>
1028 void AdapterFactory::InitializeVectorData(const RCP<T> &data,
1029  vector<const zscalar_t *> &coords,
1030  vector<int> & strides,
1031  int stride)
1032 {
1033  // set up adapter data
1034  size_t localCount = data->getLocalLength();
1035  size_t nvecs = data->getNumVectors();
1036  size_t vecsize = data->getNumVectors() * data->getLocalLength();
1037 // printf("Number of vectors by data: %zu\n", nvecs);
1038  // printf("Size of data: %zu\n", vecsize);
1039 
1040  ArrayRCP<zscalar_t> *petravectors = new ArrayRCP<zscalar_t>[nvecs];
1041 
1042  // printf("Getting t-petra vectors...\n");
1043  for (size_t i = 0; i < nvecs; i++)
1044  petravectors[i] = data->getDataNonConst(i);
1045 
1046  // debugging
1047  // for (size_t i = 0; i < nvecs; i++){
1048  // printf("Tpetra vector %zu: {",i);
1049  //
1050  // for (size_t j = 0; j < localCount; j++)
1051  // {
1052  // printf("%1.2g ",petravectors[i][j]);
1053  // }
1054  // printf("}\n");
1055  // }
1056 
1057  size_t idx = 0;
1058  zscalar_t *coordarr = new zscalar_t[vecsize];
1059 
1060  if(stride == 1 || stride != (int)nvecs)
1061  {
1062  for (size_t i = 0; i < nvecs; i++) {
1063  for (size_t j = 0; j < localCount; j++) {
1064  coordarr[idx++] = petravectors[i][j];
1065  }
1066  }
1067  }else
1068  {
1069  for (size_t j = 0; j < localCount; j++) {
1070  for (size_t i = 0; i < nvecs; i++) {
1071  coordarr[idx++] = petravectors[i][j];
1072  }
1073  }
1074  }
1075 
1076  // debugging
1077  // printf("Made coordarr : {");
1078  // for (zlno_t i = 0; i < vecsize; i++){
1079  // printf("%1.2g ",coordarr[i]);
1080  // }
1081  // printf("}\n");
1082 
1083  // always build for dim 3
1084  coords = std::vector<const zscalar_t *>(nvecs);
1085  strides = std::vector<int>(nvecs);
1086 
1087  for (size_t i = 0; i < nvecs; i++) {
1088  if(stride == 1)
1089  coords[i] = &coordarr[i*localCount];
1090  else
1091  coords[i] = &coordarr[i];
1092 
1093  strides[i] = stride;
1094  }
1095 
1096  // debugging
1097  // printf("Made coords...\n");
1098  // for (size_t i = 0; i < nvecs; i++){
1099  // const zscalar_t * tmp = coords[i];
1100  // printf("coord %zu: {",i);
1101  // for(size_t j = 0; j < localCount; j++)
1102  // {
1103  // printf("%1.2g ", tmp[j]);
1104  // }
1105  // printf("}\n");
1106  // }
1107 
1108  // printf("clean up coordarr and tpetravectors...\n\n\n");
1109  delete [] petravectors;
1110 }
1111 
1112 #ifdef HAVE_EPETRA_DATA_TYPES
1113 
1114 template <typename T>
1115 void AdapterFactory::InitializeEpetraVectorData(const RCP<T> &data,
1116  vector<const zscalar_t *> &coords,
1117  vector<int> & strides,
1118  int stride){
1119  size_t localCount = data->MyLength();
1120  size_t nvecs = data->NumVectors();
1121  size_t vecsize = nvecs * localCount;
1122 
1123  // printf("Number of vectors by data: %zu\n", nvecs);
1124  // printf("Size of data: %zu\n", vecsize);
1125 
1126  vector<zscalar_t *> epetravectors(nvecs);
1127  zscalar_t ** arr;
1128  // printf("get data from epetra vector..\n");
1129  data->ExtractView(&arr);
1130 
1131  for(size_t k = 0; k < nvecs; k++)
1132  {
1133  epetravectors[k] = arr[k];
1134  }
1135 
1136  size_t idx = 0;
1137  basic_vector_adapter::scalar_t *coordarr =
1138  new basic_vector_adapter::scalar_t[vecsize];
1139 
1140  if(stride == 1 || stride != (int)nvecs)
1141  {
1142  for (size_t i = 0; i < nvecs; i++) {
1143  for (size_t j = 0; j < localCount; j++) {
1144  coordarr[idx++] = epetravectors[i][j];
1145  }
1146  }
1147  }else
1148  {
1149  for (size_t j = 0; j < localCount; j++) {
1150  for (size_t i = 0; i < nvecs; i++) {
1151  coordarr[idx++] = epetravectors[i][j];
1152  }
1153  }
1154  }
1155 
1156  // debugging
1157 // printf("Made coordarr : {");
1158 // for (zlno_t i = 0; i < vecsize; i++){
1159 // printf("%1.2g ",coordarr[i]);
1160 // }
1161 // printf("}\n");
1162 
1163  coords = std::vector<const zscalar_t *>(nvecs);
1164  strides = std::vector<int>(nvecs);
1165 
1166  for (size_t i = 0; i < nvecs; i++) {
1167  if(stride == 1)
1168  coords[i] = &coordarr[i*localCount];
1169  else
1170  coords[i] = &coordarr[i];
1171 
1172  strides[i] = stride;
1173  }
1174 
1175 // printf("Made coords...\n");
1176 // for (size_t i = 0; i < nvecs; i++){
1177 // const zscalar_t * tmp = coords[i];
1178 // printf("coord %zu: {",i);
1179 // for(size_t j = 0; j < localCount; j++)
1180 // {
1181 // printf("%1.2g ", tmp[j]);
1182 // }
1183 // printf("}\n");
1184 // }
1185 
1186 }
1187 #endif
1188 
1189 
1190 // pamgen adapter
1192 AdapterFactory::getPamgenMeshAdapterForInput(UserInputForTests *uinput,
1193  const ParameterList &pList,
1194  const RCP<const Comm<int> > &comm)
1195 {
1196  AdapterWithTemplateName result;
1197 
1198 #ifdef HAVE_ZOLTAN2_PAMGEN
1199  if(uinput->hasPamgenMesh())
1200  {
1201  if(uinput->hasPamgenMesh())
1202  {
1203 // if(comm->getRank() == 0) cout << "Have pamgen mesh, constructing adapter...." << endl;
1204  result.adapter =
1205  new pamgen_adapter_t(*(comm.get()), "region");
1207 // if(comm->getRank() == 0)
1208 // ia->print(0);
1209  }
1210  }else{
1211  std::cerr << "Pamgen mesh is unavailable for PamgenMeshAdapter!"
1212  << std::endl;
1213  }
1214 
1215  return result;
1216 #else
1217  throw std::runtime_error("Pamgen input requested but Trilinos is not "
1218  "built with Pamgen");
1219 #endif
1220 }
1221 #endif
1222 
1223 
InputTraits< User >::scalar_t scalar_t
RCP< tMVector_t > getUITpetraMultiVector(int nvec)
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tcrsMatrix_t
Definition: GraphModel.cpp:84
Generate input for testing purposes.
Zoltan2::XpetraCrsGraphAdapter< tcrsGraph_t, tMVector_t > xCG_tCG_t
RCP< tcrsMatrix_t > getUITpetraCrsMatrix()
Defines Parameter related enumerators, declares functions.
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
RCP< tVector_t > getUITpetraVector()
double zscalar_t
Zoltan2::BaseAdapterRoot * adapter
Zoltan2::XpetraCrsGraphAdapter< xcrsGraph_t, tMVector_t > xCG_xCG_t
Provides access for Zoltan2 to Xpetra::CrsGraph data.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Defines the PamgenMeshAdapter class.
int zlno_t
#define SET_COORDS_INPUT_2(adapterClass)
Zoltan2::BasicVectorAdapter< tMVector_t > xMV_eMV_t
RCP< tMVector_t > getUIWeights()
#define Z2_TEST_UPCAST_COORDS(adptr, TEMPLATE_ACTION)
#define SET_COORDS_INPUT_1(adapterClass)
Zoltan2::BasicVectorAdapter< tMVector_t > basic_vector_adapter
EAdapterType getCoordinateAdapterType() const
Defines the XpetraMultiVectorAdapter.
RCP< xcrsMatrix_t > getUIXpetraCrsMatrix()
Zoltan2::BasicVectorAdapter< tMVector_t > xCG_eCG_t
Defines XpetraCrsGraphAdapter class.
AdapterFactory(UserInputForTests *uinput, const ParameterList &pList, const RCP< const Comm< int > > &comm)
A class method for constructing an input adapter defind in a parameter list.
RCP< xVector_t > getUIXpetraVector()
Defines the XpetraCrsMatrixAdapter class.
bool hasInputDataType(const string &input_type)
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xMVector_t
RCP< tMVector_t > getUIEdgeWeights()
EAdapterType getMainAdapterType() const
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tcrsGraph_t
Definition: GraphModel.cpp:85
Defines the EvaluatePartition class.
Xpetra::CrsGraph< zlno_t, zgno_t, znode_t > xcrsGraph_t
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xcrsMatrix_t
BaseAdapter defines methods required by all Adapters.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
Zoltan2::XpetraMultiVectorAdapter< tMVector_t > xMV_tMV_t
An adapter for Xpetra::MultiVector.
Zoltan2::XpetraCrsMatrixAdapter< tcrsMatrix_t, tMVector_t > xCM_tCM_t
int zgno_t
Defines the BasicIdentifierAdapter class.
RCP< xMVector_t > getUIXpetraMultiVector(int nvec)
Zoltan2::BaseAdapterRoot * getMainAdapter() const
Tpetra::MultiVector< double, int, int > tMVector_t
Zoltan2::XpetraMultiVectorAdapter< xMVector_t > xMV_xMV_t
RCP< tMVector_t > getUICoordinates()
Zoltan2::BasicVectorAdapter< tMVector_t > xCM_eCM_t
Zoltan2::BaseAdapterRoot * getCoordinateAdapter() const
Defines the PartitioningProblem class.
RCP< xcrsGraph_t > getUIXpetraCrsGraph()
Zoltan2::XpetraCrsMatrixAdapter< xcrsMatrix_t, tMVector_t > xCM_xCM_t
Defines the BasicVectorAdapter class.
Zoltan2::BasicIdentifierAdapter< userTypes_t > basic_id_t
Zoltan2::BasicVectorAdapter< tMVector_t > pamgen_adapter_t
RCP< tcrsGraph_t > getUITpetraCrsGraph()