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