xref: /petsc/src/mat/tests/ex96.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests sequential and parallel DMCreateMatrix(), MatMatMult() and MatPtAP()\n\
2c4762a1bSJed Brown   -Mx <xg>, where <xg> = number of coarse grid points in the x-direction\n\
3c4762a1bSJed Brown   -My <yg>, where <yg> = number of coarse grid points in the y-direction\n\
4c4762a1bSJed Brown   -Mz <zg>, where <zg> = number of coarse grid points in the z-direction\n\
5c4762a1bSJed Brown   -Npx <npx>, where <npx> = number of processors in the x-direction\n\
6c4762a1bSJed Brown   -Npy <npy>, where <npy> = number of processors in the y-direction\n\
7c4762a1bSJed Brown   -Npz <npz>, where <npz> = number of processors in the z-direction\n\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /*
10c4762a1bSJed Brown     This test is modified from ~src/ksp/tests/ex19.c.
11c4762a1bSJed Brown     Example of usage: mpiexec -n 3 ./ex96 -Mx 10 -My 10 -Mz 10
12c4762a1bSJed Brown */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscdm.h>
15c4762a1bSJed Brown #include <petscdmda.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /* User-defined application contexts */
18c4762a1bSJed Brown typedef struct {
19c4762a1bSJed Brown   PetscInt mx, my, mz;     /* number grid points in x, y and z direction */
20c4762a1bSJed Brown   Vec      localX, localF; /* local vectors with ghost region */
21c4762a1bSJed Brown   DM       da;
22c4762a1bSJed Brown   Vec      x, b, r; /* global vectors */
23c4762a1bSJed Brown   Mat      J;       /* Jacobian on grid */
24c4762a1bSJed Brown } GridCtx;
25c4762a1bSJed Brown typedef struct {
26c4762a1bSJed Brown   GridCtx  fine;
27c4762a1bSJed Brown   GridCtx  coarse;
28c4762a1bSJed Brown   PetscInt ratio;
29c4762a1bSJed Brown   Mat      Ii; /* interpolation from coarse to fine */
30c4762a1bSJed Brown } AppCtx;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown #define COARSE_LEVEL 0
33c4762a1bSJed Brown #define FINE_LEVEL   1
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown       Mm_ratio - ration of grid lines between fine and coarse grids.
37c4762a1bSJed Brown */
main(int argc,char ** argv)38d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
39d71ae5a4SJacob Faibussowitsch {
40c4762a1bSJed Brown   AppCtx          user;
41c4762a1bSJed Brown   PetscInt        Npx = PETSC_DECIDE, Npy = PETSC_DECIDE, Npz = PETSC_DECIDE;
42c4762a1bSJed Brown   PetscMPIInt     size, rank;
43c4762a1bSJed Brown   PetscInt        m, n, M, N, i, nrows;
44c4762a1bSJed Brown   PetscScalar     one  = 1.0;
45c4762a1bSJed Brown   PetscReal       fill = 2.0;
46c4762a1bSJed Brown   Mat             A, A_tmp, P, C, C1, C2;
47c4762a1bSJed Brown   PetscScalar    *array, none = -1.0, alpha;
48c4762a1bSJed Brown   Vec             x, v1, v2, v3, v4;
49c4762a1bSJed Brown   PetscReal       norm, norm_tmp, norm_tmp1, tol = 100. * PETSC_MACHINE_EPSILON;
50c4762a1bSJed Brown   PetscRandom     rdm;
51c4762a1bSJed Brown   PetscBool       Test_MatMatMult = PETSC_TRUE, Test_MatPtAP = PETSC_TRUE, Test_3D = PETSC_TRUE, flg;
52c4762a1bSJed Brown   const PetscInt *ia, *ja;
53c4762a1bSJed Brown 
54327415f7SBarry Smith   PetscFunctionBeginUser;
559566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   user.ratio     = 2;
599371c9d4SSatish Balay   user.coarse.mx = 20;
609371c9d4SSatish Balay   user.coarse.my = 20;
619371c9d4SSatish Balay   user.coarse.mz = 20;
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-Mx", &user.coarse.mx, NULL));
649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-My", &user.coarse.my, NULL));
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-Mz", &user.coarse.mz, NULL));
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ratio", &user.ratio, NULL));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   if (user.coarse.mz) Test_3D = PETSC_TRUE;
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-Npx", &Npx, NULL));
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-Npy", &Npy, NULL));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-Npz", &Npz, NULL));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Set up distributed array for fine grid */
77c4762a1bSJed Brown   if (!Test_3D) {
789566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.coarse.mx, user.coarse.my, Npx, Npy, 1, 1, NULL, NULL, &user.coarse.da));
79c4762a1bSJed Brown   } else {
809566063dSJacob Faibussowitsch     PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.coarse.mx, user.coarse.my, user.coarse.mz, Npx, Npy, Npz, 1, 1, NULL, NULL, NULL, &user.coarse.da));
81c4762a1bSJed Brown   }
829566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.coarse.da));
839566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.coarse.da));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* This makes sure the coarse DMDA has the same partition as the fine DMDA */
869566063dSJacob Faibussowitsch   PetscCall(DMRefine(user.coarse.da, PetscObjectComm((PetscObject)user.coarse.da), &user.fine.da));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Test DMCreateMatrix()                                         */
89c4762a1bSJed Brown   /*------------------------------------------------------------*/
909566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(user.fine.da, MATAIJ));
919566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.fine.da, &A));
929566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(user.fine.da, MATBAIJ));
939566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.fine.da, &C));
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(MatConvert(C, MATAIJ, MAT_INITIAL_MATRIX, &A_tmp)); /* not work for mpisbaij matrix! */
969566063dSJacob Faibussowitsch   PetscCall(MatEqual(A, A_tmp, &flg));
9728b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "A != C");
989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A_tmp));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*------------------------------------------------------------*/
102c4762a1bSJed Brown 
1039566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
1049566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
105dd400576SPatrick Sanan   /* if (rank == 0) printf("A %d, %d\n",M,N); */
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* set val=one to A */
108c4762a1bSJed Brown   if (size == 1) {
1099566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
110c4762a1bSJed Brown     if (flg) {
1119566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(A, &array));
112c4762a1bSJed Brown       for (i = 0; i < ia[nrows]; i++) array[i] = one;
1139566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(A, &array));
114c4762a1bSJed Brown     }
1159566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
116c4762a1bSJed Brown   } else {
117c4762a1bSJed Brown     Mat AA, AB;
1189566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, NULL));
1199566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
120c4762a1bSJed Brown     if (flg) {
1219566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AA, &array));
122c4762a1bSJed Brown       for (i = 0; i < ia[nrows]; i++) array[i] = one;
1239566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AA, &array));
124c4762a1bSJed Brown     }
1259566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
1269566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
127c4762a1bSJed Brown     if (flg) {
1289566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AB, &array));
129c4762a1bSJed Brown       for (i = 0; i < ia[nrows]; i++) array[i] = one;
1309566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AB, &array));
131c4762a1bSJed Brown     }
1329566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg));
133c4762a1bSJed Brown   }
1349566063dSJacob Faibussowitsch   /* PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); */
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* Create interpolation between the fine and coarse grids */
1379566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(user.coarse.da, user.fine.da, &P, NULL));
1389566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(P, &m, &n));
1399566063dSJacob Faibussowitsch   PetscCall(MatGetSize(P, &M, &N));
140dd400576SPatrick Sanan   /* if (rank == 0) printf("P %d, %d\n",M,N); */
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* Create vectors v1 and v2 that are compatible with A */
1439566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &v1));
1449566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
1459566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(v1, m, PETSC_DECIDE));
1469566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(v1));
1479566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v1, &v2));
1489566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
1499566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Test MatMatMult(): C = A*P */
152c4762a1bSJed Brown   /*----------------------------*/
153c4762a1bSJed Brown   if (Test_MatMatMult) {
1549566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A_tmp));
1559566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A_tmp, P, MAT_INITIAL_MATRIX, fill, &C));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX - reuse symbolic C */
158c4762a1bSJed Brown     alpha = 1.0;
159c4762a1bSJed Brown     for (i = 0; i < 2; i++) {
160c4762a1bSJed Brown       alpha -= 0.1;
1619566063dSJacob Faibussowitsch       PetscCall(MatScale(A_tmp, alpha));
1629566063dSJacob Faibussowitsch       PetscCall(MatMatMult(A_tmp, P, MAT_REUSE_MATRIX, fill, &C));
163c4762a1bSJed Brown     }
164c4762a1bSJed Brown     /* Free intermediate data structures created for reuse of C=Pt*A*P */
1659566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown     /* Test MatDuplicate()        */
168c4762a1bSJed Brown     /*----------------------------*/
1699566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &C1));
1709566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C1, MAT_COPY_VALUES, &C2));
1719566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
1729566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C2));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown     /* Create vector x that is compatible with P */
1759566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
1769566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(P, NULL, &n));
1779566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
1789566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(x));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown     norm = 0.0;
181c4762a1bSJed Brown     for (i = 0; i < 10; i++) {
1829566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x, rdm));
1839566063dSJacob Faibussowitsch       PetscCall(MatMult(P, x, v1));
1849566063dSJacob Faibussowitsch       PetscCall(MatMult(A_tmp, v1, v2)); /* v2 = A*P*x */
1859566063dSJacob Faibussowitsch       PetscCall(MatMult(C, x, v1));      /* v1 = C*x   */
1869566063dSJacob Faibussowitsch       PetscCall(VecAXPY(v1, none, v2));
1879566063dSJacob Faibussowitsch       PetscCall(VecNorm(v1, NORM_1, &norm_tmp));
1889566063dSJacob Faibussowitsch       PetscCall(VecNorm(v2, NORM_1, &norm_tmp1));
189c4762a1bSJed Brown       norm_tmp /= norm_tmp1;
190c4762a1bSJed Brown       if (norm_tmp > norm) norm = norm_tmp;
191c4762a1bSJed Brown     }
19248a46eb9SPierre Jolivet     if (norm >= tol && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMatMult(), |v1 - v2|/|v2|: %g\n", (double)norm));
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
1959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1969566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_tmp));
197c4762a1bSJed Brown   }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* Test P^T * A * P - MatPtAP() */
200c4762a1bSJed Brown   /*------------------------------*/
201c4762a1bSJed Brown   if (Test_MatPtAP) {
2029566063dSJacob Faibussowitsch     PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, fill, &C));
2039566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(C, &m, &n));
204c4762a1bSJed Brown 
205c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX - reuse symbolic C */
206c4762a1bSJed Brown     alpha = 1.0;
207c4762a1bSJed Brown     for (i = 0; i < 1; i++) {
208c4762a1bSJed Brown       alpha -= 0.1;
2099566063dSJacob Faibussowitsch       PetscCall(MatScale(A, alpha));
2109566063dSJacob Faibussowitsch       PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, fill, &C));
211c4762a1bSJed Brown     }
212c4762a1bSJed Brown 
213c4762a1bSJed Brown     /* Free intermediate data structures created for reuse of C=Pt*A*P */
2149566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown     /* Test MatDuplicate()        */
217c4762a1bSJed Brown     /*----------------------------*/
2189566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &C1));
2199566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C1, MAT_COPY_VALUES, &C2));
2209566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
2219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C2));
222c4762a1bSJed Brown 
223c4762a1bSJed Brown     /* Create vector x that is compatible with P */
2249566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
2259566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(P, &m, &n));
2269566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
2279566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(x));
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD, &v3));
2309566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(v3, n, PETSC_DECIDE));
2319566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(v3));
2329566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(v3, &v4));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown     norm = 0.0;
235c4762a1bSJed Brown     for (i = 0; i < 10; i++) {
2369566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x, rdm));
2379566063dSJacob Faibussowitsch       PetscCall(MatMult(P, x, v1));
2389566063dSJacob Faibussowitsch       PetscCall(MatMult(A, v1, v2)); /* v2 = A*P*x */
239c4762a1bSJed Brown 
2409566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(P, v2, v3)); /* v3 = Pt*A*P*x */
2419566063dSJacob Faibussowitsch       PetscCall(MatMult(C, x, v4));           /* v3 = C*x   */
2429566063dSJacob Faibussowitsch       PetscCall(VecAXPY(v4, none, v3));
2439566063dSJacob Faibussowitsch       PetscCall(VecNorm(v4, NORM_1, &norm_tmp));
2449566063dSJacob Faibussowitsch       PetscCall(VecNorm(v3, NORM_1, &norm_tmp1));
245c4762a1bSJed Brown 
246c4762a1bSJed Brown       norm_tmp /= norm_tmp1;
247c4762a1bSJed Brown       if (norm_tmp > norm) norm = norm_tmp;
248c4762a1bSJed Brown     }
24948a46eb9SPierre Jolivet     if (norm >= tol && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatPtAP(), |v3 - v4|/|v3|: %g\n", (double)norm));
2509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
2519566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v3));
2529566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v4));
2539566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
254c4762a1bSJed Brown   }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Clean up */
2579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2589566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
2599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v1));
2609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v2));
2619566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.fine.da));
2629566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.coarse.da));
2639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
2649566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
265b122ec5aSJacob Faibussowitsch   return 0;
266c4762a1bSJed Brown }
267c4762a1bSJed Brown 
268c4762a1bSJed Brown /*TEST
269c4762a1bSJed Brown 
270c4762a1bSJed Brown    test:
271c4762a1bSJed Brown       args: -Mx 10 -My 5 -Mz 10
272*3886731fSPierre Jolivet       output_file: output/empty.out
273c4762a1bSJed Brown 
274c4762a1bSJed Brown    test:
275c4762a1bSJed Brown       suffix: nonscalable
276c4762a1bSJed Brown       nsize: 3
277c4762a1bSJed Brown       args: -Mx 10 -My 5 -Mz 10
278*3886731fSPierre Jolivet       output_file: output/empty.out
279c4762a1bSJed Brown 
280c4762a1bSJed Brown    test:
281c4762a1bSJed Brown       suffix: scalable
282c4762a1bSJed Brown       nsize: 3
283c4762a1bSJed Brown       args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable
284*3886731fSPierre Jolivet       output_file: output/empty.out
285c4762a1bSJed Brown 
286c4762a1bSJed Brown    test:
287c4762a1bSJed Brown      suffix: seq_scalable
288c4762a1bSJed Brown      nsize: 3
2893e662e0bSHong Zhang      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm scalable -inner_offdiag_mat_product_algorithm scalable
290*3886731fSPierre Jolivet      output_file: output/empty.out
291c4762a1bSJed Brown 
292c4762a1bSJed Brown    test:
293c4762a1bSJed Brown      suffix: seq_sorted
294c4762a1bSJed Brown      nsize: 3
2953e662e0bSHong Zhang      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm sorted -inner_offdiag_mat_product_algorithm sorted
296*3886731fSPierre Jolivet      output_file: output/empty.out
297c4762a1bSJed Brown 
298c4762a1bSJed Brown    test:
299c4762a1bSJed Brown      suffix: seq_scalable_fast
300c4762a1bSJed Brown      nsize: 3
3013e662e0bSHong Zhang      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm scalable_fast -inner_offdiag_mat_product_algorithm scalable_fast
302*3886731fSPierre Jolivet      output_file: output/empty.out
303c4762a1bSJed Brown 
304c4762a1bSJed Brown    test:
305c4762a1bSJed Brown      suffix: seq_heap
306c4762a1bSJed Brown      nsize: 3
3073e662e0bSHong Zhang      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm heap -inner_offdiag_mat_product_algorithm heap
308*3886731fSPierre Jolivet      output_file: output/empty.out
309c4762a1bSJed Brown 
310c4762a1bSJed Brown    test:
311c4762a1bSJed Brown      suffix: seq_btheap
312c4762a1bSJed Brown      nsize: 3
3133e662e0bSHong Zhang      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm btheap -inner_offdiag_mat_product_algorithm btheap
314*3886731fSPierre Jolivet      output_file: output/empty.out
315c4762a1bSJed Brown 
316c4762a1bSJed Brown    test:
317c4762a1bSJed Brown      suffix: seq_llcondensed
318c4762a1bSJed Brown      nsize: 3
3193e662e0bSHong Zhang      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm llcondensed -inner_offdiag_mat_product_algorithm llcondensed
320*3886731fSPierre Jolivet      output_file: output/empty.out
321c4762a1bSJed Brown 
322c4762a1bSJed Brown    test:
323c4762a1bSJed Brown      suffix: seq_rowmerge
324c4762a1bSJed Brown      nsize: 3
3253e662e0bSHong Zhang      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm rowmerge -inner_offdiag_mat_product_algorithm rowmerge
326*3886731fSPierre Jolivet      output_file: output/empty.out
327c4762a1bSJed Brown 
328c4762a1bSJed Brown    test:
329c4762a1bSJed Brown      suffix: allatonce
330c4762a1bSJed Brown      nsize: 3
331c4762a1bSJed Brown      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce
332*3886731fSPierre Jolivet      output_file: output/empty.out
333c4762a1bSJed Brown 
334c4762a1bSJed Brown    test:
335c4762a1bSJed Brown      suffix: allatonce_merged
336c4762a1bSJed Brown      nsize: 3
337c4762a1bSJed Brown      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce_merged
338*3886731fSPierre Jolivet      output_file: output/empty.out
339c4762a1bSJed Brown 
340c4762a1bSJed Brown TEST*/
341