xref: /petsc/src/mat/tests/ex76.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
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 
23327415f7SBarry Smith   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
26be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-reorder", &reorder, NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-testaij", &TestAIJ, NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-testShift", &TestShift, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-displ", &displ, NULL));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   n = mbs * bs;
35c4762a1bSJed Brown   if (TestAIJ) { /* A is in aij format */
369566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, n, n, nz, NULL, &A));
37c4762a1bSJed Brown     TestBAIJ = PETSC_FALSE;
38c4762a1bSJed Brown   } else { /* A is in baij format */
399566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD, bs, n, n, nz, NULL, &A));
40c4762a1bSJed Brown     TestAIJ = PETSC_FALSE;
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Assemble matrix */
44c4762a1bSJed Brown   if (bs == 1) {
459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_problem", &prob, NULL));
46c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
479371c9d4SSatish Balay       value[0] = -1.0;
489371c9d4SSatish Balay       value[1] = 2.0;
499371c9d4SSatish Balay       value[2] = -1.0;
50c4762a1bSJed Brown       for (i = 1; i < n - 1; i++) {
519371c9d4SSatish Balay         col[0] = i - 1;
529371c9d4SSatish Balay         col[1] = i;
539371c9d4SSatish Balay         col[2] = i + 1;
549566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
55c4762a1bSJed Brown       }
569371c9d4SSatish Balay       i      = n - 1;
579371c9d4SSatish Balay       col[0] = 0;
589371c9d4SSatish Balay       col[1] = n - 2;
599371c9d4SSatish Balay       col[2] = n - 1;
60c4762a1bSJed Brown 
619371c9d4SSatish Balay       value[0] = 0.1;
629371c9d4SSatish Balay       value[1] = -1;
639371c9d4SSatish Balay       value[2] = 2;
649566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
65c4762a1bSJed Brown 
669371c9d4SSatish Balay       i      = 0;
679371c9d4SSatish Balay       col[0] = 0;
689371c9d4SSatish Balay       col[1] = 1;
699371c9d4SSatish Balay       col[2] = n - 1;
70c4762a1bSJed Brown 
719371c9d4SSatish Balay       value[0] = 2.0;
729371c9d4SSatish Balay       value[1] = -1.0;
739371c9d4SSatish Balay       value[2] = 0.1;
749566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
75c4762a1bSJed Brown     } else if (prob == 2) { /* matrix for the five point stencil */
76c4762a1bSJed Brown       n1 = (PetscInt)(PetscSqrtReal((PetscReal)n) + 0.001);
77bc30f867SBarry Smith       PetscCheck(n1 * n1 == n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "sqrt(n) must be a positive integer!");
78c4762a1bSJed Brown       for (i = 0; i < n1; i++) {
79c4762a1bSJed Brown         for (j = 0; j < n1; j++) {
80c4762a1bSJed Brown           Ii = j + n1 * i;
81c4762a1bSJed Brown           if (i > 0) {
82c4762a1bSJed Brown             J = Ii - n1;
839566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
84c4762a1bSJed Brown           }
85c4762a1bSJed Brown           if (i < n1 - 1) {
86c4762a1bSJed Brown             J = Ii + n1;
879566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
88c4762a1bSJed Brown           }
89c4762a1bSJed Brown           if (j > 0) {
90c4762a1bSJed Brown             J = Ii - 1;
919566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
92c4762a1bSJed Brown           }
93c4762a1bSJed Brown           if (j < n1 - 1) {
94c4762a1bSJed Brown             J = Ii + 1;
959566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
96c4762a1bSJed Brown           }
979566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &four, INSERT_VALUES));
98c4762a1bSJed Brown         }
99c4762a1bSJed Brown       }
100c4762a1bSJed Brown     }
101c4762a1bSJed Brown   } else { /* bs > 1 */
102c4762a1bSJed Brown     for (block = 0; block < n / bs; block++) {
103c4762a1bSJed Brown       /* diagonal blocks */
1049371c9d4SSatish Balay       value[0] = -1.0;
1059371c9d4SSatish Balay       value[1] = 4.0;
1069371c9d4SSatish Balay       value[2] = -1.0;
107c4762a1bSJed Brown       for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) {
1089371c9d4SSatish Balay         col[0] = i - 1;
1099371c9d4SSatish Balay         col[1] = i;
1109371c9d4SSatish Balay         col[2] = i + 1;
1119566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
112c4762a1bSJed Brown       }
1139371c9d4SSatish Balay       i      = bs - 1 + block * bs;
1149371c9d4SSatish Balay       col[0] = bs - 2 + block * bs;
1159371c9d4SSatish Balay       col[1] = bs - 1 + block * bs;
116c4762a1bSJed Brown 
1179371c9d4SSatish Balay       value[0] = -1.0;
1189371c9d4SSatish Balay       value[1] = 4.0;
1199566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
120c4762a1bSJed Brown 
1219371c9d4SSatish Balay       i      = 0 + block * bs;
1229371c9d4SSatish Balay       col[0] = 0 + block * bs;
1239371c9d4SSatish Balay       col[1] = 1 + block * bs;
124c4762a1bSJed Brown 
1259371c9d4SSatish Balay       value[0] = 4.0;
1269371c9d4SSatish Balay       value[1] = -1.0;
1279566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown     /* off-diagonal blocks */
130c4762a1bSJed Brown     value[0] = -1.0;
131c4762a1bSJed Brown     for (i = 0; i < (n / bs - 1) * bs; i++) {
132c4762a1bSJed Brown       col[0] = i + bs;
1339566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 1, col, value, INSERT_VALUES));
1349371c9d4SSatish Balay       col[0] = i;
1359371c9d4SSatish Balay       row    = i + bs;
1369566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES));
137c4762a1bSJed Brown     }
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   if (TestShift) {
141c4762a1bSJed Brown     /* set diagonals in the 0-th block as 0 for testing shift numerical factor */
142c4762a1bSJed Brown     for (i = 0; i < bs; i++) {
1439371c9d4SSatish Balay       row      = i;
1449371c9d4SSatish Balay       col[0]   = i;
1459371c9d4SSatish Balay       value[0] = 0.0;
1469566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES));
147c4762a1bSJed Brown     }
148c4762a1bSJed Brown   }
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* Test MatConvert */
1549566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
1559566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &sA));
1569566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A, sA, 20, &equal));
15728b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_USER, "A != sA");
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
1609566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &Ii, &J));
1619566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(sA, &i, &j));
162bc30f867SBarry Smith   PetscCheck(i == Ii && j == J, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetOwnershipRange() in MatSBAIJ format");
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* Vectors */
1659566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
1669566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
1679566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x));
1689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &b));
1699566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
1709566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x, rdm));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* Test MatReordering() - not work on sbaij matrix */
173c4762a1bSJed Brown   if (reorder) {
1749566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A, MATORDERINGRCM, &perm, &cperm));
175c4762a1bSJed Brown   } else {
1769566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A, MATORDERINGNATURAL, &perm, &cperm));
177c4762a1bSJed Brown   }
1789566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /* initialize factinfo */
1819566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&factinfo));
182c4762a1bSJed Brown   if (TestShift == 1) {
183c4762a1bSJed Brown     factinfo.shifttype   = (PetscReal)MAT_SHIFT_NONZERO;
184c4762a1bSJed Brown     factinfo.shiftamount = 0.1;
185c4762a1bSJed Brown   } else if (TestShift == 2) {
186c4762a1bSJed Brown     factinfo.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
187c4762a1bSJed Brown   }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Test MatCholeskyFactor(), MatICCFactor() */
190c4762a1bSJed Brown   /*------------------------------------------*/
191c4762a1bSJed Brown   /* Test aij matrix A */
192c4762a1bSJed Brown   if (TestAIJ) {
19348a46eb9SPierre Jolivet     if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "AIJ: \n"));
194c4762a1bSJed Brown     i = 0;
195c4762a1bSJed Brown     for (lvl = -1; lvl < 10; lvl++) {
196c4762a1bSJed Brown       if (lvl == -1) { /* Cholesky factor */
197c4762a1bSJed Brown         factinfo.fill = 5.0;
198c4762a1bSJed Brown 
1999566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sC));
2009566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(sC, A, perm, &factinfo));
201c4762a1bSJed Brown       } else { /* incomplete Cholesky factor */
202c4762a1bSJed Brown         factinfo.fill   = 5.0;
203c4762a1bSJed Brown         factinfo.levels = lvl;
204c4762a1bSJed Brown 
2059566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_ICC, &sC));
2069566063dSJacob Faibussowitsch         PetscCall(MatICCFactorSymbolic(sC, A, perm, &factinfo));
207c4762a1bSJed Brown       }
2089566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(sC, A, &factinfo));
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch       PetscCall(MatMult(A, x, b));
2119566063dSJacob Faibussowitsch       PetscCall(MatSolve(sC, b, y));
2129566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sC));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown       /* Check the residual */
2159566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, neg_one, x));
2169566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &norm2));
217c4762a1bSJed Brown 
21848a46eb9SPierre Jolivet       if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  lvl: %" PetscInt_FMT ", residual: %g\n", lvl, (double)norm2));
219c4762a1bSJed Brown     }
220c4762a1bSJed Brown   }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /* Test baij matrix A */
223c4762a1bSJed Brown   if (TestBAIJ) {
22448a46eb9SPierre Jolivet     if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BAIJ: \n"));
225c4762a1bSJed Brown     i = 0;
226c4762a1bSJed Brown     for (lvl = -1; lvl < 10; lvl++) {
227c4762a1bSJed Brown       if (lvl == -1) { /* Cholesky factor */
228c4762a1bSJed Brown         factinfo.fill = 5.0;
229c4762a1bSJed Brown 
2309566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sC));
2319566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(sC, A, perm, &factinfo));
232c4762a1bSJed Brown       } else { /* incomplete Cholesky factor */
233c4762a1bSJed Brown         factinfo.fill   = 5.0;
234c4762a1bSJed Brown         factinfo.levels = lvl;
235c4762a1bSJed Brown 
2369566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_ICC, &sC));
2379566063dSJacob Faibussowitsch         PetscCall(MatICCFactorSymbolic(sC, A, perm, &factinfo));
238c4762a1bSJed Brown       }
2399566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(sC, A, &factinfo));
240c4762a1bSJed Brown 
2419566063dSJacob Faibussowitsch       PetscCall(MatMult(A, x, b));
2429566063dSJacob Faibussowitsch       PetscCall(MatSolve(sC, b, y));
2439566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sC));
244c4762a1bSJed Brown 
245c4762a1bSJed Brown       /* Check the residual */
2469566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, neg_one, x));
2479566063dSJacob Faibussowitsch       PetscCall(VecNorm(y, NORM_2, &norm2));
24848a46eb9SPierre Jolivet       if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  lvl: %" PetscInt_FMT ", residual: %g\n", lvl, (double)norm2));
249c4762a1bSJed Brown     }
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* Test sbaij matrix sA */
25348a46eb9SPierre Jolivet   if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SBAIJ: \n"));
254c4762a1bSJed Brown   i = 0;
255c4762a1bSJed Brown   for (lvl = -1; lvl < 10; lvl++) {
256c4762a1bSJed Brown     if (lvl == -1) { /* Cholesky factor */
257c4762a1bSJed Brown       factinfo.fill = 5.0;
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sA, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sC));
2609566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(sC, sA, perm, &factinfo));
261c4762a1bSJed Brown     } else { /* incomplete Cholesky factor */
262c4762a1bSJed Brown       factinfo.fill   = 5.0;
263c4762a1bSJed Brown       factinfo.levels = lvl;
264c4762a1bSJed Brown 
2659566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sA, MATSOLVERPETSC, MAT_FACTOR_ICC, &sC));
2669566063dSJacob Faibussowitsch       PetscCall(MatICCFactorSymbolic(sC, sA, perm, &factinfo));
267c4762a1bSJed Brown     }
2689566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(sC, sA, &factinfo));
269c4762a1bSJed Brown 
270c4762a1bSJed Brown     if (lvl == 0 && bs == 1) { /* Test inplace ICC(0) for sbaij sA - does not work for new datastructure */
271c4762a1bSJed Brown       /*
272c4762a1bSJed Brown         Mat B;
2739566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&B));
2749566063dSJacob Faibussowitsch         PetscCall(MatICCFactor(B,perm,&factinfo));
2759566063dSJacob Faibussowitsch         PetscCall(MatEqual(sC,B,&equal));
276c4762a1bSJed Brown         if (!equal) {
277c4762a1bSJed Brown           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor");
278c4762a1bSJed Brown         }
2799566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&B));
280c4762a1bSJed Brown       */
281c4762a1bSJed Brown     }
282c4762a1bSJed Brown 
2839566063dSJacob Faibussowitsch     PetscCall(MatMult(sA, x, b));
2849566063dSJacob Faibussowitsch     PetscCall(MatSolve(sC, b, y));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown     /* Test MatSolves() */
287c4762a1bSJed Brown     if (bs == 1) {
288c4762a1bSJed Brown       Vecs xx, bb;
2899566063dSJacob Faibussowitsch       PetscCall(VecsCreateSeq(PETSC_COMM_SELF, n, 4, &xx));
2909566063dSJacob Faibussowitsch       PetscCall(VecsDuplicate(xx, &bb));
2919566063dSJacob Faibussowitsch       PetscCall(MatSolves(sC, bb, xx));
2929566063dSJacob Faibussowitsch       PetscCall(VecsDestroy(xx));
2939566063dSJacob Faibussowitsch       PetscCall(VecsDestroy(bb));
294c4762a1bSJed Brown     }
2959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&sC));
296c4762a1bSJed Brown 
297c4762a1bSJed Brown     /* Check the residual */
2989566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, neg_one, x));
2999566063dSJacob Faibussowitsch     PetscCall(VecNorm(y, NORM_2, &norm2));
30048a46eb9SPierre Jolivet     if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  lvl: %" PetscInt_FMT ", residual: %g\n", lvl, (double)norm2));
301c4762a1bSJed Brown   }
302c4762a1bSJed Brown 
3039566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
3049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
3069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
3089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
3099566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
310c4762a1bSJed Brown 
3119566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
312b122ec5aSJacob Faibussowitsch   return 0;
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown /*TEST
316c4762a1bSJed Brown 
317c4762a1bSJed Brown    test:
318c4762a1bSJed Brown       args: -bs {{1 2 3 4 5 6 7 8}}
319c4762a1bSJed Brown 
320c4762a1bSJed Brown    test:
321c4762a1bSJed Brown       suffix: 3
322c4762a1bSJed Brown       args: -testaij
323c4762a1bSJed Brown       output_file: output/ex76_1.out
324c4762a1bSJed Brown 
325c4762a1bSJed Brown TEST*/
326