xref: /petsc/src/mat/tests/ex96.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
42c4762a1bSJed Brown   AppCtx         user;
43c4762a1bSJed Brown   PetscInt       Npx=PETSC_DECIDE,Npy=PETSC_DECIDE,Npz=PETSC_DECIDE;
44c4762a1bSJed Brown   PetscMPIInt    size,rank;
45c4762a1bSJed Brown   PetscInt       m,n,M,N,i,nrows;
46c4762a1bSJed Brown   PetscScalar    one = 1.0;
47c4762a1bSJed Brown   PetscReal      fill=2.0;
48c4762a1bSJed Brown   Mat            A,A_tmp,P,C,C1,C2;
49c4762a1bSJed Brown   PetscScalar    *array,none = -1.0,alpha;
50c4762a1bSJed Brown   Vec            x,v1,v2,v3,v4;
51c4762a1bSJed Brown   PetscReal      norm,norm_tmp,norm_tmp1,tol=100.*PETSC_MACHINE_EPSILON;
52c4762a1bSJed Brown   PetscRandom    rdm;
53c4762a1bSJed Brown   PetscBool      Test_MatMatMult=PETSC_TRUE,Test_MatPtAP=PETSC_TRUE,Test_3D=PETSC_TRUE,flg;
54c4762a1bSJed Brown   const PetscInt *ia,*ja;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
57c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL);CHKERRQ(ierr);
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   user.ratio     = 2;
60c4762a1bSJed Brown   user.coarse.mx = 20; user.coarse.my = 20; user.coarse.mz = 20;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL);CHKERRQ(ierr);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   if (user.coarse.mz) Test_3D = PETSC_TRUE;
68c4762a1bSJed Brown 
69ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
70ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
71c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-Npx",&Npx,NULL);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-Npy",&Npy,NULL);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-Npz",&Npz,NULL);CHKERRQ(ierr);
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* Set up distributed array for fine grid */
76c4762a1bSJed Brown   if (!Test_3D) {
77c4762a1bSJed Brown     ierr = 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);CHKERRQ(ierr);
78c4762a1bSJed Brown   } else {
79c4762a1bSJed Brown     ierr = 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);CHKERRQ(ierr);
80c4762a1bSJed Brown   }
81c4762a1bSJed Brown   ierr = DMSetFromOptions(user.coarse.da);CHKERRQ(ierr);
82c4762a1bSJed Brown   ierr = DMSetUp(user.coarse.da);CHKERRQ(ierr);
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* This makes sure the coarse DMDA has the same partition as the fine DMDA */
85c4762a1bSJed Brown   ierr = DMRefine(user.coarse.da,PetscObjectComm((PetscObject)user.coarse.da),&user.fine.da);CHKERRQ(ierr);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Test DMCreateMatrix()                                         */
88c4762a1bSJed Brown   /*------------------------------------------------------------*/
89c4762a1bSJed Brown   ierr = DMSetMatType(user.fine.da,MATAIJ);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = DMCreateMatrix(user.fine.da,&A);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = DMSetMatType(user.fine.da,MATBAIJ);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = DMCreateMatrix(user.fine.da,&C);CHKERRQ(ierr);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   ierr = MatConvert(C,MATAIJ,MAT_INITIAL_MATRIX,&A_tmp);CHKERRQ(ierr); /* not work for mpisbaij matrix! */
95c4762a1bSJed Brown   ierr = MatEqual(A,A_tmp,&flg);CHKERRQ(ierr);
96*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != C");
97c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = MatDestroy(&A_tmp);CHKERRQ(ierr);
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*------------------------------------------------------------*/
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
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) {
108c4762a1bSJed Brown     ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
109c4762a1bSJed Brown     if (flg) {
110c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr);
111c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
112c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(A,&array);CHKERRQ(ierr);
113c4762a1bSJed Brown     }
114c4762a1bSJed Brown     ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
115c4762a1bSJed Brown   } else {
116c4762a1bSJed Brown     Mat AA,AB;
117c4762a1bSJed Brown     ierr = MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL);CHKERRQ(ierr);
118c4762a1bSJed Brown     ierr = MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
119c4762a1bSJed Brown     if (flg) {
120c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(AA,&array);CHKERRQ(ierr);
121c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
122c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(AA,&array);CHKERRQ(ierr);
123c4762a1bSJed Brown     }
124c4762a1bSJed Brown     ierr = MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
125c4762a1bSJed Brown     ierr = MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
126c4762a1bSJed Brown     if (flg) {
127c4762a1bSJed Brown       ierr = MatSeqAIJGetArray(AB,&array);CHKERRQ(ierr);
128c4762a1bSJed Brown       for (i=0; i<ia[nrows]; i++) array[i] = one;
129c4762a1bSJed Brown       ierr = MatSeqAIJRestoreArray(AB,&array);CHKERRQ(ierr);
130c4762a1bSJed Brown     }
131c4762a1bSJed Brown     ierr = MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr);
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown   /* ierr = MatView(A, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* Create interpolation between the fine and coarse grids */
136c4762a1bSJed Brown   ierr = DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL);CHKERRQ(ierr);
137c4762a1bSJed Brown   ierr = MatGetLocalSize(P,&m,&n);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = MatGetSize(P,&M,&N);CHKERRQ(ierr);
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 */
142c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&v1);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = VecSetSizes(v1,m,PETSC_DECIDE);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = VecSetFromOptions(v1);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = VecDuplicate(v1,&v2);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Test MatMatMult(): C = A*P */
151c4762a1bSJed Brown   /*----------------------------*/
152c4762a1bSJed Brown   if (Test_MatMatMult) {
153c4762a1bSJed Brown     ierr = MatDuplicate(A,MAT_COPY_VALUES,&A_tmp);CHKERRQ(ierr);
154c4762a1bSJed Brown     ierr = MatMatMult(A_tmp,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr);
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;
160c4762a1bSJed Brown       ierr   = MatScale(A_tmp,alpha);CHKERRQ(ierr);
161c4762a1bSJed Brown       ierr   = MatMatMult(A_tmp,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr);
162c4762a1bSJed Brown     }
163c4762a1bSJed Brown     /* Free intermediate data structures created for reuse of C=Pt*A*P */
1646718818eSStefano Zampini     ierr = MatProductClear(C);CHKERRQ(ierr);
165c4762a1bSJed Brown 
166c4762a1bSJed Brown     /* Test MatDuplicate()        */
167c4762a1bSJed Brown     /*----------------------------*/
168c4762a1bSJed Brown     ierr = MatDuplicate(C,MAT_COPY_VALUES,&C1);CHKERRQ(ierr);
169c4762a1bSJed Brown     ierr = MatDuplicate(C1,MAT_COPY_VALUES,&C2);CHKERRQ(ierr);
170c4762a1bSJed Brown     ierr = MatDestroy(&C1);CHKERRQ(ierr);
171c4762a1bSJed Brown     ierr = MatDestroy(&C2);CHKERRQ(ierr);
172c4762a1bSJed Brown 
173c4762a1bSJed Brown     /* Create vector x that is compatible with P */
174c4762a1bSJed Brown     ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
175c4762a1bSJed Brown     ierr = MatGetLocalSize(P,NULL,&n);CHKERRQ(ierr);
176c4762a1bSJed Brown     ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = VecSetFromOptions(x);CHKERRQ(ierr);
178c4762a1bSJed Brown 
179c4762a1bSJed Brown     norm = 0.0;
180c4762a1bSJed Brown     for (i=0; i<10; i++) {
181c4762a1bSJed Brown       ierr      = VecSetRandom(x,rdm);CHKERRQ(ierr);
182c4762a1bSJed Brown       ierr      = MatMult(P,x,v1);CHKERRQ(ierr);
183c4762a1bSJed Brown       ierr      = MatMult(A_tmp,v1,v2);CHKERRQ(ierr); /* v2 = A*P*x */
184c4762a1bSJed Brown       ierr      = MatMult(C,x,v1);CHKERRQ(ierr);  /* v1 = C*x   */
185c4762a1bSJed Brown       ierr      = VecAXPY(v1,none,v2);CHKERRQ(ierr);
186c4762a1bSJed Brown       ierr      = VecNorm(v1,NORM_1,&norm_tmp);CHKERRQ(ierr);
187c4762a1bSJed Brown       ierr      = VecNorm(v2,NORM_1,&norm_tmp1);CHKERRQ(ierr);
188c4762a1bSJed Brown       norm_tmp /= norm_tmp1;
189c4762a1bSJed Brown       if (norm_tmp > norm) norm = norm_tmp;
190c4762a1bSJed Brown     }
191dd400576SPatrick Sanan     if (norm >= tol && rank == 0) {
192c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult(), |v1 - v2|/|v2|: %g\n",(double)norm);CHKERRQ(ierr);
193c4762a1bSJed Brown     }
194c4762a1bSJed Brown 
195c4762a1bSJed Brown     ierr = VecDestroy(&x);CHKERRQ(ierr);
196c4762a1bSJed Brown     ierr = MatDestroy(&C);CHKERRQ(ierr);
197c4762a1bSJed Brown     ierr = MatDestroy(&A_tmp);CHKERRQ(ierr);
198c4762a1bSJed Brown   }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* Test P^T * A * P - MatPtAP() */
201c4762a1bSJed Brown   /*------------------------------*/
202c4762a1bSJed Brown   if (Test_MatPtAP) {
203c4762a1bSJed Brown     ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr);
204c4762a1bSJed Brown     ierr = MatGetLocalSize(C,&m,&n);CHKERRQ(ierr);
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;
210c4762a1bSJed Brown       ierr   = MatScale(A,alpha);CHKERRQ(ierr);
211c4762a1bSJed Brown       ierr   = MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr);
212c4762a1bSJed Brown     }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown     /* Free intermediate data structures created for reuse of C=Pt*A*P */
2156718818eSStefano Zampini     ierr = MatProductClear(C);CHKERRQ(ierr);
216c4762a1bSJed Brown 
217c4762a1bSJed Brown     /* Test MatDuplicate()        */
218c4762a1bSJed Brown     /*----------------------------*/
219c4762a1bSJed Brown     ierr = MatDuplicate(C,MAT_COPY_VALUES,&C1);CHKERRQ(ierr);
220c4762a1bSJed Brown     ierr = MatDuplicate(C1,MAT_COPY_VALUES,&C2);CHKERRQ(ierr);
221c4762a1bSJed Brown     ierr = MatDestroy(&C1);CHKERRQ(ierr);
222c4762a1bSJed Brown     ierr = MatDestroy(&C2);CHKERRQ(ierr);
223c4762a1bSJed Brown 
224c4762a1bSJed Brown     /* Create vector x that is compatible with P */
225c4762a1bSJed Brown     ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
226c4762a1bSJed Brown     ierr = MatGetLocalSize(P,&m,&n);CHKERRQ(ierr);
227c4762a1bSJed Brown     ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr);
228c4762a1bSJed Brown     ierr = VecSetFromOptions(x);CHKERRQ(ierr);
229c4762a1bSJed Brown 
230c4762a1bSJed Brown     ierr = VecCreate(PETSC_COMM_WORLD,&v3);CHKERRQ(ierr);
231c4762a1bSJed Brown     ierr = VecSetSizes(v3,n,PETSC_DECIDE);CHKERRQ(ierr);
232c4762a1bSJed Brown     ierr = VecSetFromOptions(v3);CHKERRQ(ierr);
233c4762a1bSJed Brown     ierr = VecDuplicate(v3,&v4);CHKERRQ(ierr);
234c4762a1bSJed Brown 
235c4762a1bSJed Brown     norm = 0.0;
236c4762a1bSJed Brown     for (i=0; i<10; i++) {
237c4762a1bSJed Brown       ierr = VecSetRandom(x,rdm);CHKERRQ(ierr);
238c4762a1bSJed Brown       ierr = MatMult(P,x,v1);CHKERRQ(ierr);
239c4762a1bSJed Brown       ierr = MatMult(A,v1,v2);CHKERRQ(ierr);  /* v2 = A*P*x */
240c4762a1bSJed Brown 
241c4762a1bSJed Brown       ierr = MatMultTranspose(P,v2,v3);CHKERRQ(ierr); /* v3 = Pt*A*P*x */
242c4762a1bSJed Brown       ierr = MatMult(C,x,v4);CHKERRQ(ierr);           /* v3 = C*x   */
243c4762a1bSJed Brown       ierr = VecAXPY(v4,none,v3);CHKERRQ(ierr);
244c4762a1bSJed Brown       ierr = VecNorm(v4,NORM_1,&norm_tmp);CHKERRQ(ierr);
245c4762a1bSJed Brown       ierr = VecNorm(v3,NORM_1,&norm_tmp1);CHKERRQ(ierr);
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) {
251c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatPtAP(), |v3 - v4|/|v3|: %g\n",(double)norm);CHKERRQ(ierr);
252c4762a1bSJed Brown     }
253c4762a1bSJed Brown     ierr = MatDestroy(&C);CHKERRQ(ierr);
254c4762a1bSJed Brown     ierr = VecDestroy(&v3);CHKERRQ(ierr);
255c4762a1bSJed Brown     ierr = VecDestroy(&v4);CHKERRQ(ierr);
256c4762a1bSJed Brown     ierr = VecDestroy(&x);CHKERRQ(ierr);
257c4762a1bSJed Brown   }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Clean up */
260c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
261c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = VecDestroy(&v1);CHKERRQ(ierr);
263c4762a1bSJed Brown   ierr = VecDestroy(&v2);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = DMDestroy(&user.fine.da);CHKERRQ(ierr);
265c4762a1bSJed Brown   ierr = DMDestroy(&user.coarse.da);CHKERRQ(ierr);
266c4762a1bSJed Brown   ierr = MatDestroy(&P);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = PetscFinalize();
268c4762a1bSJed Brown   return ierr;
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