14905a7bcSToby Isaac static char help[] = "Tests LU, Cholesky, and QR factorization and MatMatSolve() for a sequential dense matrix. \n\
2c4762a1bSJed Brown For MATSEQDENSE matrix, the factorization is just a thin wrapper to LAPACK. \n\
3c4762a1bSJed Brown For MATSEQDENSECUDA, it uses cusolverDn routines \n\n";
4c4762a1bSJed Brown
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown
createMatsAndVecs(PetscInt m,PetscInt n,PetscInt nrhs,PetscBool full,Mat * _mat,Mat * _RHS,Mat * _SOLU,Vec * _x,Vec * _y,Vec * _b)7d71ae5a4SJacob Faibussowitsch static PetscErrorCode createMatsAndVecs(PetscInt m, PetscInt n, PetscInt nrhs, PetscBool full, Mat *_mat, Mat *_RHS, Mat *_SOLU, Vec *_x, Vec *_y, Vec *_b)
8d71ae5a4SJacob Faibussowitsch {
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
66fb842aefSJose E. Roman PetscCall(MatMatTransposeMult(mat, mat, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &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;
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
844905a7bcSToby Isaac }
854905a7bcSToby Isaac
main(int argc,char ** argv)86d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
87d71ae5a4SJacob Faibussowitsch {
884905a7bcSToby Isaac Mat mat, F, RHS, SOLU;
894905a7bcSToby Isaac MatInfo info;
904905a7bcSToby Isaac PetscInt m = 15, n = 10, i, j, nrhs = 2;
914905a7bcSToby Isaac Vec x, y, b, ytmp;
924905a7bcSToby Isaac IS perm;
934905a7bcSToby Isaac PetscReal norm, tol = PETSC_SMALL;
944905a7bcSToby Isaac PetscMPIInt size;
954905a7bcSToby Isaac char solver[64];
967275615fSToby Isaac PetscBool inplace, full = PETSC_FALSE, ldl = PETSC_TRUE, qr = PETSC_TRUE;
974905a7bcSToby Isaac
98327415f7SBarry Smith PetscFunctionBeginUser;
99*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
101be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
102c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERPETSC, sizeof(solver)));
1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ldl", &ldl, NULL));
1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-qr", &qr, NULL));
1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-full", &full, NULL));
1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-solver_type", solver, sizeof(solver), NULL));
1104905a7bcSToby Isaac
1119566063dSJacob Faibussowitsch PetscCall(createMatsAndVecs(n, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b));
1129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ytmp));
113c4762a1bSJed Brown
114c4762a1bSJed Brown /* Only SeqDense* support in-place factorizations and NULL permutations */
1159566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQDENSE, &inplace));
1169566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mat, &i, NULL));
1179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &j, NULL));
1189566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, i, j, 1, &perm));
119c4762a1bSJed Brown
1209566063dSJacob Faibussowitsch PetscCall(MatGetInfo(mat, MAT_LOCAL, &info));
121d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD, "matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated));
1229566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b));
123c4762a1bSJed Brown
124c4762a1bSJed Brown /* Cholesky factorization - perm and factinfo are ignored by LAPACK */
125c4762a1bSJed Brown /* in-place Cholesky */
126c4762a1bSJed Brown if (inplace) {
127c4762a1bSJed Brown Mat RHS2;
128c4762a1bSJed Brown
1296f65b6c4SPierre Jolivet if (!ldl) PetscCall(MatSetOption(mat, MAT_SPD, PETSC_TRUE));
1309566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F));
1319566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(F, perm, 0));
1329566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y));
1339566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x));
1349566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm));
13548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place Cholesky %g\n", (double)norm));
136c4762a1bSJed Brown
1379566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU));
138fb842aefSJose E. Roman PetscCall(MatMatMult(mat, SOLU, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &RHS2));
1399566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS, -1.0, RHS2, SAME_NONZERO_PATTERN));
1409566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm));
14148a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of residual for in-place Cholesky (MatMatSolve) %g\n", (double)norm));
1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS2));
144c4762a1bSJed Brown }
145c4762a1bSJed Brown
146c4762a1bSJed Brown /* out-of-place Cholesky */
1476f65b6c4SPierre Jolivet if (!ldl) PetscCall(MatSetOption(mat, MAT_SPD, PETSC_TRUE));
1489566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_CHOLESKY, &F));
1499566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, mat, perm, 0));
1509566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, mat, 0));
1519566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y));
1529566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x));
1539566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm));
15448a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place Cholesky %g\n", (double)norm));
1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
156c4762a1bSJed Brown
157c4762a1bSJed Brown /* LU factorization - perms and factinfo are ignored by LAPACK */
158c4762a1bSJed Brown i = n - 1;
1599566063dSJacob Faibussowitsch PetscCall(MatZeroRows(mat, 1, &i, -1.0, NULL, NULL));
1609566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b));
161c4762a1bSJed Brown
162c4762a1bSJed Brown /* in-place LU */
163c4762a1bSJed Brown if (inplace) {
164c4762a1bSJed Brown Mat RHS2;
165c4762a1bSJed Brown
1669566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F));
1679566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F, perm, perm, 0));
1689566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y));
1699566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x));
1709566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm));
17148a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place LU %g\n", (double)norm));
1729566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU));
173fb842aefSJose E. Roman PetscCall(MatMatMult(mat, SOLU, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &RHS2));
1749566063dSJacob Faibussowitsch PetscCall(MatAXPY(RHS, -1.0, RHS2, SAME_NONZERO_PATTERN));
1759566063dSJacob Faibussowitsch PetscCall(MatNorm(RHS, NORM_FROBENIUS, &norm));
17648a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of residual for in-place LU (MatMatSolve) %g\n", (double)norm));
1779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
1789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS2));
179c4762a1bSJed Brown }
180c4762a1bSJed Brown
181c4762a1bSJed Brown /* out-of-place LU */
1829566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_LU, &F));
1839566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, mat, perm, perm, 0));
1849566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, mat, 0));
1859566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y));
1869566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x));
1879566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm));
18848a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place LU %g\n", (double)norm));
189c4762a1bSJed Brown
190c4762a1bSJed Brown /* free space */
1919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm));
1929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
1939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
1949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS));
1959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU));
1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
1999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp));
2004905a7bcSToby Isaac
2017275615fSToby Isaac if (qr) {
202da81f932SPierre Jolivet /* setup rectangular */
2039566063dSJacob Faibussowitsch PetscCall(createMatsAndVecs(m, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b));
2049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &ytmp));
2054905a7bcSToby Isaac
2064905a7bcSToby Isaac /* QR factorization - perms and factinfo are ignored by LAPACK */
2079566063dSJacob Faibussowitsch PetscCall(MatMult(mat, x, b));
2084905a7bcSToby Isaac
2094905a7bcSToby Isaac /* in-place QR */
2104905a7bcSToby Isaac if (inplace) {
2114905a7bcSToby Isaac Mat SOLU2;
2124905a7bcSToby Isaac
2139566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, &F));
2149566063dSJacob Faibussowitsch PetscCall(MatQRFactor(F, NULL, 0));
2159566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y));
2169566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x));
2179566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm));
21848a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for in-place QR %g\n", (double)norm));
219fb842aefSJose E. Roman PetscCall(MatMatMult(mat, SOLU, MAT_REUSE_MATRIX, PETSC_DETERMINE, &RHS));
2209566063dSJacob Faibussowitsch PetscCall(MatDuplicate(SOLU, MAT_DO_NOT_COPY_VALUES, &SOLU2));
2219566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, SOLU2));
2229566063dSJacob Faibussowitsch PetscCall(MatAXPY(SOLU2, -1.0, SOLU, SAME_NONZERO_PATTERN));
2239566063dSJacob Faibussowitsch PetscCall(MatNorm(SOLU2, NORM_FROBENIUS, &norm));
22448a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error: Norm of error for in-place QR (MatMatSolve) %g\n", (double)norm));
2259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
2269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU2));
2274905a7bcSToby Isaac }
2284905a7bcSToby Isaac
2294905a7bcSToby Isaac /* out-of-place QR */
2309566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat, solver, MAT_FACTOR_QR, &F));
2319566063dSJacob Faibussowitsch PetscCall(MatQRFactorSymbolic(F, mat, NULL, NULL));
2329566063dSJacob Faibussowitsch PetscCall(MatQRFactorNumeric(F, mat, NULL));
2339566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, y));
2349566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x));
2359566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm));
23648a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place QR %g\n", (double)norm));
2374905a7bcSToby Isaac
2384905a7bcSToby Isaac if (m == n) {
2394905a7bcSToby Isaac /* out-of-place MatSolveTranspose */
2409566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat, x, b));
2419566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F, b, y));
2429566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x));
2439566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm));
24448a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Norm of error for out-of-place QR %g\n", (double)norm));
2454905a7bcSToby Isaac }
2464905a7bcSToby Isaac
2474905a7bcSToby Isaac /* free space */
2489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
2499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
2509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS));
2519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&SOLU));
2529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
2539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
2549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp));
2567275615fSToby Isaac }
2579566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
258b122ec5aSJacob Faibussowitsch return 0;
259c4762a1bSJed Brown }
260c4762a1bSJed Brown
261c4762a1bSJed Brown /*TEST
262c4762a1bSJed Brown
263c4762a1bSJed Brown test:
264c4762a1bSJed Brown
265c4762a1bSJed Brown test:
266c4762a1bSJed Brown requires: cuda
267c4762a1bSJed Brown suffix: seqdensecuda
268db0483a4SToby Isaac args: -mat_type seqdensecuda -rhs_mat_type seqdensecuda -ldl 0 -solver_type {{petsc cuda}}
2697275615fSToby Isaac output_file: output/ex1_1.out
2707275615fSToby Isaac
2717275615fSToby Isaac test:
2727275615fSToby Isaac requires: cuda
2737275615fSToby Isaac suffix: seqdensecuda_2
274db0483a4SToby Isaac args: -ldl 0 -solver_type cuda
275c4762a1bSJed Brown output_file: output/ex1_1.out
276c4762a1bSJed Brown
277c4762a1bSJed Brown test:
278c4762a1bSJed Brown requires: cuda
279c4762a1bSJed Brown suffix: seqdensecuda_seqaijcusparse
2807275615fSToby Isaac args: -mat_type seqaijcusparse -rhs_mat_type seqdensecuda -qr 0
281c4762a1bSJed Brown output_file: output/ex1_2.out
282c4762a1bSJed Brown
283c4762a1bSJed Brown test:
284c4762a1bSJed Brown requires: cuda viennacl
285c4762a1bSJed Brown suffix: seqdensecuda_seqaijviennacl
2867275615fSToby Isaac args: -mat_type seqaijviennacl -rhs_mat_type seqdensecuda -qr 0
287c4762a1bSJed Brown output_file: output/ex1_2.out
288c4762a1bSJed Brown
2894905a7bcSToby Isaac test:
2904905a7bcSToby Isaac suffix: 4
2914905a7bcSToby Isaac args: -m 10 -n 10
2924905a7bcSToby Isaac output_file: output/ex1_1.out
2934905a7bcSToby Isaac
294c4762a1bSJed Brown TEST*/
295