Actual source code: ex70.c
1: #include <petscmat.h>
3: static char help[] = "Tests MatMat operations with MAT_REUSE_MATRIX and already allocated dense result.\n\n";
5: static PetscScalar MAGIC_NUMBER = 12345;
7: static PetscErrorCode CheckLocal(Mat A, Mat B, PetscScalar *a, PetscScalar *b)
8: {
9: PetscBool wA = PETSC_FALSE, wB = PETSC_FALSE;
10: PetscBool wAv = PETSC_FALSE, wBv = PETSC_FALSE;
11: PetscInt lda, i, j, m, n;
13: PetscFunctionBegin;
14: if (a) {
15: const PetscScalar *Aa;
16: PetscCall(MatDenseGetArrayRead(A, &Aa));
17: wA = (PetscBool)(a != Aa);
18: PetscCall(MatDenseGetLDA(A, &lda));
19: PetscCall(MatGetLocalSize(A, &m, &n));
20: for (j = 0; j < n; j++) {
21: for (i = m; i < lda; i++) {
22: if (Aa[j * lda + i] != MAGIC_NUMBER) wAv = PETSC_TRUE;
23: }
24: }
25: PetscCall(MatDenseRestoreArrayRead(A, &Aa));
26: }
27: if (b) {
28: const PetscScalar *Bb;
29: PetscCall(MatDenseGetArrayRead(B, &Bb));
30: wB = (PetscBool)(b != Bb);
31: PetscCall(MatDenseGetLDA(B, &lda));
32: PetscCall(MatGetLocalSize(B, &m, &n));
33: for (j = 0; j < n; j++) {
34: for (i = m; i < lda; i++) {
35: if (Bb[j * lda + i] != MAGIC_NUMBER) wBv = PETSC_TRUE;
36: }
37: }
38: PetscCall(MatDenseRestoreArrayRead(B, &Bb));
39: }
40: PetscCheck(!wA && !wB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong array in first Mat? %d, Wrong array in second Mat? %d", wA, wB);
41: PetscCheck(!wAv && !wBv, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong data in first Mat? %d, Wrong data in second Mat? %d", wAv, wBv);
42: PetscFunctionReturn(PETSC_SUCCESS);
43: }
45: typedef struct {
46: Mat A;
47: Mat P;
48: Mat R;
49: } proj_data;
51: PetscErrorCode proj_destroy(PetscCtxRt ctx)
52: {
53: proj_data *userdata = *(proj_data **)ctx;
55: PetscFunctionBegin;
56: PetscCheck(userdata, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing userdata");
57: PetscCall(MatDestroy(&userdata->A));
58: PetscCall(MatDestroy(&userdata->P));
59: PetscCall(MatDestroy(&userdata->R));
60: PetscCall(PetscFree(userdata));
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
64: PetscErrorCode proj_mult(Mat S, Vec X, Vec Y)
65: {
66: Mat A, R, P;
67: Vec Ax, Ay;
68: Vec Px, Py;
69: proj_data *userdata;
71: PetscFunctionBegin;
72: PetscCall(MatShellGetContext(S, &userdata));
73: PetscCheck(userdata, PetscObjectComm((PetscObject)S), PETSC_ERR_PLIB, "Missing userdata");
74: A = userdata->A;
75: R = userdata->R;
76: P = userdata->P;
77: PetscCheck(A, PetscObjectComm((PetscObject)S), PETSC_ERR_PLIB, "Missing matrix");
78: PetscCheck(R || P, PetscObjectComm((PetscObject)S), PETSC_ERR_PLIB, "Missing projectors");
79: PetscCheck(!R || !P, PetscObjectComm((PetscObject)S), PETSC_ERR_PLIB, "Both projectors");
80: PetscCall(MatCreateVecs(A, &Ax, &Ay));
81: if (R) PetscCall(MatCreateVecs(R, &Py, &Px));
82: else PetscCall(MatCreateVecs(P, &Px, &Py));
83: PetscCall(VecCopy(X, Px));
84: if (P) PetscCall(MatMult(P, Px, Py));
85: else PetscCall(MatMultTranspose(R, Px, Py));
86: PetscCall(VecCopy(Py, Ax));
87: PetscCall(MatMult(A, Ax, Ay));
88: PetscCall(VecCopy(Ay, Py));
89: if (P) PetscCall(MatMultTranspose(P, Py, Px));
90: else PetscCall(MatMult(R, Py, Px));
91: PetscCall(VecCopy(Px, Y));
92: PetscCall(VecDestroy(&Px));
93: PetscCall(VecDestroy(&Py));
94: PetscCall(VecDestroy(&Ax));
95: PetscCall(VecDestroy(&Ay));
96: PetscFunctionReturn(PETSC_SUCCESS);
97: }
99: PetscErrorCode MyPtShellPMultSymbolic(Mat S, Mat P, Mat PtAP, void **ctx)
100: {
101: proj_data *userdata;
103: PetscFunctionBegin;
104: PetscCall(PetscNew(&userdata));
105: PetscCall(MatShellSetContext(PtAP, userdata));
106: *ctx = (void *)userdata;
107: PetscFunctionReturn(PETSC_SUCCESS);
108: }
110: PetscErrorCode MyPtShellPMultNumeric(Mat S, Mat P, Mat PtAP, PetscCtx ctx)
111: {
112: Mat A;
113: proj_data *userdata = (proj_data *)ctx;
115: PetscFunctionBegin;
116: PetscCall(MatShellGetContext(S, &A));
117: PetscCall(PetscObjectReference((PetscObject)A));
118: PetscCall(PetscObjectReference((PetscObject)P));
119: PetscCall(MatDestroy(&userdata->A));
120: PetscCall(MatDestroy(&userdata->P));
121: PetscCall(MatDestroy(&userdata->R));
122: userdata->A = A;
123: userdata->P = P;
124: PetscCall(MatShellSetOperation(PtAP, MATOP_MULT, (PetscErrorCodeFn *)proj_mult));
125: PetscCall(MatSetUp(PtAP));
126: PetscCall(MatAssemblyBegin(PtAP, MAT_FINAL_ASSEMBLY));
127: PetscCall(MatAssemblyEnd(PtAP, MAT_FINAL_ASSEMBLY));
128: PetscFunctionReturn(PETSC_SUCCESS);
129: }
131: PetscErrorCode MyRShellRtMultSymbolic(Mat S, Mat R, Mat RARt, void **ctx)
132: {
133: proj_data *userdata;
135: PetscFunctionBegin;
136: PetscCall(PetscNew(&userdata));
137: PetscCall(MatShellSetContext(RARt, userdata));
138: *ctx = (void *)userdata;
139: PetscFunctionReturn(PETSC_SUCCESS);
140: }
142: PetscErrorCode MyRShellRtMultNumeric(Mat S, Mat R, Mat RARt, PetscCtx ctx)
143: {
144: Mat A;
145: proj_data *userdata = (proj_data *)ctx;
147: PetscFunctionBegin;
148: PetscCall(MatShellGetContext(S, &A));
149: PetscCall(PetscObjectReference((PetscObject)A));
150: PetscCall(PetscObjectReference((PetscObject)R));
151: PetscCall(MatDestroy(&userdata->A));
152: PetscCall(MatDestroy(&userdata->P));
153: PetscCall(MatDestroy(&userdata->R));
154: userdata->A = A;
155: userdata->R = R;
156: PetscCall(MatShellSetOperation(RARt, MATOP_MULT, (PetscErrorCodeFn *)proj_mult));
157: PetscCall(MatSetUp(RARt));
158: PetscCall(MatAssemblyBegin(RARt, MAT_FINAL_ASSEMBLY));
159: PetscCall(MatAssemblyEnd(RARt, MAT_FINAL_ASSEMBLY));
160: PetscFunctionReturn(PETSC_SUCCESS);
161: }
163: PetscErrorCode MyMatShellMatMultNumeric(Mat S, Mat B, Mat C, PetscCtx ctx)
164: {
165: Mat A;
167: PetscFunctionBegin;
168: PetscCall(MatShellGetContext(S, &A));
169: PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &C));
170: PetscFunctionReturn(PETSC_SUCCESS);
171: }
173: PetscErrorCode MyMatTransposeShellMatMultNumeric(Mat S, Mat B, Mat C, PetscCtx ctx)
174: {
175: Mat A;
177: PetscFunctionBegin;
178: PetscCall(MatShellGetContext(S, &A));
179: PetscCall(MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &C));
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: PetscErrorCode MyMatShellMatTransposeMultNumeric(Mat S, Mat B, Mat C, PetscCtx ctx)
184: {
185: Mat A;
187: PetscFunctionBegin;
188: PetscCall(MatShellGetContext(S, &A));
189: PetscCall(MatMatTransposeMult(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &C));
190: PetscFunctionReturn(PETSC_SUCCESS);
191: }
193: int main(int argc, char **args)
194: {
195: Mat X, B, A, Bt, T, T2, PtAP = NULL, RARt = NULL, R = NULL;
196: Vec r, l, rs, ls;
197: PetscInt m, n, k, M = 10, N = 10, K = 5, ldx = 3, ldb = 5, ldr = 4;
198: const char *deft = MATAIJ;
199: char mattype[256];
200: PetscBool flg, symm = PETSC_FALSE, testtt = PETSC_TRUE, testnest = PETSC_TRUE, testtranspose = PETSC_TRUE, testcircular = PETSC_FALSE, local = PETSC_TRUE;
201: PetscBool testhtranspose = PETSC_FALSE; /* Hermitian transpose is not handled correctly and generates an error */
202: PetscBool xgpu = PETSC_FALSE, bgpu = PETSC_FALSE, testshellops = PETSC_FALSE, testproj = PETSC_TRUE, testrart = PETSC_TRUE, testmatmatt = PETSC_TRUE, testmattmat = PETSC_TRUE;
203: PetscScalar *dataX = NULL, *dataB = NULL, *dataR = NULL, *dataBt = NULL;
204: PetscScalar *aX, *aB, *aBt;
205: PetscReal err;
207: PetscFunctionBeginUser;
208: PetscCall(PetscInitialize(&argc, &args, NULL, help));
209: PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
210: PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
211: PetscCall(PetscOptionsGetInt(NULL, NULL, "-K", &K, NULL));
212: PetscCall(PetscOptionsGetBool(NULL, NULL, "-symm", &symm, NULL));
213: PetscCall(PetscOptionsGetBool(NULL, NULL, "-local", &local, NULL));
214: PetscCall(PetscOptionsGetInt(NULL, NULL, "-ldx", &ldx, NULL));
215: PetscCall(PetscOptionsGetInt(NULL, NULL, "-ldb", &ldb, NULL));
216: PetscCall(PetscOptionsGetInt(NULL, NULL, "-ldr", &ldr, NULL));
217: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testtranspose", &testtranspose, NULL));
218: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testnest", &testnest, NULL));
219: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testtt", &testtt, NULL));
220: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testcircular", &testcircular, NULL));
221: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testshellops", &testshellops, NULL));
222: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testproj", &testproj, NULL));
223: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testrart", &testrart, NULL));
224: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testmatmatt", &testmatmatt, NULL));
225: PetscCall(PetscOptionsGetBool(NULL, NULL, "-testmattmat", &testmattmat, NULL));
226: PetscCall(PetscOptionsGetBool(NULL, NULL, "-xgpu", &xgpu, NULL));
227: PetscCall(PetscOptionsGetBool(NULL, NULL, "-bgpu", &bgpu, NULL));
228: PetscCall(PetscOptionsGetScalar(NULL, NULL, "-magic_number", &MAGIC_NUMBER, NULL));
229: if (M != N) testproj = PETSC_FALSE;
231: PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
232: PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
233: PetscCall(MatSetType(A, MATAIJ));
234: PetscCall(MatSeqAIJSetPreallocation(A, PETSC_DEFAULT, NULL));
235: PetscCall(MatMPIAIJSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
236: PetscCall(MatSetUp(A));
237: PetscCall(MatSetRandom(A, NULL));
238: if (M == N && symm) {
239: Mat AT;
241: PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &AT));
242: PetscCall(MatAXPY(A, 1.0, AT, DIFFERENT_NONZERO_PATTERN));
243: PetscCall(MatDestroy(&AT));
244: PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
245: }
246: PetscCall(MatViewFromOptions(A, NULL, "-A_init_view"));
247: PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "");
248: PetscCall(PetscOptionsFList("-A_mat_type", "Matrix type", "MatSetType", MatList, deft, mattype, 256, &flg));
249: PetscOptionsEnd();
250: if (flg) {
251: Mat A2;
253: PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A2));
254: PetscCall(MatConvert(A, mattype, MAT_INPLACE_MATRIX, &A));
255: PetscCall(MatMultEqual(A, A2, 10, &flg));
256: if (!flg) {
257: Mat AE, A2E;
259: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with convert\n"));
260: PetscCall(MatComputeOperator(A, MATDENSE, &AE));
261: PetscCall(MatComputeOperator(A2, MATDENSE, &A2E));
262: PetscCall(MatView(AE, NULL));
263: PetscCall(MatView(A2E, NULL));
264: PetscCall(MatAXPY(A2E, -1.0, A, SAME_NONZERO_PATTERN));
265: PetscCall(MatView(A2E, NULL));
266: PetscCall(MatDestroy(&A2E));
267: PetscCall(MatDestroy(&AE));
268: }
269: PetscCall(MatDestroy(&A2));
270: }
271: PetscCall(MatViewFromOptions(A, NULL, "-A_view"));
273: PetscCall(MatGetLocalSize(A, &m, &n));
274: if (local) {
275: PetscInt i;
277: PetscCall(PetscMalloc1((m + ldx) * K, &dataX));
278: PetscCall(PetscMalloc1((n + ldb) * K, &dataB));
279: for (i = 0; i < (m + ldx) * K; i++) dataX[i] = MAGIC_NUMBER;
280: for (i = 0; i < (n + ldb) * K; i++) dataB[i] = MAGIC_NUMBER;
281: }
282: PetscCall(MatCreateDense(PETSC_COMM_WORLD, n, PETSC_DECIDE, N, K, dataB, &B));
283: PetscCall(MatCreateDense(PETSC_COMM_WORLD, m, PETSC_DECIDE, M, K, dataX, &X));
284: if (local) {
285: PetscCall(MatDenseSetLDA(X, m + ldx));
286: PetscCall(MatDenseSetLDA(B, n + ldb));
287: }
288: PetscCall(MatGetLocalSize(B, NULL, &k));
289: if (local) {
290: PetscInt i;
292: PetscCall(PetscMalloc1((k + ldr) * N, &dataBt));
293: for (i = 0; i < (k + ldr) * N; i++) dataBt[i] = MAGIC_NUMBER;
294: }
295: PetscCall(MatCreateDense(PETSC_COMM_WORLD, k, n, K, N, dataBt, &Bt));
296: if (local) PetscCall(MatDenseSetLDA(Bt, k + ldr));
298: /* store pointer to dense data for testing */
299: PetscCall(MatDenseGetArrayRead(B, (const PetscScalar **)&dataB));
300: PetscCall(MatDenseGetArrayRead(X, (const PetscScalar **)&dataX));
301: PetscCall(MatDenseGetArrayRead(Bt, (const PetscScalar **)&dataBt));
302: aX = dataX;
303: aB = dataB;
304: aBt = dataBt;
305: PetscCall(MatDenseRestoreArrayRead(Bt, (const PetscScalar **)&dataBt));
306: PetscCall(MatDenseRestoreArrayRead(B, (const PetscScalar **)&dataB));
307: PetscCall(MatDenseRestoreArrayRead(X, (const PetscScalar **)&dataX));
308: if (local) {
309: dataX = aX;
310: dataB = aB;
311: dataBt = aBt;
312: }
314: PetscCall(MatSetRandom(X, NULL));
315: PetscCall(MatSetRandom(B, NULL));
316: PetscCall(MatSetRandom(Bt, NULL));
317: PetscCall(CheckLocal(X, NULL, aX, NULL));
318: PetscCall(CheckLocal(Bt, B, aBt, aB));
320: /* convert to CUDA if needed */
321: if (bgpu) {
322: PetscCall(MatConvert(B, MATDENSECUDA, MAT_INPLACE_MATRIX, &B));
323: PetscCall(MatConvert(Bt, MATDENSECUDA, MAT_INPLACE_MATRIX, &Bt));
324: }
325: if (xgpu) PetscCall(MatConvert(X, MATDENSECUDA, MAT_INPLACE_MATRIX, &X));
326: PetscCall(CheckLocal(B, X, aB, aX));
328: /* Test MatDenseGetSubMatrix */
329: {
330: Mat B2, T3, T4;
332: PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
333: PetscCall(MatDuplicate(B, MAT_DO_NOT_COPY_VALUES, &T4));
334: PetscCall(MatSetRandom(T4, NULL));
335: PetscCall(MatAXPY(B2, 1.0, T4, SAME_NONZERO_PATTERN));
336: PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, PETSC_DECIDE, PetscMin(1, K - 1), PetscMin(2, K), &T));
337: PetscCall(MatDenseGetSubMatrix(T4, PETSC_DECIDE, PETSC_DECIDE, PetscMin(1, K - 1), PetscMin(2, K), &T2));
338: PetscCall(MatDenseGetSubMatrix(B2, PETSC_DECIDE, PETSC_DECIDE, PetscMin(1, K - 1), PetscMin(2, K), &T3));
339: PetscCall(MatAXPY(T, 1.0, T2, SAME_NONZERO_PATTERN));
340: PetscCall(MatAXPY(T3, -1.0, T, SAME_NONZERO_PATTERN));
341: PetscCall(MatNorm(T3, NORM_FROBENIUS, &err));
342: if (err > PETSC_SMALL) {
343: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with MatDenseGetSubMatrix\n"));
344: PetscCall(MatView(T3, NULL));
345: }
346: PetscCall(MatDenseRestoreSubMatrix(B, &T));
347: PetscCall(MatDenseRestoreSubMatrix(T4, &T2));
348: PetscCall(MatDenseRestoreSubMatrix(B2, &T3));
349: PetscCall(CheckLocal(B, NULL, aB, NULL));
350: PetscCall(MatDestroy(&B2));
351: PetscCall(MatDestroy(&T4));
352: if (N >= 2) {
353: PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
354: PetscCall(MatDuplicate(B, MAT_DO_NOT_COPY_VALUES, &T4));
355: PetscCall(MatSetRandom(T4, NULL));
356: PetscCall(MatAXPY(B2, 1.0, T4, SAME_NONZERO_PATTERN));
357: PetscCall(MatDenseGetSubMatrix(B, N - 2, PETSC_DECIDE, PetscMin(1, K - 1), PetscMin(2, K), &T));
358: PetscCall(MatDenseGetSubMatrix(T4, N - 2, PETSC_DECIDE, PetscMin(1, K - 1), PetscMin(2, K), &T2));
359: PetscCall(MatDenseGetSubMatrix(B2, N - 2, PETSC_DECIDE, PetscMin(1, K - 1), PetscMin(2, K), &T3));
360: PetscCall(MatAXPY(T, 1.0, T2, SAME_NONZERO_PATTERN));
361: PetscCall(MatAXPY(T3, -1.0, T, SAME_NONZERO_PATTERN));
362: PetscCall(MatNorm(T3, NORM_FROBENIUS, &err));
363: if (err > PETSC_SMALL) {
364: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with MatDenseGetSubMatrix\n"));
365: PetscCall(MatView(T3, NULL));
366: }
367: PetscCall(MatDenseRestoreSubMatrix(B, &T));
368: PetscCall(MatDenseRestoreSubMatrix(T4, &T2));
369: PetscCall(MatDenseRestoreSubMatrix(B2, &T3));
370: PetscCall(CheckLocal(B, NULL, aB, NULL));
371: PetscCall(MatDestroy(&B2));
372: PetscCall(MatDestroy(&T4));
373: PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &B2));
374: PetscCall(MatDuplicate(B, MAT_DO_NOT_COPY_VALUES, &T4));
375: PetscCall(MatSetRandom(T4, NULL));
376: PetscCall(MatAXPY(B2, 1.0, T4, SAME_NONZERO_PATTERN));
377: PetscCall(MatDenseGetSubMatrix(B, PETSC_DECIDE, 2, PetscMin(1, K - 1), PetscMin(2, K), &T));
378: PetscCall(MatDenseGetSubMatrix(T4, PETSC_DECIDE, 2, PetscMin(1, K - 1), PetscMin(2, K), &T2));
379: PetscCall(MatDenseGetSubMatrix(B2, PETSC_DECIDE, 2, PetscMin(1, K - 1), PetscMin(2, K), &T3));
380: PetscCall(MatAXPY(T, 1.0, T2, SAME_NONZERO_PATTERN));
381: PetscCall(MatAXPY(T3, -1.0, T, SAME_NONZERO_PATTERN));
382: PetscCall(MatNorm(T3, NORM_FROBENIUS, &err));
383: if (err > PETSC_SMALL) {
384: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with MatDenseGetSubMatrix\n"));
385: PetscCall(MatView(T3, NULL));
386: }
387: PetscCall(MatDenseRestoreSubMatrix(B, &T));
388: PetscCall(MatDenseRestoreSubMatrix(T4, &T2));
389: PetscCall(MatDenseRestoreSubMatrix(B2, &T3));
390: PetscCall(CheckLocal(B, NULL, aB, NULL));
391: PetscCall(MatDestroy(&B2));
392: PetscCall(MatDestroy(&T4));
393: }
394: }
396: /* Test reusing a previously allocated dense buffer */
397: PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
398: PetscCall(CheckLocal(B, X, aB, aX));
399: PetscCall(MatMatMultEqual(A, B, X, 10, &flg));
400: if (!flg) {
401: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with reusage\n"));
402: PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
403: PetscCall(MatAXPY(T, -1.0, X, SAME_NONZERO_PATTERN));
404: PetscCall(MatView(T, NULL));
405: PetscCall(MatDestroy(&T));
406: }
408: /* Test MatTransposeMat and MatMatTranspose */
409: if (testmattmat) {
410: PetscCall(MatTransposeMatMult(A, X, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B));
411: PetscCall(CheckLocal(B, X, aB, aX));
412: PetscCall(MatTransposeMatMultEqual(A, X, B, 10, &flg));
413: if (!flg) {
414: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with reusage (MatTransposeMat)\n"));
415: PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &B));
416: PetscCall(MatAXPY(T, -1.0, B, SAME_NONZERO_PATTERN));
417: PetscCall(MatView(T, NULL));
418: PetscCall(MatDestroy(&T));
419: }
420: }
421: if (testmatmatt) {
422: PetscCall(MatMatTransposeMult(A, Bt, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
423: PetscCall(CheckLocal(Bt, X, aBt, aX));
424: PetscCall(MatMatTransposeMultEqual(A, Bt, X, 10, &flg));
425: if (!flg) {
426: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with reusage (MatMatTranspose)\n"));
427: PetscCall(MatMatTransposeMult(A, Bt, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
428: PetscCall(MatAXPY(T, -1.0, X, SAME_NONZERO_PATTERN));
429: PetscCall(MatView(T, NULL));
430: PetscCall(MatDestroy(&T));
431: }
432: }
434: /* Test projection operations (PtAP and RARt) */
435: if (testproj) {
436: PetscCall(MatPtAP(A, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &PtAP));
437: PetscCall(MatPtAPMultEqual(A, B, PtAP, 10, &flg));
438: if (!flg) {
439: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with PtAP\n"));
440: PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
441: PetscCall(MatTransposeMatMult(B, T, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T2));
442: PetscCall(MatAXPY(T2, -1.0, PtAP, SAME_NONZERO_PATTERN));
443: PetscCall(MatView(T2, NULL));
444: PetscCall(MatDestroy(&T2));
445: PetscCall(MatDestroy(&T));
446: }
447: PetscCall(PetscMalloc1((k + ldr) * M, &dataR));
448: PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, m, K, M, dataR, &R));
449: PetscCall(MatDenseSetLDA(R, k + ldr));
450: PetscCall(MatSetRandom(R, NULL));
451: if (testrart) { /* fails for AIJCUSPARSE because RA operation is not defined */
452: PetscCall(MatRARt(A, R, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &RARt));
453: PetscCall(MatRARtMultEqual(A, R, RARt, 10, &flg));
454: if (!flg) {
455: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with RARt\n"));
456: PetscCall(MatMatTransposeMult(A, R, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
457: PetscCall(MatMatMult(R, T, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T2));
458: PetscCall(MatAXPY(T2, -1.0, RARt, SAME_NONZERO_PATTERN));
459: PetscCall(MatView(T2, NULL));
460: PetscCall(MatDestroy(&T2));
461: PetscCall(MatDestroy(&T));
462: }
463: }
464: }
466: /* Test MatDenseGetColumnVec and friends */
467: PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
468: PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
469: PetscCall(MatDuplicate(T, MAT_DO_NOT_COPY_VALUES, &T2));
470: for (k = 0; k < K; k++) {
471: Vec Xv, Tv, T2v;
473: PetscCall(MatDenseGetColumnVecRead(X, k, &Xv));
474: PetscCall(MatDenseGetColumnVec(T, k, &Tv));
475: PetscCall(MatDenseGetColumnVecWrite(T2, k, &T2v));
476: PetscCall(VecCopy(Xv, T2v));
477: PetscCall(VecAXPY(Tv, -1., Xv));
478: PetscCall(MatDenseRestoreColumnVecRead(X, k, &Xv));
479: PetscCall(MatDenseRestoreColumnVec(T, k, &Tv));
480: PetscCall(MatDenseRestoreColumnVecWrite(T2, k, &T2v));
481: }
482: PetscCall(MatNorm(T, NORM_FROBENIUS, &err));
483: if (err > PETSC_SMALL) {
484: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with MatDenseGetColumnVec\n"));
485: PetscCall(MatView(T, NULL));
486: }
487: PetscCall(MatAXPY(T2, -1., X, SAME_NONZERO_PATTERN));
488: PetscCall(MatNorm(T2, NORM_FROBENIUS, &err));
489: if (err > PETSC_SMALL) {
490: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with MatDenseGetColumnVecWrite\n"));
491: PetscCall(MatView(T2, NULL));
492: }
493: PetscCall(MatDestroy(&T));
494: PetscCall(MatDestroy(&T2));
496: /* Test with MatShell */
497: PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &T));
498: PetscCall(MatConvert(T, MATSHELL, MAT_INITIAL_MATRIX, &T2));
499: PetscCall(MatDestroy(&T));
501: /* scale matrix */
502: PetscCall(MatScale(A, 2.0));
503: PetscCall(MatScale(T2, 2.0));
504: PetscCall(MatCreateVecs(A, &r, &l));
505: PetscCall(VecSetRandom(r, NULL));
506: PetscCall(VecSetRandom(l, NULL));
507: PetscCall(MatCreateVecs(T2, &rs, &ls));
508: PetscCall(VecCopy(r, rs));
509: PetscCall(VecCopy(l, ls));
510: if (testproj) {
511: PetscCall(MatDiagonalScale(A, r, r));
512: PetscCall(MatDiagonalScale(T2, rs, rs));
513: } else {
514: PetscCall(MatDiagonalScale(A, l, r));
515: PetscCall(MatDiagonalScale(T2, ls, rs));
516: }
517: PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &T));
518: PetscCall(MatAXPY(A, 4.5, T, SAME_NONZERO_PATTERN));
519: PetscCall(MatAXPY(T2, 4.5, T, DIFFERENT_NONZERO_PATTERN));
520: PetscCall(MatMultEqual(T2, A, 10, &flg));
521: if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with MATSHELL (MatMult)\n"));
522: PetscCall(MatMultTransposeEqual(T2, A, 10, &flg));
523: if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with MATSHELL (MatMultTranspose)\n"));
524: PetscCall(MatDestroy(&T));
525: PetscCall(VecDestroy(&ls));
526: PetscCall(VecDestroy(&rs));
527: PetscCall(VecDestroy(&l));
528: PetscCall(VecDestroy(&r));
530: /* recompute projections, test reusage */
531: if (PtAP) PetscCall(MatPtAP(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &PtAP));
532: if (RARt) PetscCall(MatRARt(A, R, MAT_REUSE_MATRIX, PETSC_DETERMINE, &RARt));
533: if (testshellops) { /* test callbacks for user defined MatProducts */
534: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_AB, NULL, MyMatShellMatMultNumeric, NULL, MATDENSE, MATDENSE));
535: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_AB, NULL, MyMatShellMatMultNumeric, NULL, MATDENSECUDA, MATDENSECUDA));
536: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_AtB, NULL, MyMatTransposeShellMatMultNumeric, NULL, MATDENSE, MATDENSE));
537: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_AtB, NULL, MyMatTransposeShellMatMultNumeric, NULL, MATDENSECUDA, MATDENSECUDA));
538: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_ABt, NULL, MyMatShellMatTransposeMultNumeric, NULL, MATDENSE, MATDENSE));
539: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_ABt, NULL, MyMatShellMatTransposeMultNumeric, NULL, MATDENSECUDA, MATDENSECUDA));
540: if (testproj) {
541: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_PtAP, MyPtShellPMultSymbolic, MyPtShellPMultNumeric, proj_destroy, MATDENSE, MATSHELL));
542: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_PtAP, MyPtShellPMultSymbolic, MyPtShellPMultNumeric, proj_destroy, MATDENSECUDA, MATSHELL));
543: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_RARt, MyRShellRtMultSymbolic, MyRShellRtMultNumeric, proj_destroy, MATDENSE, MATSHELL));
544: PetscCall(MatShellSetMatProductOperation(T2, MATPRODUCT_RARt, MyRShellRtMultSymbolic, MyRShellRtMultNumeric, proj_destroy, MATDENSECUDA, MATSHELL));
545: }
546: }
547: PetscCall(CheckLocal(B, X, aB, aX));
548: /* we either use the shell operations or the loop over columns code, applying the operator */
549: PetscCall(MatMatMult(T2, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
550: PetscCall(CheckLocal(B, X, aB, aX));
551: PetscCall(MatMatMultEqual(T2, B, X, 10, &flg));
552: if (!flg) {
553: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with reusage (MATSHELL)\n"));
554: PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
555: PetscCall(MatAXPY(T, -1.0, X, SAME_NONZERO_PATTERN));
556: PetscCall(MatView(T, NULL));
557: PetscCall(MatDestroy(&T));
558: }
559: if (testproj) {
560: PetscCall(MatPtAPMultEqual(T2, B, PtAP, 10, &flg));
561: if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with PtAP (MATSHELL)\n"));
562: if (testshellops) { /* projections fail if the product operations are not specified */
563: PetscCall(MatPtAP(T2, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
564: PetscCall(MatPtAP(T2, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &T));
565: PetscCall(MatPtAPMultEqual(T2, B, T, 10, &flg));
566: if (!flg) {
567: Mat TE;
569: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with PtAP (MATSHELL user defined)\n"));
570: PetscCall(MatComputeOperator(T, MATDENSE, &TE));
571: PetscCall(MatView(TE, NULL));
572: PetscCall(MatView(PtAP, NULL));
573: PetscCall(MatAXPY(TE, -1.0, PtAP, SAME_NONZERO_PATTERN));
574: PetscCall(MatView(TE, NULL));
575: PetscCall(MatDestroy(&TE));
576: }
577: PetscCall(MatDestroy(&T));
578: }
579: if (RARt) {
580: PetscCall(MatRARtMultEqual(T2, R, RARt, 10, &flg));
581: if (!flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with RARt (MATSHELL)\n"));
582: }
583: if (testshellops) {
584: PetscCall(MatRARt(T2, R, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
585: PetscCall(MatRARt(T2, R, MAT_REUSE_MATRIX, PETSC_DETERMINE, &T));
586: PetscCall(MatRARtMultEqual(T2, R, T, 10, &flg));
587: if (!flg) {
588: Mat TE;
590: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with RARt (MATSHELL user defined)\n"));
591: PetscCall(MatComputeOperator(T, MATDENSE, &TE));
592: PetscCall(MatView(TE, NULL));
593: if (RARt) {
594: PetscCall(MatView(RARt, NULL));
595: PetscCall(MatAXPY(TE, -1.0, RARt, SAME_NONZERO_PATTERN));
596: PetscCall(MatView(TE, NULL));
597: }
598: PetscCall(MatDestroy(&TE));
599: }
600: PetscCall(MatDestroy(&T));
601: }
602: }
604: if (testmattmat) { /* we either use the shell operations or the loop over columns code applying the transposed operator */
605: PetscCall(MatTransposeMatMult(T2, X, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B));
606: PetscCall(CheckLocal(B, X, aB, aX));
607: PetscCall(MatTransposeMatMultEqual(T2, X, B, 10, &flg));
608: if (!flg) {
609: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with reusage (MatTranspose, MATSHELL)\n"));
610: PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
611: PetscCall(MatAXPY(T, -1.0, B, SAME_NONZERO_PATTERN));
612: PetscCall(MatView(T, NULL));
613: PetscCall(MatDestroy(&T));
614: }
615: }
616: if (testmatmatt && testshellops) { /* only when shell operations are set */
617: PetscCall(MatMatTransposeMult(T2, Bt, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
618: PetscCall(CheckLocal(Bt, X, aBt, aX));
619: PetscCall(MatMatTransposeMultEqual(T2, Bt, X, 10, &flg));
620: if (!flg) {
621: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with reusage (MatMatTranspose, MATSHELL)\n"));
622: PetscCall(MatMatTransposeMult(A, Bt, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
623: PetscCall(MatAXPY(T, -1.0, X, SAME_NONZERO_PATTERN));
624: PetscCall(MatView(T, NULL));
625: PetscCall(MatDestroy(&T));
626: }
627: }
628: PetscCall(MatDestroy(&T2));
630: if (testnest) { /* test with MatNest */
631: Mat NA;
633: PetscCall(MatCreateNest(PETSC_COMM_WORLD, 1, NULL, 1, NULL, &A, &NA));
634: PetscCall(MatViewFromOptions(NA, NULL, "-NA_view"));
635: PetscCall(MatMatMult(NA, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
636: PetscCall(CheckLocal(B, X, aB, aX));
637: PetscCall(MatMatMultEqual(NA, B, X, 10, &flg));
638: if (!flg) {
639: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with Nest\n"));
640: PetscCall(MatMatMult(NA, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
641: PetscCall(MatAXPY(T, -1.0, X, SAME_NONZERO_PATTERN));
642: PetscCall(MatView(T, NULL));
643: PetscCall(MatDestroy(&T));
644: }
645: PetscCall(MatDestroy(&NA));
646: }
648: if (testtranspose) { /* test with Transpose */
649: Mat TA;
651: PetscCall(MatCreateTranspose(A, &TA));
652: PetscCall(MatMatMult(TA, X, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B));
653: PetscCall(CheckLocal(B, X, aB, aX));
654: PetscCall(MatMatMultEqual(TA, X, B, 10, &flg));
655: if (!flg) {
656: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with Transpose\n"));
657: PetscCall(MatMatMult(TA, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
658: PetscCall(MatAXPY(T, -1.0, B, SAME_NONZERO_PATTERN));
659: PetscCall(MatView(T, NULL));
660: PetscCall(MatDestroy(&T));
661: }
662: PetscCall(MatDestroy(&TA));
663: }
665: if (testhtranspose) { /* test with Hermitian Transpose */
666: Mat TA;
668: PetscCall(MatCreateHermitianTranspose(A, &TA));
669: PetscCall(MatMatMult(TA, X, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B));
670: PetscCall(CheckLocal(B, X, aB, aX));
671: PetscCall(MatMatMultEqual(TA, X, B, 10, &flg));
672: if (!flg) {
673: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with Transpose\n"));
674: PetscCall(MatMatMult(TA, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
675: PetscCall(MatAXPY(T, -1.0, B, SAME_NONZERO_PATTERN));
676: PetscCall(MatView(T, NULL));
677: PetscCall(MatDestroy(&T));
678: }
679: PetscCall(MatDestroy(&TA));
680: }
682: if (testtt) { /* test with Transpose(Transpose) */
683: Mat TA, TTA;
685: PetscCall(MatCreateTranspose(A, &TA));
686: PetscCall(MatCreateTranspose(TA, &TTA));
687: PetscCall(MatMatMult(TTA, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
688: PetscCall(CheckLocal(B, X, aB, aX));
689: PetscCall(MatMatMultEqual(TTA, B, X, 10, &flg));
690: if (!flg) {
691: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with Transpose(Transpose)\n"));
692: PetscCall(MatMatMult(TTA, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
693: PetscCall(MatAXPY(T, -1.0, X, SAME_NONZERO_PATTERN));
694: PetscCall(MatView(T, NULL));
695: PetscCall(MatDestroy(&T));
696: }
697: PetscCall(MatDestroy(&TA));
698: PetscCall(MatDestroy(&TTA));
699: }
701: if (testcircular) { /* test circular */
702: Mat AB;
704: PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &AB));
705: PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &X));
706: PetscCall(CheckLocal(B, X, aB, aX));
707: if (M == N && N == K) PetscCall(MatMatMult(A, X, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B));
708: else PetscCall(MatTransposeMatMult(A, X, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B));
709: PetscCall(CheckLocal(B, X, aB, aX));
710: PetscCall(MatDestroy(&AB));
711: }
713: /* Test by Pierre Jolivet */
714: {
715: Mat C, D, D2, AtA;
716: PetscCall(MatCreateNormal(A, &AtA));
717: PetscCall(MatDuplicate(X, MAT_DO_NOT_COPY_VALUES, &C));
718: PetscCall(MatDuplicate(B, MAT_DO_NOT_COPY_VALUES, &D));
719: PetscCall(MatDuplicate(B, MAT_DO_NOT_COPY_VALUES, &D2));
720: PetscCall(MatSetRandom(B, NULL));
721: PetscCall(MatSetRandom(C, NULL));
722: PetscCall(MatSetRandom(D, NULL));
723: PetscCall(MatSetRandom(D2, NULL));
724: PetscCall(MatProductCreateWithMat(A, B, NULL, C));
725: PetscCall(MatProductSetType(C, MATPRODUCT_AB));
726: PetscCall(MatProductSetFromOptions(C));
727: PetscCall(MatProductSymbolic(C));
728: PetscCall(MatProductCreateWithMat(A, C, NULL, D));
729: PetscCall(MatProductSetType(D, MATPRODUCT_AtB));
730: PetscCall(MatProductSetFromOptions(D));
731: PetscCall(MatProductSymbolic(D));
732: PetscCall(MatProductNumeric(C));
733: PetscCall(MatMatMultEqual(A, B, C, 10, &flg));
734: if (!flg) {
735: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with Normal (AB != C)\n"));
736: PetscCall(MatView(A, NULL));
737: PetscCall(MatView(B, NULL));
738: PetscCall(MatView(C, NULL));
739: }
740: PetscCall(MatProductNumeric(D));
741: PetscCall(MatMatMultEqual(AtA, B, D, 10, &flg));
742: if (!flg) {
743: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error with Normal (2)\n"));
744: PetscCall(MatMatMult(AtA, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &T));
745: PetscCall(MatView(D, NULL));
746: PetscCall(MatView(T, NULL));
747: PetscCall(MatAXPY(T, -1.0, D, SAME_NONZERO_PATTERN));
748: PetscCall(MatView(T, NULL));
749: PetscCall(MatDestroy(&T));
750: }
751: PetscCall(MatDestroy(&C));
752: PetscCall(MatDestroy(&D));
753: PetscCall(MatDestroy(&D2));
754: PetscCall(MatDestroy(&AtA));
755: }
757: PetscCall(MatDestroy(&X));
758: PetscCall(MatDestroy(&Bt));
759: PetscCall(MatDestroy(&B));
760: PetscCall(MatDestroy(&A));
761: PetscCall(MatDestroy(&R));
762: PetscCall(MatDestroy(&PtAP));
763: PetscCall(MatDestroy(&RARt));
764: PetscCall(PetscFree(dataX));
765: PetscCall(PetscFree(dataB));
766: PetscCall(PetscFree(dataR));
767: PetscCall(PetscFree(dataBt));
768: PetscCall(PetscFinalize());
769: return 0;
770: }
772: /*TEST
774: test:
775: output_file: output/empty.out
776: suffix: 1
777: args: -local {{0 1}} -testshellops
779: test:
780: output_file: output/empty.out
781: requires: cuda
782: suffix: 1_cuda
783: args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{seqaijcusparse seqaij}} -testshellops {{0 1}}
785: test:
786: output_file: output/empty.out
787: nsize: 2
788: suffix: 1_par
789: args: -local {{0 1}} -testmatmatt 0
791: test:
792: output_file: output/empty.out
793: requires: cuda
794: nsize: 2
795: suffix: 1_par_cuda
796: args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{mpiaijcusparse mpiaij}} -testnest 0 -testmatmatt 0 -matproduct_batch_size 3
798: test:
799: output_file: output/empty.out
800: suffix: 2
801: nsize: 1
802: args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}}
804: testset:
805: requires: cuda
806: output_file: output/empty.out
807: nsize: 1
808: args: -M 7 -N 9 -K 2 -local {{0 1}} -testnest 0 -A_mat_type {{seqdensecuda seqdense}} -xgpu {{0 1}} -bgpu {{0 1}}
809: test:
810: requires: !complex
811: suffix: 2_cuda_real
812: test:
813: # complex+single gives a little bigger error in the MatDenseGetColumnVec test
814: requires: complex !single
815: suffix: 2_cuda_complex
817: test:
818: output_file: output/empty.out
819: suffix: 2_par
820: nsize: 2
821: args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} -testcircular -testmatmatt 0
823: test:
824: requires: cuda
825: output_file: output/empty.out
826: suffix: 2_par_cuda
827: nsize: 2
828: args: -M 11 -N 9 -K 1 -local {{0 1}} -testcircular 0 -A_mat_type mpiaijcusparse -xgpu -bgpu -testnest 0 -testmatmatt 0
830: test:
831: output_file: output/empty.out
832: suffix: 3
833: nsize: {{1 3}}
834: args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type sbaij -symm -testproj 0 -testmatmatt 0
836: test:
837: output_file: output/empty.out
838: suffix: 4
839: nsize: 1
840: args: -M 3 -N 3 -K 3 -local {{0 1}} -testcircular
842: test:
843: output_file: output/empty.out
844: suffix: 5
845: nsize: {{2 4}}
846: args: -M 3 -N 3 -K 3 -local {{0 1}} -testcircular -testmatmatt 0
848: test:
849: output_file: output/empty.out
850: suffix: 6
851: nsize: 1
852: args: -M {{1 3}} -N {{2 5}} -K {{1 2}} -local {{0 1}} -testcircular
854: test:
855: output_file: output/empty.out
856: suffix: 7
857: nsize: 1
858: args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type dense -testnest -testcircular
860: TEST*/