Actual source code: ex69.c
1: static char help[] = "Tests MatCreateDenseCUDA(), MatDenseCUDAPlaceArray(), MatDenseCUDAReplaceArray(), MatDenseCUDAResetArray()\n";
3: #include <petscmat.h>
5: static PetscErrorCode MatMult_S(Mat S, Vec x, Vec y)
6: {
7: Mat A;
9: PetscFunctionBeginUser;
10: PetscCall(MatShellGetContext(S, &A));
11: PetscCall(MatMult(A, x, y));
12: PetscFunctionReturn(PETSC_SUCCESS);
13: }
15: static PetscBool test_cusparse_transgen = PETSC_FALSE;
17: static PetscErrorCode MatMultTranspose_S(Mat S, Vec x, Vec y)
18: {
19: Mat A;
21: PetscFunctionBeginUser;
22: PetscCall(MatShellGetContext(S, &A));
23: PetscCall(MatMultTranspose(A, x, y));
25: /* alternate transgen true and false to test code logic */
26: PetscCall(MatSetOption(A, MAT_FORM_EXPLICIT_TRANSPOSE, test_cusparse_transgen));
27: test_cusparse_transgen = (PetscBool)!test_cusparse_transgen;
28: PetscFunctionReturn(PETSC_SUCCESS);
29: }
31: int main(int argc, char **argv)
32: {
33: Mat A, B, C, S;
34: Vec t, v;
35: PetscScalar *vv, *aa;
36: PetscInt n = 30, k = 6, l = 0, i, Istart, Iend, nloc, bs, test = 1;
37: PetscBool flg, reset, use_shell = PETSC_FALSE;
38: VecType vtype;
40: PetscFunctionBeginUser;
41: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
42: PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
43: PetscCall(PetscOptionsGetInt(NULL, NULL, "-k", &k, NULL));
44: PetscCall(PetscOptionsGetInt(NULL, NULL, "-l", &l, NULL));
45: PetscCall(PetscOptionsGetInt(NULL, NULL, "-test", &test, NULL));
46: PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_shell", &use_shell, NULL));
47: PetscCheck(k >= 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "k %" PetscInt_FMT " must be positive", k);
48: PetscCheck(l >= 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "l %" PetscInt_FMT " must be positive", l);
49: PetscCheck(l <= k, PETSC_COMM_WORLD, PETSC_ERR_USER, "l %" PetscInt_FMT " must be smaller or equal than k %" PetscInt_FMT, l, k);
51: /* sparse matrix */
52: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
53: PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n));
54: PetscCall(MatSetType(A, MATAIJCUSPARSE));
55: PetscCall(MatSetOptionsPrefix(A, "A_"));
56: PetscCall(MatSetFromOptions(A));
57: PetscCall(MatSetUp(A));
59: /* test special case for SeqAIJCUSPARSE to generate explicit transpose (not default) */
60: PetscCall(MatSetOption(A, MAT_FORM_EXPLICIT_TRANSPOSE, PETSC_TRUE));
62: PetscCall(MatGetOwnershipRange(A, &Istart, &Iend));
63: for (i = Istart; i < Iend; i++) {
64: if (i > 0) PetscCall(MatSetValue(A, i, i - 1, -1.0, INSERT_VALUES));
65: if (i < n - 1) PetscCall(MatSetValue(A, i, i + 1, -1.0, INSERT_VALUES));
66: PetscCall(MatSetValue(A, i, i, 2.0, INSERT_VALUES));
67: }
68: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
71: /* template vector */
72: PetscCall(MatCreateVecs(A, NULL, &t));
73: PetscCall(VecGetType(t, &vtype));
75: /* long vector, contains the stacked columns of an nxk dense matrix */
76: PetscCall(VecGetLocalSize(t, &nloc));
77: PetscCall(VecGetBlockSize(t, &bs));
78: PetscCall(VecCreate(PetscObjectComm((PetscObject)t), &v));
79: PetscCall(VecSetType(v, vtype));
80: PetscCall(VecSetSizes(v, k * nloc, k * n));
81: PetscCall(VecSetBlockSize(v, bs));
82: PetscCall(VecSetRandom(v, NULL));
84: /* dense matrix that contains the columns of v */
85: PetscCall(VecCUDAGetArray(v, &vv));
87: /* test few cases for MatDenseCUDA handling pointers */
88: switch (test) {
89: case 1:
90: PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, vv, &B)); /* pass a pointer to avoid allocation of storage */
91: PetscCall(MatDenseCUDAReplaceArray(B, NULL)); /* replace with a null pointer, the value after BVRestoreMat */
92: PetscCall(MatDenseCUDAPlaceArray(B, vv + l * nloc)); /* set the actual pointer */
93: reset = PETSC_TRUE;
94: break;
95: case 2:
96: PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, NULL, &B));
97: PetscCall(MatDenseCUDAPlaceArray(B, vv + l * nloc)); /* set the actual pointer */
98: reset = PETSC_TRUE;
99: break;
100: default:
101: PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, vv + l * nloc, &B));
102: reset = PETSC_FALSE;
103: break;
104: }
106: /* Test MatMatMult */
107: if (use_shell) {
108: /* we could have called the general convertor below, but we explicit set the operations
109: ourselves to test MatProductSymbolic_X_Dense, MatProductNumeric_X_Dense code */
110: /* PetscCall(MatConvert(A,MATSHELL,MAT_INITIAL_MATRIX,&S)); */
111: PetscCall(MatCreateShell(PetscObjectComm((PetscObject)v), nloc, nloc, n, n, A, &S));
112: PetscCall(MatShellSetOperation(S, MATOP_MULT, (void (*)(void))MatMult_S));
113: PetscCall(MatShellSetOperation(S, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_S));
114: PetscCall(MatShellSetVecType(S, vtype));
115: } else {
116: PetscCall(PetscObjectReference((PetscObject)A));
117: S = A;
118: }
120: PetscCall(MatCreateDenseCUDA(PetscObjectComm((PetscObject)v), nloc, PETSC_DECIDE, n, k - l, NULL, &C));
122: /* test MatMatMult */
123: PetscCall(MatProductCreateWithMat(S, B, NULL, C));
124: PetscCall(MatProductSetType(C, MATPRODUCT_AB));
125: PetscCall(MatProductSetFromOptions(C));
126: PetscCall(MatProductSymbolic(C));
127: PetscCall(MatProductNumeric(C));
128: PetscCall(MatMatMultEqual(S, B, C, 10, &flg));
129: if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatMatMult\n"));
131: /* test MatTransposeMatMult */
132: PetscCall(MatProductCreateWithMat(S, B, NULL, C));
133: PetscCall(MatProductSetType(C, MATPRODUCT_AtB));
134: PetscCall(MatProductSetFromOptions(C));
135: PetscCall(MatProductSymbolic(C));
136: PetscCall(MatProductNumeric(C));
137: PetscCall(MatTransposeMatMultEqual(S, B, C, 10, &flg));
138: if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error MatTransposeMatMult\n"));
140: PetscCall(MatDestroy(&C));
141: PetscCall(MatDestroy(&S));
143: /* finished using B */
144: PetscCall(MatDenseGetArrayAndMemType(B, &aa, NULL));
145: PetscCheck(vv == aa - l * nloc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong array");
146: PetscCall(MatDenseRestoreArrayAndMemType(B, &aa));
147: if (reset) PetscCall(MatDenseCUDAResetArray(B));
148: PetscCall(VecCUDARestoreArray(v, &vv));
150: if (test == 1) {
151: PetscCall(MatDenseGetArrayAndMemType(B, &aa, NULL));
152: PetscCheck(!aa, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected a null pointer");
153: PetscCall(MatDenseRestoreArrayAndMemType(B, &aa));
154: }
156: /* free work space */
157: PetscCall(MatDestroy(&B));
158: PetscCall(MatDestroy(&A));
159: PetscCall(VecDestroy(&t));
160: PetscCall(VecDestroy(&v));
161: PetscCall(PetscFinalize());
162: return 0;
163: }
165: /*TEST
167: build:
168: requires: cuda
170: test:
171: requires: cuda
172: suffix: 1
173: nsize: {{1 2}}
174: args: -A_mat_type {{aij aijcusparse}} -test {{0 1 2}} -k 6 -l {{0 5}} -use_shell {{0 1}}
176: TEST*/