xref: /petsc/src/mat/tests/ex96.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Tests sequential and parallel DMCreateMatrix(), MatMatMult() and MatPtAP()\n\
3c4762a1bSJed Brown   -Mx <xg>, where <xg> = number of coarse grid points in the x-direction\n\
4c4762a1bSJed Brown   -My <yg>, where <yg> = number of coarse grid points in the y-direction\n\
5c4762a1bSJed Brown   -Mz <zg>, where <zg> = number of coarse grid points in the z-direction\n\
6c4762a1bSJed Brown   -Npx <npx>, where <npx> = number of processors in the x-direction\n\
7c4762a1bSJed Brown   -Npy <npy>, where <npy> = number of processors in the y-direction\n\
8c4762a1bSJed Brown   -Npz <npz>, where <npz> = number of processors in the z-direction\n\n";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown     This test is modified from ~src/ksp/tests/ex19.c.
12c4762a1bSJed Brown     Example of usage: mpiexec -n 3 ./ex96 -Mx 10 -My 10 -Mz 10
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown #include <petscdm.h>
16c4762a1bSJed Brown #include <petscdmda.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /* User-defined application contexts */
19c4762a1bSJed Brown typedef struct {
20c4762a1bSJed Brown   PetscInt mx,my,mz;            /* number grid points in x, y and z direction */
21c4762a1bSJed Brown   Vec      localX,localF;       /* local vectors with ghost region */
22c4762a1bSJed Brown   DM       da;
23c4762a1bSJed Brown   Vec      x,b,r;               /* global vectors */
24c4762a1bSJed Brown   Mat      J;                   /* Jacobian on grid */
25c4762a1bSJed Brown } GridCtx;
26c4762a1bSJed Brown typedef struct {
27c4762a1bSJed Brown   GridCtx  fine;
28c4762a1bSJed Brown   GridCtx  coarse;
29c4762a1bSJed Brown   PetscInt ratio;
30c4762a1bSJed Brown   Mat      Ii;                  /* interpolation from coarse to fine */
31c4762a1bSJed Brown } AppCtx;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown #define COARSE_LEVEL 0
34c4762a1bSJed Brown #define FINE_LEVEL   1
35c4762a1bSJed Brown 
36c4762a1bSJed Brown /*
37c4762a1bSJed Brown       Mm_ratio - ration of grid lines between fine and coarse grids.
38c4762a1bSJed Brown */
39c4762a1bSJed Brown int main(int argc,char **argv)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   AppCtx         user;
42c4762a1bSJed Brown   PetscInt       Npx=PETSC_DECIDE,Npy=PETSC_DECIDE,Npz=PETSC_DECIDE;
43c4762a1bSJed Brown   PetscMPIInt    size,rank;
44c4762a1bSJed Brown   PetscInt       m,n,M,N,i,nrows;
45c4762a1bSJed Brown   PetscScalar    one = 1.0;
46c4762a1bSJed Brown   PetscReal      fill=2.0;
47c4762a1bSJed Brown   Mat            A,A_tmp,P,C,C1,C2;
48c4762a1bSJed Brown   PetscScalar    *array,none = -1.0,alpha;
49c4762a1bSJed Brown   Vec            x,v1,v2,v3,v4;
50c4762a1bSJed Brown   PetscReal      norm,norm_tmp,norm_tmp1,tol=100.*PETSC_MACHINE_EPSILON;
51c4762a1bSJed Brown   PetscRandom    rdm;
52c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_TRUE,Test_MatPtAP=PETSC_TRUE,Test_3D=PETSC_TRUE,flg;
53c4762a1bSJed Brown   const PetscInt *ia,*ja;
54c4762a1bSJed Brown 
55*327415f7SBarry Smith   PetscFunctionBeginUser;
569566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   user.ratio     = 2;
60c4762a1bSJed Brown   user.coarse.mx = 20; user.coarse.my = 20; user.coarse.mz = 20;
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL));
649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL));
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   if (user.coarse.mz) Test_3D = PETSC_TRUE;
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Npx",&Npx,NULL));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Npy",&Npy,NULL));
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-Npz",&Npz,NULL));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* Set up distributed array for fine grid */
76c4762a1bSJed Brown   if (!Test_3D) {
779566063dSJacob 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));
78c4762a1bSJed Brown   } else {
799566063dSJacob 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));
80c4762a1bSJed Brown   }
819566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.coarse.da));
829566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.coarse.da));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* This makes sure the coarse DMDA has the same partition as the fine DMDA */
859566063dSJacob Faibussowitsch   PetscCall(DMRefine(user.coarse.da,PetscObjectComm((PetscObject)user.coarse.da),&user.fine.da));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Test DMCreateMatrix()                                         */
88c4762a1bSJed Brown   /*------------------------------------------------------------*/
899566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(user.fine.da,MATAIJ));
909566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.fine.da,&A));
919566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(user.fine.da,MATBAIJ));
929566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.fine.da,&C));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(MatConvert(C,MATAIJ,MAT_INITIAL_MATRIX,&A_tmp)); /* not work for mpisbaij matrix! */
959566063dSJacob Faibussowitsch   PetscCall(MatEqual(A,A_tmp,&flg));
9628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != C");
979566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A_tmp));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*------------------------------------------------------------*/
101c4762a1bSJed Brown 
1029566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
1039566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&M,&N));
104dd400576SPatrick Sanan   /* if (rank == 0) printf("A %d, %d\n",M,N); */
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* set val=one to A */
107c4762a1bSJed Brown   if (size == 1) {
1089566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
109c4762a1bSJed Brown     if (flg) {
1109566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(A,&array));
111c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
1129566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(A,&array));
113c4762a1bSJed Brown     }
1149566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
115c4762a1bSJed Brown   } else {
116c4762a1bSJed Brown     Mat AA,AB;
1179566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL));
1189566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
119c4762a1bSJed Brown     if (flg) {
1209566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AA,&array));
121c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
1229566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AA,&array));
123c4762a1bSJed Brown     }
1249566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
1259566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
126c4762a1bSJed Brown     if (flg) {
1279566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetArray(AB,&array));
128c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
1299566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJRestoreArray(AB,&array));
130c4762a1bSJed Brown     }
1319566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg));
132c4762a1bSJed Brown   }
1339566063dSJacob Faibussowitsch   /* PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); */
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* Create interpolation between the fine and coarse grids */
1369566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL));
1379566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(P,&m,&n));
1389566063dSJacob Faibussowitsch   PetscCall(MatGetSize(P,&M,&N));
139dd400576SPatrick Sanan   /* if (rank == 0) printf("P %d, %d\n",M,N); */
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Create vectors v1 and v2 that are compatible with A */
1429566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&v1));
1439566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,NULL));
1449566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(v1,m,PETSC_DECIDE));
1459566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(v1));
1469566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v1,&v2));
1479566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
1489566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Test MatMatMult(): C = A*P */
151c4762a1bSJed Brown   /*----------------------------*/
152c4762a1bSJed Brown   if (Test_MatMatMult) {
1539566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A_tmp));
1549566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A_tmp,P,MAT_INITIAL_MATRIX,fill,&C));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX - reuse symbolic C */
157c4762a1bSJed Brown     alpha=1.0;
158c4762a1bSJed Brown     for (i=0; i<2; i++) {
159c4762a1bSJed Brown       alpha -= 0.1;
1609566063dSJacob Faibussowitsch       PetscCall(MatScale(A_tmp,alpha));
1619566063dSJacob Faibussowitsch       PetscCall(MatMatMult(A_tmp,P,MAT_REUSE_MATRIX,fill,&C));
162c4762a1bSJed Brown     }
163c4762a1bSJed Brown     /* Free intermediate data structures created for reuse of C=Pt*A*P */
1649566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown     /* Test MatDuplicate()        */
167c4762a1bSJed Brown     /*----------------------------*/
1689566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&C1));
1699566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C1,MAT_COPY_VALUES,&C2));
1709566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
1719566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C2));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown     /* Create vector x that is compatible with P */
1749566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
1759566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(P,NULL,&n));
1769566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
1779566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(x));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown     norm = 0.0;
180c4762a1bSJed Brown     for (i=0; i<10; i++) {
1819566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x,rdm));
1829566063dSJacob Faibussowitsch       PetscCall(MatMult(P,x,v1));
1839566063dSJacob Faibussowitsch       PetscCall(MatMult(A_tmp,v1,v2)); /* v2 = A*P*x */
1849566063dSJacob Faibussowitsch       PetscCall(MatMult(C,x,v1));  /* v1 = C*x   */
1859566063dSJacob Faibussowitsch       PetscCall(VecAXPY(v1,none,v2));
1869566063dSJacob Faibussowitsch       PetscCall(VecNorm(v1,NORM_1,&norm_tmp));
1879566063dSJacob Faibussowitsch       PetscCall(VecNorm(v2,NORM_1,&norm_tmp1));
188c4762a1bSJed Brown       norm_tmp /= norm_tmp1;
189c4762a1bSJed Brown       if (norm_tmp > norm) norm = norm_tmp;
190c4762a1bSJed Brown     }
191dd400576SPatrick Sanan     if (norm >= tol && rank == 0) {
1929566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult(), |v1 - v2|/|v2|: %g\n",(double)norm));
193c4762a1bSJed Brown     }
194c4762a1bSJed Brown 
1959566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
1969566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1979566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_tmp));
198c4762a1bSJed Brown   }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* Test P^T * A * P - MatPtAP() */
201c4762a1bSJed Brown   /*------------------------------*/
202c4762a1bSJed Brown   if (Test_MatPtAP) {
2039566063dSJacob Faibussowitsch     PetscCall(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C));
2049566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(C,&m,&n));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX - reuse symbolic C */
207c4762a1bSJed Brown     alpha=1.0;
208c4762a1bSJed Brown     for (i=0; i<1; i++) {
209c4762a1bSJed Brown       alpha -= 0.1;
2109566063dSJacob Faibussowitsch       PetscCall(MatScale(A,alpha));
2119566063dSJacob Faibussowitsch       PetscCall(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C));
212c4762a1bSJed Brown     }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown     /* Free intermediate data structures created for reuse of C=Pt*A*P */
2159566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown     /* Test MatDuplicate()        */
218c4762a1bSJed Brown     /*----------------------------*/
2199566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&C1));
2209566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(C1,MAT_COPY_VALUES,&C2));
2219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C1));
2229566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C2));
223c4762a1bSJed Brown 
224c4762a1bSJed Brown     /* Create vector x that is compatible with P */
2259566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
2269566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(P,&m,&n));
2279566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
2289566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(x));
229c4762a1bSJed Brown 
2309566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD,&v3));
2319566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(v3,n,PETSC_DECIDE));
2329566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(v3));
2339566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(v3,&v4));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown     norm = 0.0;
236c4762a1bSJed Brown     for (i=0; i<10; i++) {
2379566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x,rdm));
2389566063dSJacob Faibussowitsch       PetscCall(MatMult(P,x,v1));
2399566063dSJacob Faibussowitsch       PetscCall(MatMult(A,v1,v2));  /* v2 = A*P*x */
240c4762a1bSJed Brown 
2419566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(P,v2,v3)); /* v3 = Pt*A*P*x */
2429566063dSJacob Faibussowitsch       PetscCall(MatMult(C,x,v4));           /* v3 = C*x   */
2439566063dSJacob Faibussowitsch       PetscCall(VecAXPY(v4,none,v3));
2449566063dSJacob Faibussowitsch       PetscCall(VecNorm(v4,NORM_1,&norm_tmp));
2459566063dSJacob Faibussowitsch       PetscCall(VecNorm(v3,NORM_1,&norm_tmp1));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown       norm_tmp /= norm_tmp1;
248c4762a1bSJed Brown       if (norm_tmp > norm) norm = norm_tmp;
249c4762a1bSJed Brown     }
250dd400576SPatrick Sanan     if (norm >= tol && rank == 0) {
2519566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatPtAP(), |v3 - v4|/|v3|: %g\n",(double)norm));
252c4762a1bSJed Brown     }
2539566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
2549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v3));
2559566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v4));
2569566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
257c4762a1bSJed Brown   }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Clean up */
2609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2619566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
2629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v1));
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v2));
2649566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.fine.da));
2659566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.coarse.da));
2669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
2679566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
268b122ec5aSJacob Faibussowitsch   return 0;
269c4762a1bSJed Brown }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown /*TEST
272c4762a1bSJed Brown 
273c4762a1bSJed Brown    test:
274c4762a1bSJed Brown       args: -Mx 10 -My 5 -Mz 10
275c4762a1bSJed Brown       output_file: output/ex96_1.out
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    test:
278c4762a1bSJed Brown       suffix: nonscalable
279c4762a1bSJed Brown       nsize: 3
280c4762a1bSJed Brown       args: -Mx 10 -My 5 -Mz 10
281c4762a1bSJed Brown       output_file: output/ex96_1.out
282c4762a1bSJed Brown 
283c4762a1bSJed Brown    test:
284c4762a1bSJed Brown       suffix: scalable
285c4762a1bSJed Brown       nsize: 3
286c4762a1bSJed Brown       args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable
287c4762a1bSJed Brown       output_file: output/ex96_1.out
288c4762a1bSJed Brown 
289c4762a1bSJed Brown    test:
290c4762a1bSJed Brown      suffix: seq_scalable
291c4762a1bSJed Brown      nsize: 3
2923e662e0bSHong 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
293c4762a1bSJed Brown      output_file: output/ex96_1.out
294c4762a1bSJed Brown 
295c4762a1bSJed Brown    test:
296c4762a1bSJed Brown      suffix: seq_sorted
297c4762a1bSJed Brown      nsize: 3
2983e662e0bSHong 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
299c4762a1bSJed Brown      output_file: output/ex96_1.out
300c4762a1bSJed Brown 
301c4762a1bSJed Brown    test:
302c4762a1bSJed Brown      suffix: seq_scalable_fast
303c4762a1bSJed Brown      nsize: 3
3043e662e0bSHong 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
305c4762a1bSJed Brown      output_file: output/ex96_1.out
306c4762a1bSJed Brown 
307c4762a1bSJed Brown    test:
308c4762a1bSJed Brown      suffix: seq_heap
309c4762a1bSJed Brown      nsize: 3
3103e662e0bSHong 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
311c4762a1bSJed Brown      output_file: output/ex96_1.out
312c4762a1bSJed Brown 
313c4762a1bSJed Brown    test:
314c4762a1bSJed Brown      suffix: seq_btheap
315c4762a1bSJed Brown      nsize: 3
3163e662e0bSHong 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
317c4762a1bSJed Brown      output_file: output/ex96_1.out
318c4762a1bSJed Brown 
319c4762a1bSJed Brown    test:
320c4762a1bSJed Brown      suffix: seq_llcondensed
321c4762a1bSJed Brown      nsize: 3
3223e662e0bSHong 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
323c4762a1bSJed Brown      output_file: output/ex96_1.out
324c4762a1bSJed Brown 
325c4762a1bSJed Brown    test:
326c4762a1bSJed Brown      suffix: seq_rowmerge
327c4762a1bSJed Brown      nsize: 3
3283e662e0bSHong 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
329c4762a1bSJed Brown      output_file: output/ex96_1.out
330c4762a1bSJed Brown 
331c4762a1bSJed Brown    test:
332c4762a1bSJed Brown      suffix: allatonce
333c4762a1bSJed Brown      nsize: 3
334c4762a1bSJed Brown      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce
335c4762a1bSJed Brown      output_file: output/ex96_1.out
336c4762a1bSJed Brown 
337c4762a1bSJed Brown    test:
338c4762a1bSJed Brown      suffix: allatonce_merged
339c4762a1bSJed Brown      nsize: 3
340c4762a1bSJed Brown      args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce_merged
341c4762a1bSJed Brown      output_file: output/ex96_1.out
342c4762a1bSJed Brown 
343c4762a1bSJed Brown TEST*/
344