xref: /petsc/src/mat/tests/ex76.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests cholesky, icc factorization and solve on sequential aij, baij and sbaij matrices. \n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Vec            x,y,b;
9c4762a1bSJed Brown   Mat            A;           /* linear system matrix */
10c4762a1bSJed Brown   Mat            sA,sC;       /* symmetric part of the matrices */
11c4762a1bSJed Brown   PetscInt       n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],block, row,Ii,J,n1,lvl;
12c4762a1bSJed Brown   PetscMPIInt    size;
13c4762a1bSJed Brown   PetscReal      norm2;
14c4762a1bSJed Brown   PetscScalar    neg_one = -1.0,four=4.0,value[3];
15c4762a1bSJed Brown   IS             perm,cperm;
16c4762a1bSJed Brown   PetscRandom    rdm;
17c4762a1bSJed Brown   PetscBool      reorder = PETSC_FALSE,displ = PETSC_FALSE;
18c4762a1bSJed Brown   MatFactorInfo  factinfo;
19c4762a1bSJed Brown   PetscBool      equal;
20c4762a1bSJed Brown   PetscBool      TestAIJ = PETSC_FALSE,TestBAIJ = PETSC_TRUE;
21c4762a1bSJed Brown   PetscInt       TestShift=0;
22c4762a1bSJed Brown 
23*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
24*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
252c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
26*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
27*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
28*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-reorder",&reorder,NULL));
29*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testaij",&TestAIJ,NULL));
30*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-testShift",&TestShift,NULL));
31*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   n = mbs*bs;
34c4762a1bSJed Brown   if (TestAIJ) { /* A is in aij format */
35*9566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,nz,NULL,&A));
36c4762a1bSJed Brown     TestBAIJ = PETSC_FALSE;
37c4762a1bSJed Brown   } else { /* A is in baij format */
38*9566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL,&A));
39c4762a1bSJed Brown     TestAIJ = PETSC_FALSE;
40c4762a1bSJed Brown   }
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* Assemble matrix */
43c4762a1bSJed Brown   if (bs == 1) {
44*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL));
45c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
46c4762a1bSJed Brown       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
47c4762a1bSJed Brown       for (i=1; i<n-1; i++) {
48c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
49*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
50c4762a1bSJed Brown       }
51c4762a1bSJed Brown       i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown       value[0]= 0.1; value[1]=-1; value[2]=2;
54*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown       i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown       value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
59*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
60c4762a1bSJed Brown     } else if (prob ==2) { /* matrix for the five point stencil */
61c4762a1bSJed Brown       n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001);
622c71b3e2SJacob Faibussowitsch       PetscCheckFalse(n1*n1 - n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!");
63c4762a1bSJed Brown       for (i=0; i<n1; i++) {
64c4762a1bSJed Brown         for (j=0; j<n1; j++) {
65c4762a1bSJed Brown           Ii = j + n1*i;
66c4762a1bSJed Brown           if (i>0) {
67c4762a1bSJed Brown             J    = Ii - n1;
68*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
69c4762a1bSJed Brown           }
70c4762a1bSJed Brown           if (i<n1-1) {
71c4762a1bSJed Brown             J    = Ii + n1;
72*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
73c4762a1bSJed Brown           }
74c4762a1bSJed Brown           if (j>0) {
75c4762a1bSJed Brown             J    = Ii - 1;
76*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
77c4762a1bSJed Brown           }
78c4762a1bSJed Brown           if (j<n1-1) {
79c4762a1bSJed Brown             J    = Ii + 1;
80*9566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
81c4762a1bSJed Brown           }
82*9566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES));
83c4762a1bSJed Brown         }
84c4762a1bSJed Brown       }
85c4762a1bSJed Brown     }
86c4762a1bSJed Brown   } else { /* bs > 1 */
87c4762a1bSJed Brown     for (block=0; block<n/bs; block++) {
88c4762a1bSJed Brown       /* diagonal blocks */
89c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
90c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
91c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
92*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
93c4762a1bSJed Brown       }
94c4762a1bSJed Brown       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
95c4762a1bSJed Brown 
96c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
97*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
100c4762a1bSJed Brown 
101c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
102*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
103c4762a1bSJed Brown     }
104c4762a1bSJed Brown     /* off-diagonal blocks */
105c4762a1bSJed Brown     value[0]=-1.0;
106c4762a1bSJed Brown     for (i=0; i<(n/bs-1)*bs; i++) {
107c4762a1bSJed Brown       col[0]=i+bs;
108*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES));
109c4762a1bSJed Brown       col[0]=i; row=i+bs;
110*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES));
111c4762a1bSJed Brown     }
112c4762a1bSJed Brown   }
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   if (TestShift) {
115c4762a1bSJed Brown     /* set diagonals in the 0-th block as 0 for testing shift numerical factor */
116c4762a1bSJed Brown     for (i=0; i<bs; i++) {
117c4762a1bSJed Brown       row  = i; col[0] = i; value[0] = 0.0;
118*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES));
119c4762a1bSJed Brown     }
120c4762a1bSJed Brown   }
121c4762a1bSJed Brown 
122*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
123*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Test MatConvert */
126*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
127*9566063dSJacob Faibussowitsch   PetscCall(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA));
128*9566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A,sA,20,&equal));
12928b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_USER,"A != sA");
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
132*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&Ii,&J));
133*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(sA,&i,&j));
1342c71b3e2SJacob Faibussowitsch   PetscCheckFalse(i-Ii || j-J,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetOwnershipRange() in MatSBAIJ format");
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* Vectors */
137*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
138*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
139*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x));
140*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&b));
141*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&y));
142*9566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x,rdm));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Test MatReordering() - not work on sbaij matrix */
145c4762a1bSJed Brown   if (reorder) {
146*9566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A,MATORDERINGRCM,&perm,&cperm));
147c4762a1bSJed Brown   } else {
148*9566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&cperm));
149c4762a1bSJed Brown   }
150*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* initialize factinfo */
153*9566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&factinfo));
154c4762a1bSJed Brown   if (TestShift == 1) {
155c4762a1bSJed Brown     factinfo.shifttype   = (PetscReal)MAT_SHIFT_NONZERO;
156c4762a1bSJed Brown     factinfo.shiftamount = 0.1;
157c4762a1bSJed Brown   } else if (TestShift == 2) {
158c4762a1bSJed Brown     factinfo.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* Test MatCholeskyFactor(), MatICCFactor() */
162c4762a1bSJed Brown   /*------------------------------------------*/
163c4762a1bSJed Brown   /* Test aij matrix A */
164c4762a1bSJed Brown   if (TestAIJ) {
165c4762a1bSJed Brown     if (displ) {
166*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"AIJ: \n"));
167c4762a1bSJed Brown     }
168c4762a1bSJed Brown     i = 0;
169c4762a1bSJed Brown     for (lvl=-1; lvl<10; lvl++) {
170c4762a1bSJed Brown       if (lvl==-1) {  /* Cholesky factor */
171c4762a1bSJed Brown         factinfo.fill = 5.0;
172c4762a1bSJed Brown 
173*9566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC));
174*9566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(sC,A,perm,&factinfo));
175c4762a1bSJed Brown       } else {       /* incomplete Cholesky factor */
176c4762a1bSJed Brown         factinfo.fill   = 5.0;
177c4762a1bSJed Brown         factinfo.levels = lvl;
178c4762a1bSJed Brown 
179*9566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC));
180*9566063dSJacob Faibussowitsch         PetscCall(MatICCFactorSymbolic(sC,A,perm,&factinfo));
181c4762a1bSJed Brown       }
182*9566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(sC,A,&factinfo));
183c4762a1bSJed Brown 
184*9566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,b));
185*9566063dSJacob Faibussowitsch       PetscCall(MatSolve(sC,b,y));
186*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sC));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown       /* Check the residual */
189*9566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y,neg_one,x));
190*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(y,NORM_2,&norm2));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown       if (displ) {
193*9566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2));
194c4762a1bSJed Brown       }
195c4762a1bSJed Brown     }
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Test baij matrix A */
199c4762a1bSJed Brown   if (TestBAIJ) {
200c4762a1bSJed Brown     if (displ) {
201*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"BAIJ: \n"));
202c4762a1bSJed Brown     }
203c4762a1bSJed Brown     i = 0;
204c4762a1bSJed Brown     for (lvl=-1; lvl<10; lvl++) {
205c4762a1bSJed Brown       if (lvl==-1) {  /* Cholesky factor */
206c4762a1bSJed Brown         factinfo.fill = 5.0;
207c4762a1bSJed Brown 
208*9566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC));
209*9566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(sC,A,perm,&factinfo));
210c4762a1bSJed Brown       } else {       /* incomplete Cholesky factor */
211c4762a1bSJed Brown         factinfo.fill   = 5.0;
212c4762a1bSJed Brown         factinfo.levels = lvl;
213c4762a1bSJed Brown 
214*9566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC));
215*9566063dSJacob Faibussowitsch         PetscCall(MatICCFactorSymbolic(sC,A,perm,&factinfo));
216c4762a1bSJed Brown       }
217*9566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(sC,A,&factinfo));
218c4762a1bSJed Brown 
219*9566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,b));
220*9566063dSJacob Faibussowitsch       PetscCall(MatSolve(sC,b,y));
221*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sC));
222c4762a1bSJed Brown 
223c4762a1bSJed Brown       /* Check the residual */
224*9566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y,neg_one,x));
225*9566063dSJacob Faibussowitsch       PetscCall(VecNorm(y,NORM_2,&norm2));
226c4762a1bSJed Brown       if (displ) {
227*9566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2));
228c4762a1bSJed Brown       }
229c4762a1bSJed Brown     }
230c4762a1bSJed Brown   }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* Test sbaij matrix sA */
233c4762a1bSJed Brown   if (displ) {
234*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"SBAIJ: \n"));
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown   i = 0;
237c4762a1bSJed Brown   for (lvl=-1; lvl<10; lvl++) {
238c4762a1bSJed Brown     if (lvl==-1) {  /* Cholesky factor */
239c4762a1bSJed Brown       factinfo.fill = 5.0;
240c4762a1bSJed Brown 
241*9566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC));
242*9566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(sC,sA,perm,&factinfo));
243c4762a1bSJed Brown     } else {       /* incomplete Cholesky factor */
244c4762a1bSJed Brown       factinfo.fill   = 5.0;
245c4762a1bSJed Brown       factinfo.levels = lvl;
246c4762a1bSJed Brown 
247*9566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC));
248*9566063dSJacob Faibussowitsch       PetscCall(MatICCFactorSymbolic(sC,sA,perm,&factinfo));
249c4762a1bSJed Brown     }
250*9566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(sC,sA,&factinfo));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown     if (lvl==0 && bs==1) { /* Test inplace ICC(0) for sbaij sA - does not work for new datastructure */
253c4762a1bSJed Brown       /*
254c4762a1bSJed Brown         Mat B;
255*9566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&B));
256*9566063dSJacob Faibussowitsch         PetscCall(MatICCFactor(B,perm,&factinfo));
257*9566063dSJacob Faibussowitsch         PetscCall(MatEqual(sC,B,&equal));
258c4762a1bSJed Brown         if (!equal) {
259c4762a1bSJed Brown           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor");
260c4762a1bSJed Brown         }
261*9566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&B));
262c4762a1bSJed Brown       */
263c4762a1bSJed Brown     }
264c4762a1bSJed Brown 
265*9566063dSJacob Faibussowitsch     PetscCall(MatMult(sA,x,b));
266*9566063dSJacob Faibussowitsch     PetscCall(MatSolve(sC,b,y));
267c4762a1bSJed Brown 
268c4762a1bSJed Brown     /* Test MatSolves() */
269c4762a1bSJed Brown     if (bs == 1) {
270c4762a1bSJed Brown       Vecs xx,bb;
271*9566063dSJacob Faibussowitsch       PetscCall(VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx));
272*9566063dSJacob Faibussowitsch       PetscCall(VecsDuplicate(xx,&bb));
273*9566063dSJacob Faibussowitsch       PetscCall(MatSolves(sC,bb,xx));
274*9566063dSJacob Faibussowitsch       PetscCall(VecsDestroy(xx));
275*9566063dSJacob Faibussowitsch       PetscCall(VecsDestroy(bb));
276c4762a1bSJed Brown     }
277*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&sC));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown     /* Check the residual */
280*9566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y,neg_one,x));
281*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(y,NORM_2,&norm2));
282c4762a1bSJed Brown     if (displ) {
283*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2));
284c4762a1bSJed Brown     }
285c4762a1bSJed Brown   }
286c4762a1bSJed Brown 
287*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
288*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
289*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
290*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
291*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
292*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
293*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
294c4762a1bSJed Brown 
295*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
296b122ec5aSJacob Faibussowitsch   return 0;
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown /*TEST
300c4762a1bSJed Brown 
301c4762a1bSJed Brown    test:
302c4762a1bSJed Brown       args: -bs {{1 2 3 4 5 6 7 8}}
303c4762a1bSJed Brown 
304c4762a1bSJed Brown    test:
305c4762a1bSJed Brown       suffix: 3
306c4762a1bSJed Brown       args: -testaij
307c4762a1bSJed Brown       output_file: output/ex76_1.out
308c4762a1bSJed Brown 
309c4762a1bSJed Brown TEST*/
310