59 #include <Teuchos_GlobalMPISession.hpp> 60 #include <Teuchos_DefaultComm.hpp> 61 #include <Teuchos_RCP.hpp> 62 #include <Teuchos_Comm.hpp> 63 #include <Teuchos_CommHelpers.hpp> 68 using Teuchos::rcp_const_cast;
70 using Teuchos::DefaultComm;
72 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t>
tmatrix_t;
73 typedef Xpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t>
xmatrix_t;
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 (
zlno_t j=offsets[i]; j < offsets[i+1]; j++){
88 std::cout << colIds[j] <<
" ";
90 std::cout << std::endl;
99 template <
typename User>
103 RCP<const Comm<int> > comm = M.getComm();
104 int fail = 0, gfail=0;
109 if (M.getNodeNumRows()){
116 const zgno_t *rowIds=NULL, *colIds=NULL;
117 const zlno_t *offsets=NULL;
126 if (nrows != M.getNodeNumRows())
135 if (!fail) fail = 10;
141 int main(
int argc,
char *argv[])
143 Teuchos::GlobalMPISession session(&argc, &argv);
144 RCP<const Comm<int> > comm = DefaultComm<int>::getComm();
145 int rank = comm->getRank();
146 int fail = 0, gfail=0;
151 RCP<UserInputForTests> uinput;
158 catch(std::exception &e){
165 tM = uinput->getUITpetraCrsMatrix();
166 size_t nrows = tM->getNodeNumRows();
176 typedef adapter_t::part_t part_t;
178 part_t *p =
new part_t [nrows];
179 memset(p, 0,
sizeof(part_t) * nrows);
180 ArrayRCP<part_t> solnParts(p, 0, nrows,
true);
182 soln_t solution(env, comm, nWeights);
183 solution.setParts(solnParts);
188 RCP<const tmatrix_t> ctM = rcp_const_cast<
const tmatrix_t>(tM);
189 RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > tMInput;
195 catch (std::exception &e){
197 string(
"XpetraCrsMatrixAdapter ")+e.what(), 1);
201 std::cout <<
"Input adapter for Tpetra::CrsMatrix" << std::endl;
203 fail = verifyInputAdapter<tmatrix_t>(*tMInput, *tM);
210 tMInput->applyPartitioningSolution(*tM, mMigrate, solution);
211 newM = rcp(mMigrate);
213 catch (std::exception &e){
215 cout <<
"Error caught: " << e.what() << endl;
221 RCP<const tmatrix_t> cnewM = rcp_const_cast<
const tmatrix_t>(newM);
222 RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > newInput;
226 catch (std::exception &e){
228 string(
"XpetraCrsMatrixAdapter 2 ")+e.what(), 1);
233 "Input adapter for Tpetra::CrsMatrix migrated to proc 0" <<
236 fail = verifyInputAdapter<tmatrix_t>(*newInput, *newM);
237 if (fail) fail += 100;
249 RCP<xmatrix_t> xM = uinput->getUIXpetraCrsMatrix();
250 RCP<const xmatrix_t> cxM = rcp_const_cast<
const xmatrix_t>(xM);
251 RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > xMInput;
257 catch (std::exception &e){
259 string(
"XpetraCrsMatrixAdapter 3 ")+e.what(), 1);
263 std::cout <<
"Input adapter for Xpetra::CrsMatrix" << std::endl;
265 fail = verifyInputAdapter<xmatrix_t>(*xMInput, *tM);
272 xMInput->applyPartitioningSolution(*xM, mMigrate, solution);
274 catch (std::exception &e){
275 cout <<
"Error caught: " << e.what() << endl;
282 RCP<const xmatrix_t> cnewM(mMigrate);
283 RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > newInput;
288 catch (std::exception &e){
290 string(
"XpetraCrsMatrixAdapter 4 ")+e.what(), 1);
295 "Input adapter for Xpetra::CrsMatrix migrated to proc 0" <<
298 fail = verifyInputAdapter<xmatrix_t>(*newInput, *newM);
299 if (fail) fail += 100;
308 #ifdef HAVE_EPETRA_DATA_TYPES 312 RCP<ematrix_t> eM = uinput->getUIEpetraCrsMatrix();
313 RCP<const ematrix_t> ceM = rcp_const_cast<
const ematrix_t>(eM);
314 RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > eMInput;
320 catch (std::exception &e){
322 string(
"XpetraCrsMatrixAdapter 5 ")+e.what(), 1);
326 std::cout <<
"Input adapter for Epetra_CrsMatrix" << std::endl;
328 fail = verifyInputAdapter<ematrix_t>(*eMInput, *tM);
335 eMInput->applyPartitioningSolution(*eM, mMigrate, solution);
337 catch (std::exception &e){
338 cout <<
"Error caught: " << e.what() << endl;
345 RCP<const ematrix_t> cnewM(mMigrate,
true);
346 RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > newInput;
351 catch (std::exception &e){
353 string(
"XpetraCrsMatrixAdapter 6 ")+e.what(), 1);
358 "Input adapter for Epetra_CrsMatrix migrated to proc 0" <<
361 fail = verifyInputAdapter<ematrix_t>(*newInput, *newM);
362 if (fail) fail += 100;
376 std::cout <<
"PASS" << std::endl;
int globalFail(const RCP< const Comm< int > > &comm, int fail)
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
void getCRSView(const lno_t *&offsets, const gno_t *&colIds) const
size_t getLocalNumColumns() const
Returns the number of columns on this process.
common code used by tests
Defines the XpetraCrsMatrixAdapter class.
A PartitioningSolution is a solution to a partitioning problem.
size_t getLocalNumRows() const
Returns the number of rows on this process.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
void printFailureCode(const RCP< const Comm< int > > &comm, int fail)
void getRowIDsView(const gno_t *&rowIds) const
std::string testDataFilePath(".")