1c4762a1bSJed Brown 24905a7bcSToby Isaac static char help[] = "Tests LU, Cholesky, and QR factorization and MatMatSolve() for a sequential dense matrix. \n\ 3c4762a1bSJed Brown For MATSEQDENSE matrix, the factorization is just a thin wrapper to LAPACK. \n\ 4c4762a1bSJed Brown For MATSEQDENSECUDA, it uses cusolverDn routines \n\n"; 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petscmat.h> 7c4762a1bSJed Brown 89371c9d4SSatish Balay static PetscErrorCode createMatsAndVecs(PetscInt m, PetscInt n, PetscInt nrhs, PetscBool full, Mat *_mat, Mat *_RHS, Mat *_SOLU, Vec *_x, Vec *_y, Vec *_b) { 9c4762a1bSJed Brown PetscRandom rand; 104905a7bcSToby Isaac Mat mat, RHS, SOLU; 114905a7bcSToby Isaac PetscInt rstart, rend; 124905a7bcSToby Isaac PetscInt cstart, cend; 134905a7bcSToby Isaac PetscScalar value = 1.0; 144905a7bcSToby Isaac Vec x, y, b; 15c4762a1bSJed Brown 164905a7bcSToby Isaac PetscFunctionBegin; 17c4762a1bSJed Brown /* create multiple vectors RHS and SOLU */ 189566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &RHS)); 199566063dSJacob Faibussowitsch PetscCall(MatSetSizes(RHS, PETSC_DECIDE, PETSC_DECIDE, m, nrhs)); 209566063dSJacob Faibussowitsch PetscCall(MatSetType(RHS, MATDENSE)); 219566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(RHS, "rhs_")); 229566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(RHS)); 239566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(RHS, NULL)); 24c4762a1bSJed Brown 259566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 269566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 279566063dSJacob Faibussowitsch PetscCall(MatSetRandom(RHS, rand)); 28c4762a1bSJed Brown 294905a7bcSToby Isaac if (m == n) { 309566063dSJacob Faibussowitsch PetscCall(MatDuplicate(RHS, MAT_DO_NOT_COPY_VALUES, &SOLU)); 314905a7bcSToby Isaac } else { 329566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &SOLU)); 339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(SOLU, PETSC_DECIDE, PETSC_DECIDE, n, nrhs)); 349566063dSJacob Faibussowitsch PetscCall(MatSetType(SOLU, MATDENSE)); 359566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(SOLU, NULL)); 364905a7bcSToby Isaac } 379566063dSJacob Faibussowitsch PetscCall(MatSetRandom(SOLU, rand)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* create matrix */ 409566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &mat)); 419566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, PETSC_DECIDE, PETSC_DECIDE, m, n)); 429566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, MATDENSE)); 439566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 449566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat)); 459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &rstart, &rend)); 469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(mat, &cstart, &cend)); 47c4762a1bSJed Brown if (!full) { 484905a7bcSToby Isaac for (PetscInt i = rstart; i < rend; i++) { 494905a7bcSToby Isaac if (m == n) { 50c4762a1bSJed Brown value = (PetscReal)i + 1; 519566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &i, &value, INSERT_VALUES)); 524905a7bcSToby Isaac } else { 534905a7bcSToby Isaac for (PetscInt j = cstart; j < cend; j++) { 544905a7bcSToby Isaac value = ((PetscScalar)i + 1.) / (PetscSqr(i - j) + 1.); 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &j, &value, INSERT_VALUES)); 564905a7bcSToby Isaac } 574905a7bcSToby Isaac } 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 61c4762a1bSJed Brown } else { 629566063dSJacob Faibussowitsch PetscCall(MatSetRandom(mat, rand)); 634905a7bcSToby Isaac if (m == n) { 64c4762a1bSJed Brown Mat T; 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(mat, mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &T)); 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 68c4762a1bSJed Brown mat = T; 69c4762a1bSJed Brown } 704905a7bcSToby Isaac } 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* create single vectors */ 739566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, &x, &b)); 749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 759566063dSJacob Faibussowitsch PetscCall(VecSet(x, value)); 769566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 774905a7bcSToby Isaac *_mat = mat; 784905a7bcSToby Isaac *_RHS = RHS; 794905a7bcSToby Isaac *_SOLU = SOLU; 804905a7bcSToby Isaac *_x = x; 814905a7bcSToby Isaac *_y = y; 824905a7bcSToby Isaac *_b = b; 834905a7bcSToby Isaac PetscFunctionReturn(0); 844905a7bcSToby Isaac } 854905a7bcSToby Isaac 869371c9d4SSatish Balay int main(int argc, char **argv) { 874905a7bcSToby Isaac Mat mat, F, RHS, SOLU; 884905a7bcSToby Isaac MatInfo info; 894905a7bcSToby Isaac PetscInt m = 15, n = 10, i, j, nrhs = 2; 904905a7bcSToby Isaac Vec x, y, b, ytmp; 914905a7bcSToby Isaac IS perm; 924905a7bcSToby Isaac PetscReal norm, tol = PETSC_SMALL; 934905a7bcSToby Isaac PetscMPIInt size; 944905a7bcSToby Isaac char solver[64]; 957275615fSToby Isaac PetscBool inplace, full = PETSC_FALSE, ldl = PETSC_TRUE, qr = PETSC_TRUE; 964905a7bcSToby Isaac 97327415f7SBarry Smith PetscFunctionBeginUser; 989566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 100be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 1019566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(solver, "petsc")); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ldl", &ldl, NULL)); 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-qr", &qr, NULL)); 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-full", &full, NULL)); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-solver_type", solver, sizeof(solver), NULL)); 1094905a7bcSToby Isaac 1109566063dSJacob Faibussowitsch PetscCall(createMatsAndVecs(n, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b)); 1119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ytmp)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Only SeqDense* support in-place factorizations and NULL permutations */ 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQDENSE, &inplace)); 1159566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &i, NULL)); 1169566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &j, NULL)); 1179566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, i, j, 1, &perm)); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(MatGetInfo(mat, MAT_LOCAL, &info)); 120d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD, "matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated)); 1219566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* Cholesky factorization - perm and factinfo are ignored by LAPACK */ 124c4762a1bSJed Brown /* in-place Cholesky */ 125c4762a1bSJed Brown if (inplace) { 126c4762a1bSJed Brown Mat RHS2; 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F)); 1299566063dSJacob Faibussowitsch if (!ldl) PetscCall(MatSetOption(F, MAT_SPD, PETSC_TRUE)); 1309566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(F, perm, 0)); 1319566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1329566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1339566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 134*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place Cholesky %g\n", (double)norm)); 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU)); 1379566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, SOLU, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &RHS2)); 1389566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS, -1.0, RHS2, SAME_NONZERO_PATTERN)); 1399566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm)); 140*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of residual for in-place Cholesky (MatMatSolve) %g\n", (double)norm)); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS2)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* out-of-place Cholesky */ 1469566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_CHOLESKY, &F)); 1479566063dSJacob Faibussowitsch if (!ldl) PetscCall(MatSetOption(F, MAT_SPD, PETSC_TRUE)); 1489566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, mat, perm, 0)); 1499566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, mat, 0)); 1509566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1519566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1529566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 153*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place Cholesky %g\n", (double)norm)); 1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* LU factorization - perms and factinfo are ignored by LAPACK */ 157c4762a1bSJed Brown i = n - 1; 1589566063dSJacob Faibussowitsch PetscCall(MatZeroRows(mat, 1, &i, -1.0, NULL, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* in-place LU */ 162c4762a1bSJed Brown if (inplace) { 163c4762a1bSJed Brown Mat RHS2; 164c4762a1bSJed Brown 1659566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F)); 1669566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F, perm, perm, 0)); 1679566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1689566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1699566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 170*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place LU %g\n", (double)norm)); 1719566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU)); 1729566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, SOLU, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &RHS2)); 1739566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS, -1.0, RHS2, SAME_NONZERO_PATTERN)); 1749566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm)); 175*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of residual for in-place LU (MatMatSolve) %g\n", (double)norm)); 1769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS2)); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* out-of-place LU */ 1819566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_LU, &F)); 1829566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, mat, perm, perm, 0)); 1839566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, mat, 0)); 1849566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1859566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1869566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 187*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place LU %g\n", (double)norm)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* free space */ 1909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 1939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 1949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU)); 1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 1994905a7bcSToby Isaac 2007275615fSToby Isaac if (qr) { 2014905a7bcSToby Isaac /* setup rectanglar */ 2029566063dSJacob Faibussowitsch PetscCall(createMatsAndVecs(m, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b)); 2039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ytmp)); 2044905a7bcSToby Isaac 2054905a7bcSToby Isaac /* QR factorization - perms and factinfo are ignored by LAPACK */ 2069566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b)); 2074905a7bcSToby Isaac 2084905a7bcSToby Isaac /* in-place QR */ 2094905a7bcSToby Isaac if (inplace) { 2104905a7bcSToby Isaac Mat SOLU2; 2114905a7bcSToby Isaac 2129566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F)); 2139566063dSJacob Faibussowitsch PetscCall(MatQRFactor(F, NULL, 0)); 2149566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 2159566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 2169566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 217*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place QR %g\n", (double)norm)); 2189566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, SOLU, MAT_REUSE_MATRIX, PETSC_DEFAULT, &RHS)); 2199566063dSJacob Faibussowitsch PetscCall(MatDuplicate(SOLU, MAT_DO_NOT_COPY_VALUES, &SOLU2)); 2209566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU2)); 2219566063dSJacob Faibussowitsch PetscCall(MatAXPY(SOLU2, -1.0, SOLU, SAME_NONZERO_PATTERN)); 2229566063dSJacob Faibussowitsch PetscCall(MatNorm(SOLU2, NORM_FROBENIUS, &norm)); 223*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of error for in-place QR (MatMatSolve) %g\n", (double)norm)); 2249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 2259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU2)); 2264905a7bcSToby Isaac } 2274905a7bcSToby Isaac 2284905a7bcSToby Isaac /* out-of-place QR */ 2299566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_QR, &F)); 2309566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(F, mat, NULL, NULL)); 2319566063dSJacob Faibussowitsch PetscCall(MatQRFactorNumeric(F, mat, NULL)); 2329566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 2339566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 2349566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 235*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place QR %g\n", (double)norm)); 2364905a7bcSToby Isaac 2374905a7bcSToby Isaac if (m == n) { 2384905a7bcSToby Isaac /* out-of-place MatSolveTranspose */ 2399566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat, x, b)); 2409566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F, b, y)); 2419566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 2429566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 243*48a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place QR %g\n", (double)norm)); 2444905a7bcSToby Isaac } 2454905a7bcSToby Isaac 2464905a7bcSToby Isaac /* free space */ 2479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 2489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 2499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 2509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU)); 2519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 2557275615fSToby Isaac } 2569566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 257b122ec5aSJacob Faibussowitsch return 0; 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown /*TEST 261c4762a1bSJed Brown 262c4762a1bSJed Brown test: 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown requires: cuda 266c4762a1bSJed Brown suffix: seqdensecuda 267db0483a4SToby Isaac args: -mat_type seqdensecuda -rhs_mat_type seqdensecuda -ldl 0 -solver_type {{petsc cuda}} 2687275615fSToby Isaac output_file: output/ex1_1.out 2697275615fSToby Isaac 2707275615fSToby Isaac test: 2717275615fSToby Isaac requires: cuda 2727275615fSToby Isaac suffix: seqdensecuda_2 273db0483a4SToby Isaac args: -ldl 0 -solver_type cuda 274c4762a1bSJed Brown output_file: output/ex1_1.out 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown requires: cuda 278c4762a1bSJed Brown suffix: seqdensecuda_seqaijcusparse 2797275615fSToby Isaac args: -mat_type seqaijcusparse -rhs_mat_type seqdensecuda -qr 0 280c4762a1bSJed Brown output_file: output/ex1_2.out 281c4762a1bSJed Brown 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown requires: cuda viennacl 284c4762a1bSJed Brown suffix: seqdensecuda_seqaijviennacl 2857275615fSToby Isaac args: -mat_type seqaijviennacl -rhs_mat_type seqdensecuda -qr 0 286c4762a1bSJed Brown output_file: output/ex1_2.out 287c4762a1bSJed Brown 2884905a7bcSToby Isaac test: 2894905a7bcSToby Isaac suffix: 4 2904905a7bcSToby Isaac args: -m 10 -n 10 2914905a7bcSToby Isaac output_file: output/ex1_1.out 2924905a7bcSToby Isaac 293c4762a1bSJed Brown TEST*/ 294