Zoltan2
Loading...
Searching...
No Matches
AdapterForTests.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15#ifndef ADAPTERFORTESTS
16#define ADAPTERFORTESTS
17
19#include <UserInputForTests.hpp>
20
23
29
30#ifdef HAVE_ZOLTAN2_PAMGEN
32#endif
33
34#include <Teuchos_DefaultComm.hpp>
35#include <Teuchos_XMLObject.hpp>
36#include <Teuchos_FileInputSource.hpp>
37
38#include <Tpetra_MultiVector.hpp>
39#include <Tpetra_CrsMatrix.hpp>
40
41#include <string>
42#include <iostream>
43#include <vector>
44
45using Teuchos::RCP;
46using Teuchos::ArrayRCP;
47using Teuchos::ArrayView;
48using Teuchos::Array;
49using Teuchos::Comm;
50using Teuchos::rcp;
51using Teuchos::arcp;
52using Teuchos::rcp_const_cast;
53using Teuchos::ParameterList;
54using std::string;
55using namespace Zoltan2_TestingFramework;
56
57// helper struct to store both an adapter and the coordinate adapter
59{
60 Zoltan2::BaseAdapterRoot * adapter = nullptr; // generic base class
61 EAdapterType adapterType; // convert back to proper adapter type
62};
63
65{
66 AdapterWithTemplateName main; // the main adapter - never null
68};
69
70/* \brief A class for constructing Zoltan2 input adapters */
72public:
83 UserInputForTests *uinput, const ParameterList &pList,
84 const RCP<const Comm<int> > &comm);
85
86 ~AdapterFactory(); // handles deleting BaseAdapterRoot * data for adapter
87
89 return adaptersSet.main.adapter;
90 }
91
93 return adaptersSet.main.adapterType;
94 }
95
97 return adaptersSet.coordinate.adapter;
98 }
99
101 return adaptersSet.coordinate.adapterType;
102 }
103
104private:
106
116 getBasicIdentiferAdapterForInput(UserInputForTests *uinput,
117 const ParameterList &pList, const RCP<const Comm<int> > &comm);
118
128 getXpetraMVAdapterForInput(UserInputForTests *uinput,
129 const ParameterList &pList, const RCP<const Comm<int> > &comm);
130
140 getXpetraCrsGraphAdapterForInput(UserInputForTests *uinput,
141 const ParameterList &pList, const RCP<const Comm<int> > &comm);
142
152 getXpetraCrsMatrixAdapterForInput(UserInputForTests *uinput,
153 const ParameterList &pList, const RCP<const Comm<int> > &comm);
154
164 getBasicVectorAdapterForInput(UserInputForTests *uinput,
165 const ParameterList &pList, const RCP<const Comm<int> > &comm);
166
176 getPamgenMeshAdapterForInput(UserInputForTests *uinput,
177 const ParameterList &pList, const RCP<const Comm<int> > &comm);
178
187 template <typename T>
188 void InitializeVectorData(const RCP<T> &data,
189 std::vector<const zscalar_t *> &coords,
190 std::vector<int> & strides,
191 int stride);
192
193};
194
195
197 UserInputForTests *uinput,
198 const ParameterList &pList,
199 const RCP<const Comm<int> > &comm)
200{
201 if(!pList.isParameter("input adapter"))
202 {
203 std::cerr << "Input adapter unspecified" << std::endl;
204 return;
205 }
206
207 // pick method for chosen adapter
208 std::string input_adapter_name = pList.get<string>("input adapter");
209
210 if(input_adapter_name == "BasicIdentifier")
211 adaptersSet.main = getBasicIdentiferAdapterForInput(uinput, pList, comm);
212 else if(input_adapter_name == "XpetraMultiVector")
213 adaptersSet.main = getXpetraMVAdapterForInput(uinput, pList, comm);
214 else if(input_adapter_name == "XpetraCrsGraph")
215 adaptersSet = getXpetraCrsGraphAdapterForInput(uinput,pList, comm);
216 else if(input_adapter_name == "XpetraCrsMatrix")
217 adaptersSet = getXpetraCrsMatrixAdapterForInput(uinput,pList, comm);
218 else if(input_adapter_name == "BasicVector")
219 adaptersSet.main = getBasicVectorAdapterForInput(uinput,pList, comm);
220 else if(input_adapter_name == "PamgenMesh")
221 adaptersSet.main = getPamgenMeshAdapterForInput(uinput,pList, comm);
222
223 if(adaptersSet.main.adapter == nullptr) {
224 throw std::logic_error("AdapterFactory failed to create adapter!");
225 }
226}
227
229 if( adaptersSet.main.adapter ) {
230 delete adaptersSet.main.adapter;
231 }
232
233 if( adaptersSet.coordinate.adapter ) {
234 delete adaptersSet.coordinate.adapter;
235 }
236}
237
238
240 AdapterFactory::getBasicIdentiferAdapterForInput(UserInputForTests *uinput,
241 const ParameterList &pList,
242 const RCP<const Comm<int> > &comm)
243{
245
246 if(!pList.isParameter("data type"))
247 {
248 std::cerr << "Input data type unspecified" << std::endl;
249 return result;
250 }
251
252 string input_type = pList.get<string>("data type"); // get the input type
253
254 if (!uinput->hasInputDataType(input_type))
255 {
256 std::cerr << "Input type: " + input_type + " unavailable or misspelled."
257 << std::endl; // bad type
258 return result;
259 }
260
261 std::vector<const zscalar_t *> weights;
262 std::vector<int> weightStrides;
263 const zgno_t *globalIds = NULL;
264 size_t localCount = 0;
265
266 // get weights if any
267 // get weights if any
268 if(uinput->hasUIWeights())
269 {
270 RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
271 // copy to weight
272 size_t cols = vtx_weights->getNumVectors();
273 for (size_t i = 0; i< cols; i++) {
274 weights.push_back(vtx_weights->getData(i).getRawPtr());
275 weightStrides.push_back((int)vtx_weights->getStride());
276 }
277 }
278
279 if(input_type == "coordinates")
280 {
281 RCP<tMVector_t> data = uinput->getUICoordinates();
282 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
283 localCount = data->getLocalLength();
284 }
285 else if(input_type == "tpetra_vector")
286 {
287 RCP<tVector_t> data = uinput->getUITpetraVector();
288 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
289 localCount = data->getLocalLength();
290 }
291 else if(input_type == "tpetra_multivector")
292 {
293 int nvec = pList.get<int>("vector_dimension");
294 RCP<tMVector_t> data = uinput->getUITpetraMultiVector(nvec);
295 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
296 localCount = data->getLocalLength();
297 }
298 else if(input_type == "tpetra_crs_graph")
299 {
300 RCP<tcrsGraph_t> data = uinput->getUITpetraCrsGraph();
301 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
302 localCount = data->getLocalNumCols();
303 }
304 else if(input_type == "tpetra_crs_matrix")
305 {
306 RCP<tcrsMatrix_t> data = uinput->getUITpetraCrsMatrix();
307 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
308 localCount = data->getLocalNumCols();
309 }
310 else if(input_type == "xpetra_vector")
311 {
312 RCP<xVector_t> data = uinput->getUIXpetraVector();
313 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
314 localCount = data->getLocalLength();
315 }
316 else if(input_type == "xpetra_multivector")
317 {
318 int nvec = pList.get<int>("vector_dimension");
319 RCP<xMVector_t> data = uinput->getUIXpetraMultiVector(nvec);
320 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
321 localCount = data->getLocalLength();
322 }
323 else if(input_type == "xpetra_crs_graph")
324 {
325 RCP<xcrsGraph_t> data = uinput->getUIXpetraCrsGraph();
326 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
327 localCount = data->getLocalNumCols();
328 }
329 else if(input_type == "xpetra_crs_matrix")
330 {
331 RCP<xcrsMatrix_t> data = uinput->getUIXpetraCrsMatrix();
332 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
333 localCount = data->getLocalNumCols();
334 }
335
336 result.adapterType = AT_basic_id_t;
337 result.adapter = new Zoltan2_TestingFramework::basic_id_t(zlno_t(localCount),
338 globalIds,
339 weights,weightStrides);
340 return result;
341}
342
343
344AdapterWithTemplateName AdapterFactory::getXpetraMVAdapterForInput(
345 UserInputForTests *uinput,
346 const ParameterList &pList,
347 const RCP<const Comm<int> > &comm)
348{
350
351 if(!pList.isParameter("data type"))
352 {
353 std::cerr << "Input data type unspecified" << std::endl;
354 return result;
355 }
356
357 string input_type = pList.get<string>("data type");
358 if (!uinput->hasInputDataType(input_type))
359 {
360 std::cerr << "Input type:" + input_type + ", unavailable or misspelled."
361 << std::endl; // bad type
362 return result;
363 }
364
365 std::vector<const zscalar_t *> weights;
366 std::vector<int> weightStrides;
367
368 // get weights if any
369 if(uinput->hasUIWeights())
370 {
371 RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
372 // copy to weight
373 size_t weightsPerRow = vtx_weights->getNumVectors();
374 for (size_t i = 0; i< weightsPerRow; i++) {
375 weights.push_back(vtx_weights->getData(i).getRawPtr());
376 weightStrides.push_back(1);
377 }
378 }
379
380 // set adapter
381 if(input_type == "coordinates")
382 {
383 RCP<tMVector_t> data = uinput->getUICoordinates();
384 RCP<const tMVector_t> const_data = rcp_const_cast<const tMVector_t>(data);
385 if(weights.empty())
386 result.adapter = new xMV_tMV_t(const_data);
387 else {
388 result.adapter = new xMV_tMV_t(const_data,weights,weightStrides);
389 }
390 result.adapterType = AT_xMV_tMV_t;
391 }
392 else if(input_type == "tpetra_multivector")
393 {
394 int nvec = pList.get<int>("vector_dimension");
395 RCP<tMVector_t> data = uinput->getUITpetraMultiVector(nvec);
396 RCP<const tMVector_t> const_data = rcp_const_cast<const tMVector_t>(data);
397 if(weights.empty())
398 result.adapter = new xMV_tMV_t(const_data);
399 else
400 result.adapter = new xMV_tMV_t(const_data,weights,weightStrides);
401 result.adapterType = AT_xMV_tMV_t;
402 }
403 else if(input_type == "xpetra_multivector")
404 {
405 int nvec = pList.get<int>("vector_dimension");
406 RCP<xMVector_t> data = uinput->getUIXpetraMultiVector(nvec);
407 RCP<const xMVector_t> const_data = rcp_const_cast<const xMVector_t>(data);
408 if(weights.empty())
409 result.adapter = new xMV_xMV_t(const_data);
410 else{
411 result.adapter = new xMV_xMV_t(const_data,weights,weightStrides);
412 }
413 result.adapterType = AT_xMV_xMV_t;
414 }
415
416 if(result.adapter == nullptr)
417 std::cerr << "Input data chosen not compatible with xpetra multi-vector adapter." << std::endl;
418
419 return result;
420}
421
422
423AdapterWithOptionalCoordinateAdapter AdapterFactory::getXpetraCrsGraphAdapterForInput(
424 UserInputForTests *uinput,
425 const ParameterList &pList,
426 const RCP<const Comm<int> > &comm)
427{
428
430
431 if(!pList.isParameter("data type"))
432 {
433 std::cerr << "Input data type unspecified" << std::endl;
434 return adapters;
435 }
436
437 string input_type = pList.get<string>("data type");
438 if (!uinput->hasInputDataType(input_type))
439 {
440 std::cerr << "Input type: " + input_type + ", unavailable or misspelled."
441 << std::endl; // bad type
442 return adapters;
443 }
444
445 std::vector<const zscalar_t *> vtx_weights;
446 std::vector<const zscalar_t *> edge_weights;
447 std::vector<int> vtx_weightStride;
448 std::vector<int> edge_weightStride;
449
450 // get vtx weights if any
451 if(uinput->hasUIWeights())
452 {
453 RCP<tMVector_t> vtx_weights_tmp = uinput->getUIWeights();
454 // copy to weight
455 size_t weightsPerRow = vtx_weights_tmp->getNumVectors();
456 for (size_t i = 0; i< weightsPerRow; i++) {
457 vtx_weights.push_back(vtx_weights_tmp->getData(i).getRawPtr());
458 vtx_weightStride.push_back(1);
459 }
460 }
461
462 // get edge weights if any
463 if(uinput->hasUIEdgeWeights())
464 {
465 RCP<tMVector_t> edge_weights_tmp = uinput->getUIEdgeWeights();
466 // copy to weight
467 size_t weightsPerRow = edge_weights_tmp->getNumVectors();
468 for (size_t i = 0; i< weightsPerRow; i++) {
469 edge_weights.push_back(edge_weights_tmp->getData(i).getRawPtr());
470 edge_weightStride.push_back(1);
471 }
472 }
473
474 // make the coordinate adapter
475 // get an adapter for the coordinates
476 // need to make a copy of the plist and change the vector type
477 Teuchos::ParameterList pCopy(pList);
478 pCopy = pCopy.set<std::string>("data type","coordinates");
479
480 // for coordinate adapter
481 #define SET_COORDS_INPUT_1(adapterClass) \
482 auto * ca = dynamic_cast<adapterClass*>(adapters.coordinate.adapter); \
483 if(!ca) {throw std::logic_error( "Coordinate adapter case failed!" );} \
484 ia->setCoordinateInput(ca);
485
486 if(input_type == "tpetra_crs_graph")
487 {
488 RCP<tcrsGraph_t> data = uinput->getUITpetraCrsGraph();
489 RCP<const tcrsGraph_t> const_data = rcp_const_cast<const tcrsGraph_t>(data);
490
491 xCG_tCG_t * ia = new xCG_tCG_t(const_data,(int)vtx_weights.size(),(int)edge_weights.size());
492 adapters.main.adapterType = AT_xCG_tCG_t;
493 adapters.main.adapter = ia;
494
495 if(!vtx_weights.empty()) {
496 for(int i = 0; i < (int)vtx_weights.size(); i++)
497 ia->setVertexWeights(vtx_weights[i],vtx_weightStride[i],i);
498 }
499
500 if(!edge_weights.empty()) {
501 for(int i = 0; i < (int)edge_weights.size(); i++)
502 ia->setEdgeWeights(edge_weights[i],edge_weightStride[i],i);
503 }
504
505 if (uinput->hasUICoordinates()) {
506 adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
508 }
509 }
510 else if(input_type == "xpetra_crs_graph")
511 {
512 RCP<xcrsGraph_t> data = uinput->getUIXpetraCrsGraph();
513 RCP<const xcrsGraph_t> const_data = rcp_const_cast<const xcrsGraph_t>(data);
514
515 xCG_xCG_t * ia = new xCG_xCG_t(const_data, (int)vtx_weights.size(), (int)edge_weights.size());
516 adapters.main.adapterType = AT_xCG_xCG_t;
517 adapters.main.adapter = ia;
518 if(!vtx_weights.empty())
519 {
520 for(int i = 0; i < (int)vtx_weights.size(); i++)
521 ia->setVertexWeights(vtx_weights[i],vtx_weightStride[i],i);
522 }
523
524 if(!edge_weights.empty())
525 {
526 for(int i = 0; i < (int)edge_weights.size(); i++)
527 ia->setEdgeWeights(edge_weights[i],edge_weightStride[i],i);
528 }
529
530 if (uinput->hasUICoordinates()) {
531 adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
533 }
534 }
535
536 if(adapters.main.adapter == nullptr) {
537 std::cerr << "Input data chosen not compatible with "
538 << "XpetraCrsGraph adapter." << std::endl;
539 return adapters;
540 }
541
542 return adapters;
543}
544
545
546AdapterWithOptionalCoordinateAdapter AdapterFactory::getXpetraCrsMatrixAdapterForInput(
547 UserInputForTests *uinput,
548 const ParameterList &pList,
549 const RCP<const Comm<int> > &comm)
550{
552
553 if(!pList.isParameter("data type"))
554 {
555 std::cerr << "Input data type unspecified" << std::endl;
556 return adapters;
557 }
558
559 string input_type = pList.get<string>("data type");
560 if (!uinput->hasInputDataType(input_type))
561 {
562 std::cerr << "Input type:" + input_type + ", unavailable or misspelled."
563 << std::endl; // bad type
564 return adapters;
565 }
566
567 std::vector<const zscalar_t *> weights;
568 std::vector<int> strides;
569
570 // get weights if any
571 if(uinput->hasUIWeights())
572 {
573 if(comm->getRank() == 0) std::cout << "Have weights...." << std::endl;
574 RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
575
576 // copy to weight
577 int weightsPerRow = (int)vtx_weights->getNumVectors();
578 for (int i = 0; i< weightsPerRow; i++)
579 {
580 weights.push_back(vtx_weights->getData(i).getRawPtr());
581 strides.push_back(1);
582 }
583
584 }
585
586 // make the coordinate adapter
587 // get an adapter for the coordinates
588 // need to make a copy of the plist and change the vector type
589 Teuchos::ParameterList pCopy(pList);
590 pCopy = pCopy.set<std::string>("data type","coordinates");
591
592 // for coordinate adapter
593 #define SET_COORDS_INPUT_2(adapterClass) \
594 auto * ca = dynamic_cast<adapterClass*>(adapters.coordinate.adapter); \
595 if(!ca) {throw std::logic_error( "Coordinate adapter case failed!" );} \
596 ia->setCoordinateInput(ca);
597
598 // set adapter
599 if(input_type == "tpetra_crs_matrix")
600 {
601 if(comm->getRank() == 0) std::cout << "Make tpetra crs matrix adapter...." << std::endl;
602
603 // get pointer to data
604 RCP<tcrsMatrix_t> data = uinput->getUITpetraCrsMatrix();
605 RCP<const tcrsMatrix_t> const_data = rcp_const_cast<const tcrsMatrix_t>(data); // const cast data
606
607 // new adapter
608 xCM_tCM_t *ia = new xCM_tCM_t(const_data, (int)weights.size());
609 adapters.main.adapterType = AT_xCM_tCM_t;
610 adapters.main.adapter = ia;
611
612 // if we have weights set them
613 if(!weights.empty())
614 {
615 for(int i = 0; i < (int)weights.size(); i++)
616 ia->setWeights(weights[i],strides[i],i);
617 }
618
619 if (uinput->hasUICoordinates()) {
620 adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
622 }
623 }
624 else if(input_type == "xpetra_crs_matrix")
625 {
626 RCP<xcrsMatrix_t> data = uinput->getUIXpetraCrsMatrix();
627 RCP<const xcrsMatrix_t> const_data = rcp_const_cast<const xcrsMatrix_t>(data);
628
629 // new adapter
630 xCM_xCM_t *ia = new xCM_xCM_t(const_data, (int)weights.size());
631 adapters.main.adapterType = AT_xCM_xCM_t;
632 adapters.main.adapter = ia;
633
634 // if we have weights set them
635 if(!weights.empty())
636 {
637 for(int i = 0; i < (int)weights.size(); i++)
638 ia->setWeights(weights[i],strides[i],i);
639 }
640
641 if (uinput->hasUICoordinates()) {
642 adapters.coordinate = getXpetraMVAdapterForInput(uinput, pCopy, comm);
644 }
645 }
646
647 if(adapters.main.adapter == nullptr)
648 {
649 std::cerr << "Input data chosen not compatible with "
650 << "XpetraCrsMatrix adapter." << std::endl;
651 return adapters;
652 }
653
654 return adapters;
655}
656
657AdapterWithTemplateName AdapterFactory::getBasicVectorAdapterForInput(
658 UserInputForTests *uinput,
659 const ParameterList &pList,
660 const RCP<const Comm<int> > &comm)
661{
662
664
665 if(!pList.isParameter("data type"))
666 {
667 std::cerr << "Input data type unspecified" << std::endl;
668 return result;
669 }
670
671 string input_type = pList.get<string>("data type");
672 if (!uinput->hasInputDataType(input_type))
673 {
674 std::cerr << "Input type:" + input_type + ", unavailable or misspelled."
675 << std::endl; // bad type
676 return result;
677 }
678
679 std::vector<const zscalar_t *> weights;
680 std::vector<int> weightStrides;
681 const zgno_t * globalIds;
682 zlno_t localCount = 0;
683
684 // get weights if any
685 // get weights if any
686 if(uinput->hasUIWeights())
687 {
688 RCP<tMVector_t> vtx_weights = uinput->getUIWeights();
689 // copy to weight
690 size_t cols = vtx_weights->getNumVectors();
691 for (size_t i = 0; i< cols; i++) {
692 weights.push_back(vtx_weights->getData(i).getRawPtr());
693 weightStrides.push_back(1);
694 }
695 }
696
697 // get vector stride
698 int stride = 1;
699 if(pList.isParameter("stride"))
700 stride = pList.get<int>("stride");
701
703
704 if(input_type == "coordinates")
705 {
706 RCP<tMVector_t> data = uinput->getUICoordinates();
707 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
708 localCount = static_cast<zlno_t>(data->getLocalLength());
709
710 // get strided data
711 std::vector<const zscalar_t *> coords;
712 std::vector<int> entry_strides;
713 InitializeVectorData(data,coords,entry_strides,stride);
714
715
716
717 if (weights.empty()) {
718 size_t dim = coords.size(); //BDD add NULL for constructor call
719 size_t push_null = 3-dim;
720 for (size_t i = 0; i < push_null; i ++) coords.push_back(NULL);
722 zlno_t(localCount),
723 globalIds,
724 coords[0],
725 coords[1],coords[2],
726 stride, stride, stride);
727 } else if (weights.size() == 1) {
728 size_t dim = coords.size(); //BDD add NULL for constructor call
729 size_t push_null = 3-dim;
730 for (size_t i = 0; i < push_null; i ++) coords.push_back(NULL);
732 zlno_t(localCount),
733 globalIds,
734 coords[0],
735 coords[1],coords[2],
736 stride, stride, stride,
737 true,
738 weights[0],
739 weightStrides[0]);
740 } else { // More than one weight per ID
742 zlno_t(localCount),
743 globalIds,
744 coords, entry_strides,
745 weights, weightStrides);
746 }
747 }
748 else if(input_type == "tpetra_vector")
749 {
750 RCP<tVector_t> data = uinput->getUITpetraVector();
751 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
752 localCount = static_cast<zlno_t>(data->getLocalLength());
753
754 // get strided data
755 std::vector<const zscalar_t *> coords;
756 std::vector<int> entry_strides;
757 InitializeVectorData(data,coords,entry_strides,stride);
758
759 if(weights.empty())
760 {
761 result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
762 coords[0], entry_strides[0]);
763 }else{
764 result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
765 coords[0], entry_strides[0],
766 true,
767 weights[0],
768 weightStrides[0]);
769
770 }
771
772 }
773 else if(input_type == "tpetra_multivector")
774 {
775 int nvec = pList.get<int>("vector_dimension");
776
777 RCP<tMVector_t> data = uinput->getUITpetraMultiVector(nvec);
778 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
779 localCount = static_cast<zlno_t>(data->getLocalLength());
780
781 // get strided data
782 std::vector<const zscalar_t *> coords;
783 std::vector<int> entry_strides;
784 InitializeVectorData(data,coords,entry_strides,stride);
785
786 result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
787 coords, entry_strides,
788 weights,weightStrides);
789
790 }
791 else if(input_type == "xpetra_vector")
792 {
793 RCP<xVector_t> data = uinput->getUIXpetraVector();
794 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
795 localCount = static_cast<zlno_t>(data->getLocalLength());
796
797 // get strided data
798 std::vector<const zscalar_t *> coords;
799 std::vector<int> entry_strides;
800 InitializeVectorData(data,coords,entry_strides,stride);
801
802 if(weights.empty())
803 {
804 result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
805 coords[0], entry_strides[0]);
806 }else{
807 result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
808 coords[0], entry_strides[0],
809 true,
810 weights[0],
811 weightStrides[0]);
812
813 }
814 }
815 else if(input_type == "xpetra_multivector")
816 {
817 int nvec = pList.get<int>("vector_dimension");
818 RCP<xMVector_t> data = uinput->getUIXpetraMultiVector(nvec);
819 globalIds = (zgno_t *)data->getMap()->getLocalElementList().getRawPtr();
820 localCount = static_cast<zlno_t>(data->getLocalLength());
821
822 // get strided data
823 std::vector<const zscalar_t *> coords;
824 std::vector<int> entry_strides;
825 InitializeVectorData(data,coords,entry_strides,stride);
826 if(comm->getRank() == 0) std::cout << "size of entry strides: " << entry_strides.size() << std::endl;
827 if(comm->getRank() == 0) std::cout << "size of coords: " << coords.size() << std::endl;
828
829 // make vector!
830 result.adapter = new Zoltan2_TestingFramework::basic_vector_adapter(localCount, globalIds,
831 coords, entry_strides,
832 weights,weightStrides);
833 }
834
835 return result;
836}
837
838template <typename T>
839void AdapterFactory::InitializeVectorData(const RCP<T> &data,
840 std::vector<const zscalar_t *> &coords,
841 std::vector<int> & strides,
842 int stride)
843{
844 // set up adapter data
845 size_t localCount = data->getLocalLength();
846 size_t nvecs = data->getNumVectors();
847 size_t vecsize = data->getNumVectors() * data->getLocalLength();
848// printf("Number of vectors by data: %zu\n", nvecs);
849 // printf("Size of data: %zu\n", vecsize);
850
851 ArrayRCP<zscalar_t> *petravectors = new ArrayRCP<zscalar_t>[nvecs];
852
853 // printf("Getting t-petra vectors...\n");
854 for (size_t i = 0; i < nvecs; i++)
855 petravectors[i] = data->getDataNonConst(i);
856
857 // debugging
858 // for (size_t i = 0; i < nvecs; i++){
859 // printf("Tpetra vector %zu: {",i);
860 //
861 // for (size_t j = 0; j < localCount; j++)
862 // {
863 // printf("%1.2g ",petravectors[i][j]);
864 // }
865 // printf("}\n");
866 // }
867
868 size_t idx = 0;
869 zscalar_t *coordarr = new zscalar_t[vecsize];
870
871 if(stride == 1 || stride != (int)nvecs)
872 {
873 for (size_t i = 0; i < nvecs; i++) {
874 for (size_t j = 0; j < localCount; j++) {
875 coordarr[idx++] = petravectors[i][j];
876 }
877 }
878 }else
879 {
880 for (size_t j = 0; j < localCount; j++) {
881 for (size_t i = 0; i < nvecs; i++) {
882 coordarr[idx++] = petravectors[i][j];
883 }
884 }
885 }
886
887 // debugging
888 // printf("Made coordarr : {");
889 // for (zlno_t i = 0; i < vecsize; i++){
890 // printf("%1.2g ",coordarr[i]);
891 // }
892 // printf("}\n");
893
894 // always build for dim 3
895 coords = std::vector<const zscalar_t *>(nvecs);
896 strides = std::vector<int>(nvecs);
897
898 for (size_t i = 0; i < nvecs; i++) {
899 if(stride == 1)
900 coords[i] = &coordarr[i*localCount];
901 else
902 coords[i] = &coordarr[i];
903
904 strides[i] = stride;
905 }
906
907 // debugging
908 // printf("Made coords...\n");
909 // for (size_t i = 0; i < nvecs; i++){
910 // const zscalar_t * tmp = coords[i];
911 // printf("coord %zu: {",i);
912 // for(size_t j = 0; j < localCount; j++)
913 // {
914 // printf("%1.2g ", tmp[j]);
915 // }
916 // printf("}\n");
917 // }
918
919 // printf("clean up coordarr and tpetravectors...\n\n\n");
920 delete [] petravectors;
921}
922
923// pamgen adapter
925AdapterFactory::getPamgenMeshAdapterForInput(UserInputForTests *uinput,
926 const ParameterList &pList,
927 const RCP<const Comm<int> > &comm)
928{
930
931#ifdef HAVE_ZOLTAN2_PAMGEN
932 if(uinput->hasPamgenMesh())
933 {
934 if(uinput->hasPamgenMesh())
935 {
936// if(comm->getRank() == 0) std::cout << "Have pamgen mesh, constructing adapter...." << std::endl;
937 result.adapter =
938 new pamgen_adapter_t(*(comm.get()), "region");
940// if(comm->getRank() == 0)
941// ia->print(0);
942 }
943 }else{
944 std::cerr << "Pamgen mesh is unavailable for PamgenMeshAdapter!"
945 << std::endl;
946 }
947
948 return result;
949#else
950 throw std::runtime_error("Pamgen input requested but Trilinos is not "
951 "built with Pamgen");
952#endif
953}
954#endif
955
956
#define SET_COORDS_INPUT_1(adapterClass)
#define SET_COORDS_INPUT_2(adapterClass)
Generate input for testing purposes.
Defines the BasicIdentifierAdapter class.
Defines the BasicVectorAdapter class.
Defines the EvaluatePartition class.
Defines the PamgenMeshAdapter class.
Defines Parameter related enumerators, declares functions.
Defines the PartitioningProblem class.
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
#define Z2_TEST_UPCAST_COORDS(adptr, TEMPLATE_ACTION)
Defines XpetraCrsGraphAdapter class.
Defines the XpetraCrsMatrixAdapter class.
Defines the XpetraMultiVectorAdapter.
Zoltan2::BaseAdapterRoot * getCoordinateAdapter() const
Zoltan2::BaseAdapterRoot * getMainAdapter() const
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.
EAdapterType getCoordinateAdapterType() const
EAdapterType getMainAdapterType() const
RCP< xcrsGraph_t > getUIXpetraCrsGraph()
RCP< tVector_t > getUITpetraVector()
RCP< tMVector_t > getUIEdgeWeights()
RCP< xMVector_t > getUIXpetraMultiVector(int nvec)
bool hasInputDataType(const string &input_type)
RCP< tMVector_t > getUITpetraMultiVector(int nvec)
RCP< tMVector_t > getUICoordinates()
RCP< tcrsGraph_t > getUITpetraCrsGraph()
RCP< xVector_t > getUIXpetraVector()
RCP< xcrsMatrix_t > getUIXpetraCrsMatrix()
RCP< tcrsMatrix_t > getUITpetraCrsMatrix()
RCP< tMVector_t > getUIWeights()
BaseAdapter defines methods required by all Adapters.
Provides access for Zoltan2 to Xpetra::CrsGraph data.
void setVertexWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to vertex weights.
void setEdgeWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to edge weights.
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
void setWeights(const scalar_t *weightVal, int stride, int idx=0)
Specify a weight for each entity of the primaryEntityType.
Zoltan2::XpetraCrsMatrixAdapter< tcrsMatrix_t, tMVector_t > xCM_tCM_t
Zoltan2::XpetraMultiVectorAdapter< tMVector_t > xMV_tMV_t
Zoltan2::XpetraCrsMatrixAdapter< xcrsMatrix_t, tMVector_t > xCM_xCM_t
Zoltan2::BasicVectorAdapter< userTypes_t > pamgen_adapter_t
Zoltan2::XpetraMultiVectorAdapter< xMVector_t > xMV_xMV_t
Zoltan2::XpetraCrsGraphAdapter< tcrsGraph_t, tMVector_t > xCG_tCG_t
Zoltan2::BasicVectorAdapter< tMVector_t > basic_vector_adapter
Zoltan2::XpetraCrsGraphAdapter< xcrsGraph_t, tMVector_t > xCG_xCG_t
Zoltan2::BasicIdentifierAdapter< userTypes_t > basic_id_t
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Zoltan2::BaseAdapterRoot * adapter