59 #include <Teuchos_DefaultComm.hpp>
60 #include <Teuchos_RCP.hpp>
61 #include <Teuchos_Comm.hpp>
62 #include <Teuchos_CommHelpers.hpp>
66 using Teuchos::rcp_const_cast;
67 using Teuchos::rcp_dynamic_cast;
70 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t>
ztcrsmatrix_t;
71 typedef Tpetra::RowMatrix<zscalar_t, zlno_t, zgno_t, znode_t>
ztrowmatrix_t;
75 template<
typename offset_t>
77 const zgno_t *rowIds,
const offset_t *offsets,
const zgno_t *colIds)
79 int rank = comm->getRank();
80 int nprocs = comm->getSize();
82 for (
int p=0; p < nprocs; p++){
84 std::cout << rank <<
":" << std::endl;
85 for (
zlno_t i=0; i < nrows; i++){
86 std::cout <<
" row " << rowIds[i] <<
": ";
87 for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
88 std::cout << colIds[j] <<
" ";
90 std::cout << std::endl;
101 template <
typename User>
107 RCP<const Comm<int> > comm = M.getComm();
108 int fail = 0, gfail=0;
113 if (M.getNodeNumRows()){
120 const zgno_t *rowIds=NULL, *colIds=NULL;
121 const offset_t *offsets=NULL;
130 if (nrows != M.getNodeNumRows())
136 printMatrix<offset_t>(comm, nrows, rowIds, offsets, colIds);
139 if (!fail) fail = 10;
146 int main(
int narg,
char *arg[])
148 Tpetra::ScopeGuard tscope(&narg, &arg);
149 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
151 int rank = comm->getRank();
152 int fail = 0, gfail=0;
157 RCP<UserInputForTests> uinput;
158 Teuchos::ParameterList params;
159 params.set(
"input file",
"simple");
160 params.set(
"file type",
"Chaco");
165 catch(std::exception &e){
167 std::cout << e.what() << std::endl;
172 RCP<ztcrsmatrix_t> tM = uinput->getUITpetraCrsMatrix();
173 RCP<ztrowmatrix_t> trM = rcp_dynamic_cast<
ztrowmatrix_t>(tM);
175 RCP<ztrowmatrix_t> newM;
177 size_t nrows = trM->getNodeNumRows();
189 part_t *p =
new part_t [nrows];
190 memset(p, 0,
sizeof(part_t) * nrows);
191 ArrayRCP<part_t> solnParts(p, 0, nrows,
true);
193 soln_t solution(env, comm, nWeights);
194 solution.setParts(solnParts);
200 std::cout <<
"Input adapter for Tpetra::RowMatrix" << std::endl;
202 RCP<const ztrowmatrix_t> ctrM = rcp_const_cast<
const ztrowmatrix_t>(
204 RCP<adapter_t> trMInput;
207 trMInput = rcp(
new adapter_t(ctrM));
209 catch (std::exception &e){
211 std::cout << e.what() << std::endl;
215 fail = verifyInputAdapter<ztrowmatrix_t>(*trMInput, *trM);
220 ztrowmatrix_t *mMigrate = NULL;
222 trMInput->applyPartitioningSolution(*trM, mMigrate, solution);
223 newM = rcp(mMigrate);
225 catch (std::exception &e){
227 std::cout <<
"Error caught: " << e.what() << std::endl;
233 RCP<const ztrowmatrix_t> cnewM =
234 rcp_const_cast<
const ztrowmatrix_t>(newM);
235 RCP<adapter_t> newInput;
237 newInput = rcp(
new adapter_t(cnewM));
239 catch (std::exception &e){
241 std::cout << e.what() << std::endl;
247 "Input adapter for Tpetra::RowMatrix migrated to proc 0" <<
250 fail = verifyInputAdapter<ztrowmatrix_t>(*newInput, *newM);
251 if (fail) fail += 100;
264 std::cout <<
"PASS" << std::endl;
void printFailureCode(const Comm< int > &comm, int fail)
int main(int narg, char *arg[])
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
common code used by tests
Provides access for Zoltan2 to Tpetra::RowMatrix data.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
void getRowIDsView(const gno_t *&rowIds) const
Sets pointer to this process' rows' global IDs.
size_t getLocalNumColumns() const
Returns the number of columns on this process.
Defines the TpetraRowMatrixAdapter class.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
void getCRSView(const offset_t *&offsets, const gno_t *&colIds) const
static const std::string fail
int globalFail(const Comm< int > &comm, int fail)
size_t getLocalNumRows() const
Returns the number of rows on this process.