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 8d71ae5a4SJacob Faibussowitsch static PetscErrorCode createMatsAndVecs(PetscInt m, PetscInt n, PetscInt nrhs, PetscBool full, Mat *_mat, Mat *_RHS, Mat *_SOLU, Vec *_x, Vec *_y, Vec *_b) 9d71ae5a4SJacob Faibussowitsch { 10c4762a1bSJed Brown PetscRandom rand; 114905a7bcSToby Isaac Mat mat, RHS, SOLU; 124905a7bcSToby Isaac PetscInt rstart, rend; 134905a7bcSToby Isaac PetscInt cstart, cend; 144905a7bcSToby Isaac PetscScalar value = 1.0; 154905a7bcSToby Isaac Vec x, y, b; 16c4762a1bSJed Brown 174905a7bcSToby Isaac PetscFunctionBegin; 18c4762a1bSJed Brown /* create multiple vectors RHS and SOLU */ 199566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &RHS)); 209566063dSJacob Faibussowitsch PetscCall(MatSetSizes(RHS, PETSC_DECIDE, PETSC_DECIDE, m, nrhs)); 219566063dSJacob Faibussowitsch PetscCall(MatSetType(RHS, MATDENSE)); 229566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(RHS, "rhs_")); 239566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(RHS)); 249566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(RHS, NULL)); 25c4762a1bSJed Brown 269566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 279566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 289566063dSJacob Faibussowitsch PetscCall(MatSetRandom(RHS, rand)); 29c4762a1bSJed Brown 304905a7bcSToby Isaac if (m == n) { 319566063dSJacob Faibussowitsch PetscCall(MatDuplicate(RHS, MAT_DO_NOT_COPY_VALUES, &SOLU)); 324905a7bcSToby Isaac } else { 339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &SOLU)); 349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(SOLU, PETSC_DECIDE, PETSC_DECIDE, n, nrhs)); 359566063dSJacob Faibussowitsch PetscCall(MatSetType(SOLU, MATDENSE)); 369566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(SOLU, NULL)); 374905a7bcSToby Isaac } 389566063dSJacob Faibussowitsch PetscCall(MatSetRandom(SOLU, rand)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* create matrix */ 419566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &mat)); 429566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, PETSC_DECIDE, PETSC_DECIDE, m, n)); 439566063dSJacob Faibussowitsch PetscCall(MatSetType(mat, MATDENSE)); 449566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 459566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat)); 469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &rstart, &rend)); 479566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(mat, &cstart, &cend)); 48c4762a1bSJed Brown if (!full) { 494905a7bcSToby Isaac for (PetscInt i = rstart; i < rend; i++) { 504905a7bcSToby Isaac if (m == n) { 51c4762a1bSJed Brown value = (PetscReal)i + 1; 529566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &i, &value, INSERT_VALUES)); 534905a7bcSToby Isaac } else { 544905a7bcSToby Isaac for (PetscInt j = cstart; j < cend; j++) { 554905a7bcSToby Isaac value = ((PetscScalar)i + 1.) / (PetscSqr(i - j) + 1.); 569566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat, 1, &i, 1, &j, &value, INSERT_VALUES)); 574905a7bcSToby Isaac } 584905a7bcSToby Isaac } 59c4762a1bSJed Brown } 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 62c4762a1bSJed Brown } else { 639566063dSJacob Faibussowitsch PetscCall(MatSetRandom(mat, rand)); 644905a7bcSToby Isaac if (m == n) { 65c4762a1bSJed Brown Mat T; 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(mat, mat, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &T)); 689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 69c4762a1bSJed Brown mat = T; 70c4762a1bSJed Brown } 714905a7bcSToby Isaac } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* create single vectors */ 749566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, &x, &b)); 759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 769566063dSJacob Faibussowitsch PetscCall(VecSet(x, value)); 779566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 784905a7bcSToby Isaac *_mat = mat; 794905a7bcSToby Isaac *_RHS = RHS; 804905a7bcSToby Isaac *_SOLU = SOLU; 814905a7bcSToby Isaac *_x = x; 824905a7bcSToby Isaac *_y = y; 834905a7bcSToby Isaac *_b = b; 843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 854905a7bcSToby Isaac } 864905a7bcSToby Isaac 87d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 88d71ae5a4SJacob Faibussowitsch { 894905a7bcSToby Isaac Mat mat, F, RHS, SOLU; 904905a7bcSToby Isaac MatInfo info; 914905a7bcSToby Isaac PetscInt m = 15, n = 10, i, j, nrhs = 2; 924905a7bcSToby Isaac Vec x, y, b, ytmp; 934905a7bcSToby Isaac IS perm; 944905a7bcSToby Isaac PetscReal norm, tol = PETSC_SMALL; 954905a7bcSToby Isaac PetscMPIInt size; 964905a7bcSToby Isaac char solver[64]; 977275615fSToby Isaac PetscBool inplace, full = PETSC_FALSE, ldl = PETSC_TRUE, qr = PETSC_TRUE; 984905a7bcSToby Isaac 99327415f7SBarry Smith PetscFunctionBeginUser; 1009566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 102be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 103*c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERPETSC, sizeof(solver))); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL)); 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ldl", &ldl, NULL)); 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-qr", &qr, NULL)); 1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-full", &full, NULL)); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-solver_type", solver, sizeof(solver), NULL)); 1114905a7bcSToby Isaac 1129566063dSJacob Faibussowitsch PetscCall(createMatsAndVecs(n, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b)); 1139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ytmp)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* Only SeqDense* support in-place factorizations and NULL permutations */ 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQDENSE, &inplace)); 1179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &i, NULL)); 1189566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &j, NULL)); 1199566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, i, j, 1, &perm)); 120c4762a1bSJed Brown 1219566063dSJacob Faibussowitsch PetscCall(MatGetInfo(mat, MAT_LOCAL, &info)); 122d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD, "matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated)); 1239566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Cholesky factorization - perm and factinfo are ignored by LAPACK */ 126c4762a1bSJed Brown /* in-place Cholesky */ 127c4762a1bSJed Brown if (inplace) { 128c4762a1bSJed Brown Mat RHS2; 129c4762a1bSJed Brown 1309566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F)); 1319566063dSJacob Faibussowitsch if (!ldl) PetscCall(MatSetOption(F, MAT_SPD, PETSC_TRUE)); 1329566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(F, perm, 0)); 1339566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1349566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1359566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 13648a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place Cholesky %g\n", (double)norm)); 137c4762a1bSJed Brown 1389566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU)); 1399566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, SOLU, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &RHS2)); 1409566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS, -1.0, RHS2, SAME_NONZERO_PATTERN)); 1419566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm)); 14248a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of residual for in-place Cholesky (MatMatSolve) %g\n", (double)norm)); 1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS2)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* out-of-place Cholesky */ 1489566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_CHOLESKY, &F)); 1499566063dSJacob Faibussowitsch if (!ldl) PetscCall(MatSetOption(F, MAT_SPD, PETSC_TRUE)); 1509566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, mat, perm, 0)); 1519566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, mat, 0)); 1529566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1539566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1549566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 15548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place Cholesky %g\n", (double)norm)); 1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* LU factorization - perms and factinfo are ignored by LAPACK */ 159c4762a1bSJed Brown i = n - 1; 1609566063dSJacob Faibussowitsch PetscCall(MatZeroRows(mat, 1, &i, -1.0, NULL, NULL)); 1619566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b)); 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* in-place LU */ 164c4762a1bSJed Brown if (inplace) { 165c4762a1bSJed Brown Mat RHS2; 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F)); 1689566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F, perm, perm, 0)); 1699566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1709566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1719566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 17248a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place LU %g\n", (double)norm)); 1739566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU)); 1749566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, SOLU, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &RHS2)); 1759566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS, -1.0, RHS2, SAME_NONZERO_PATTERN)); 1769566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm)); 17748a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of residual for in-place LU (MatMatSolve) %g\n", (double)norm)); 1789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS2)); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* out-of-place LU */ 1839566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_LU, &F)); 1849566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, mat, perm, perm, 0)); 1859566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, mat, 0)); 1869566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 1879566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1889566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 18948a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place LU %g\n", (double)norm)); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* free space */ 1929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 1959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 1969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU)); 1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 2014905a7bcSToby Isaac 2027275615fSToby Isaac if (qr) { 203da81f932SPierre Jolivet /* setup rectangular */ 2049566063dSJacob Faibussowitsch PetscCall(createMatsAndVecs(m, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b)); 2059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ytmp)); 2064905a7bcSToby Isaac 2074905a7bcSToby Isaac /* QR factorization - perms and factinfo are ignored by LAPACK */ 2089566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b)); 2094905a7bcSToby Isaac 2104905a7bcSToby Isaac /* in-place QR */ 2114905a7bcSToby Isaac if (inplace) { 2124905a7bcSToby Isaac Mat SOLU2; 2134905a7bcSToby Isaac 2149566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F)); 2159566063dSJacob Faibussowitsch PetscCall(MatQRFactor(F, NULL, 0)); 2169566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 2179566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 2189566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 21948a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place QR %g\n", (double)norm)); 2209566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, SOLU, MAT_REUSE_MATRIX, PETSC_DEFAULT, &RHS)); 2219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(SOLU, MAT_DO_NOT_COPY_VALUES, &SOLU2)); 2229566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU2)); 2239566063dSJacob Faibussowitsch PetscCall(MatAXPY(SOLU2, -1.0, SOLU, SAME_NONZERO_PATTERN)); 2249566063dSJacob Faibussowitsch PetscCall(MatNorm(SOLU2, NORM_FROBENIUS, &norm)); 22548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of error for in-place QR (MatMatSolve) %g\n", (double)norm)); 2269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU2)); 2284905a7bcSToby Isaac } 2294905a7bcSToby Isaac 2304905a7bcSToby Isaac /* out-of-place QR */ 2319566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_QR, &F)); 2329566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(F, mat, NULL, NULL)); 2339566063dSJacob Faibussowitsch PetscCall(MatQRFactorNumeric(F, mat, NULL)); 2349566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y)); 2359566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 2369566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 23748a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place QR %g\n", (double)norm)); 2384905a7bcSToby Isaac 2394905a7bcSToby Isaac if (m == n) { 2404905a7bcSToby Isaac /* out-of-place MatSolveTranspose */ 2419566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat, x, b)); 2429566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F, b, y)); 2439566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 2449566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 24548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place QR %g\n", (double)norm)); 2464905a7bcSToby Isaac } 2474905a7bcSToby Isaac 2484905a7bcSToby Isaac /* free space */ 2499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 2509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 2519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 2529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU)); 2539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 2577275615fSToby Isaac } 2589566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 259b122ec5aSJacob Faibussowitsch return 0; 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262c4762a1bSJed Brown /*TEST 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown 266c4762a1bSJed Brown test: 267c4762a1bSJed Brown requires: cuda 268c4762a1bSJed Brown suffix: seqdensecuda 269db0483a4SToby Isaac args: -mat_type seqdensecuda -rhs_mat_type seqdensecuda -ldl 0 -solver_type {{petsc cuda}} 2707275615fSToby Isaac output_file: output/ex1_1.out 2717275615fSToby Isaac 2727275615fSToby Isaac test: 2737275615fSToby Isaac requires: cuda 2747275615fSToby Isaac suffix: seqdensecuda_2 275db0483a4SToby Isaac args: -ldl 0 -solver_type cuda 276c4762a1bSJed Brown output_file: output/ex1_1.out 277c4762a1bSJed Brown 278c4762a1bSJed Brown test: 279c4762a1bSJed Brown requires: cuda 280c4762a1bSJed Brown suffix: seqdensecuda_seqaijcusparse 2817275615fSToby Isaac args: -mat_type seqaijcusparse -rhs_mat_type seqdensecuda -qr 0 282c4762a1bSJed Brown output_file: output/ex1_2.out 283c4762a1bSJed Brown 284c4762a1bSJed Brown test: 285c4762a1bSJed Brown requires: cuda viennacl 286c4762a1bSJed Brown suffix: seqdensecuda_seqaijviennacl 2877275615fSToby Isaac args: -mat_type seqaijviennacl -rhs_mat_type seqdensecuda -qr 0 288c4762a1bSJed Brown output_file: output/ex1_2.out 289c4762a1bSJed Brown 2904905a7bcSToby Isaac test: 2914905a7bcSToby Isaac suffix: 4 2924905a7bcSToby Isaac args: -m 10 -n 10 2934905a7bcSToby Isaac output_file: output/ex1_1.out 2944905a7bcSToby Isaac 295c4762a1bSJed Brown TEST*/ 296