xref: /petsc/src/mat/tests/ex76.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry 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 */
47c4762a1bSJed Brown       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
48c4762a1bSJed Brown       for (i=1; i<n-1; i++) {
49c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
509566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
51c4762a1bSJed Brown       }
52c4762a1bSJed Brown       i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown       value[0]= 0.1; value[1]=-1; value[2]=2;
559566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown       i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown       value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
609566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
61c4762a1bSJed Brown     } else if (prob ==2) { /* matrix for the five point stencil */
62c4762a1bSJed Brown       n1 = (PetscInt) (PetscSqrtReal((PetscReal)n) + 0.001);
63bc30f867SBarry Smith       PetscCheck(n1*n1 == n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive integer!");
64c4762a1bSJed Brown       for (i=0; i<n1; i++) {
65c4762a1bSJed Brown         for (j=0; j<n1; j++) {
66c4762a1bSJed Brown           Ii = j + n1*i;
67c4762a1bSJed Brown           if (i>0) {
68c4762a1bSJed Brown             J    = Ii - n1;
699566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
70c4762a1bSJed Brown           }
71c4762a1bSJed Brown           if (i<n1-1) {
72c4762a1bSJed Brown             J    = Ii + n1;
739566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
74c4762a1bSJed Brown           }
75c4762a1bSJed Brown           if (j>0) {
76c4762a1bSJed Brown             J    = Ii - 1;
779566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
78c4762a1bSJed Brown           }
79c4762a1bSJed Brown           if (j<n1-1) {
80c4762a1bSJed Brown             J    = Ii + 1;
819566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));
82c4762a1bSJed Brown           }
839566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES));
84c4762a1bSJed Brown         }
85c4762a1bSJed Brown       }
86c4762a1bSJed Brown     }
87c4762a1bSJed Brown   } else { /* bs > 1 */
88c4762a1bSJed Brown     for (block=0; block<n/bs; block++) {
89c4762a1bSJed Brown       /* diagonal blocks */
90c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
91c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
92c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
939566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
94c4762a1bSJed Brown       }
95c4762a1bSJed Brown       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
989566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
1039566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
104c4762a1bSJed Brown     }
105c4762a1bSJed Brown     /* off-diagonal blocks */
106c4762a1bSJed Brown     value[0]=-1.0;
107c4762a1bSJed Brown     for (i=0; i<(n/bs-1)*bs; i++) {
108c4762a1bSJed Brown       col[0]=i+bs;
1099566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES));
110c4762a1bSJed Brown       col[0]=i; row=i+bs;
1119566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES));
112c4762a1bSJed Brown     }
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   if (TestShift) {
116c4762a1bSJed Brown     /* set diagonals in the 0-th block as 0 for testing shift numerical factor */
117c4762a1bSJed Brown     for (i=0; i<bs; i++) {
118c4762a1bSJed Brown       row  = i; col[0] = i; value[0] = 0.0;
1199566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES));
120c4762a1bSJed Brown     }
121c4762a1bSJed Brown   }
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1249566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* Test MatConvert */
1279566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
1289566063dSJacob Faibussowitsch   PetscCall(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA));
1299566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A,sA,20,&equal));
13028b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_USER,"A != sA");
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
1339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&Ii,&J));
1349566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(sA,&i,&j));
135bc30f867SBarry Smith   PetscCheck(i == Ii && j == J ,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetOwnershipRange() in MatSBAIJ format");
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* Vectors */
1389566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
1399566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
1409566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x));
1419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&b));
1429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&y));
1439566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x,rdm));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Test MatReordering() - not work on sbaij matrix */
146c4762a1bSJed Brown   if (reorder) {
1479566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A,MATORDERINGRCM,&perm,&cperm));
148c4762a1bSJed Brown   } else {
1499566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&cperm));
150c4762a1bSJed Brown   }
1519566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* initialize factinfo */
1549566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&factinfo));
155c4762a1bSJed Brown   if (TestShift == 1) {
156c4762a1bSJed Brown     factinfo.shifttype   = (PetscReal)MAT_SHIFT_NONZERO;
157c4762a1bSJed Brown     factinfo.shiftamount = 0.1;
158c4762a1bSJed Brown   } else if (TestShift == 2) {
159c4762a1bSJed Brown     factinfo.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
160c4762a1bSJed Brown   }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* Test MatCholeskyFactor(), MatICCFactor() */
163c4762a1bSJed Brown   /*------------------------------------------*/
164c4762a1bSJed Brown   /* Test aij matrix A */
165c4762a1bSJed Brown   if (TestAIJ) {
166c4762a1bSJed Brown     if (displ) {
1679566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"AIJ: \n"));
168c4762a1bSJed Brown     }
169c4762a1bSJed Brown     i = 0;
170c4762a1bSJed Brown     for (lvl=-1; lvl<10; lvl++) {
171c4762a1bSJed Brown       if (lvl==-1) {  /* Cholesky factor */
172c4762a1bSJed Brown         factinfo.fill = 5.0;
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC));
1759566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(sC,A,perm,&factinfo));
176c4762a1bSJed Brown       } else {       /* incomplete Cholesky factor */
177c4762a1bSJed Brown         factinfo.fill   = 5.0;
178c4762a1bSJed Brown         factinfo.levels = lvl;
179c4762a1bSJed Brown 
1809566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC));
1819566063dSJacob Faibussowitsch         PetscCall(MatICCFactorSymbolic(sC,A,perm,&factinfo));
182c4762a1bSJed Brown       }
1839566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(sC,A,&factinfo));
184c4762a1bSJed Brown 
1859566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,b));
1869566063dSJacob Faibussowitsch       PetscCall(MatSolve(sC,b,y));
1879566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sC));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown       /* Check the residual */
1909566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y,neg_one,x));
1919566063dSJacob Faibussowitsch       PetscCall(VecNorm(y,NORM_2,&norm2));
192c4762a1bSJed Brown 
193c4762a1bSJed Brown       if (displ) {
1949566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2));
195c4762a1bSJed Brown       }
196c4762a1bSJed Brown     }
197c4762a1bSJed Brown   }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* Test baij matrix A */
200c4762a1bSJed Brown   if (TestBAIJ) {
201c4762a1bSJed Brown     if (displ) {
2029566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"BAIJ: \n"));
203c4762a1bSJed Brown     }
204c4762a1bSJed Brown     i = 0;
205c4762a1bSJed Brown     for (lvl=-1; lvl<10; lvl++) {
206c4762a1bSJed Brown       if (lvl==-1) {  /* Cholesky factor */
207c4762a1bSJed Brown         factinfo.fill = 5.0;
208c4762a1bSJed Brown 
2099566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC));
2109566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(sC,A,perm,&factinfo));
211c4762a1bSJed Brown       } else {       /* incomplete Cholesky factor */
212c4762a1bSJed Brown         factinfo.fill   = 5.0;
213c4762a1bSJed Brown         factinfo.levels = lvl;
214c4762a1bSJed Brown 
2159566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC));
2169566063dSJacob Faibussowitsch         PetscCall(MatICCFactorSymbolic(sC,A,perm,&factinfo));
217c4762a1bSJed Brown       }
2189566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(sC,A,&factinfo));
219c4762a1bSJed Brown 
2209566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,b));
2219566063dSJacob Faibussowitsch       PetscCall(MatSolve(sC,b,y));
2229566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&sC));
223c4762a1bSJed Brown 
224c4762a1bSJed Brown       /* Check the residual */
2259566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y,neg_one,x));
2269566063dSJacob Faibussowitsch       PetscCall(VecNorm(y,NORM_2,&norm2));
227c4762a1bSJed Brown       if (displ) {
2289566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2));
229c4762a1bSJed Brown       }
230c4762a1bSJed Brown     }
231c4762a1bSJed Brown   }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* Test sbaij matrix sA */
234c4762a1bSJed Brown   if (displ) {
2359566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"SBAIJ: \n"));
236c4762a1bSJed Brown   }
237c4762a1bSJed Brown   i = 0;
238c4762a1bSJed Brown   for (lvl=-1; lvl<10; lvl++) {
239c4762a1bSJed Brown     if (lvl==-1) {  /* Cholesky factor */
240c4762a1bSJed Brown       factinfo.fill = 5.0;
241c4762a1bSJed Brown 
2429566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sC));
2439566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(sC,sA,perm,&factinfo));
244c4762a1bSJed Brown     } else {       /* incomplete Cholesky factor */
245c4762a1bSJed Brown       factinfo.fill   = 5.0;
246c4762a1bSJed Brown       factinfo.levels = lvl;
247c4762a1bSJed Brown 
2489566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sA,MATSOLVERPETSC,MAT_FACTOR_ICC,&sC));
2499566063dSJacob Faibussowitsch       PetscCall(MatICCFactorSymbolic(sC,sA,perm,&factinfo));
250c4762a1bSJed Brown     }
2519566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(sC,sA,&factinfo));
252c4762a1bSJed Brown 
253c4762a1bSJed Brown     if (lvl==0 && bs==1) { /* Test inplace ICC(0) for sbaij sA - does not work for new datastructure */
254c4762a1bSJed Brown       /*
255c4762a1bSJed Brown         Mat B;
2569566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&B));
2579566063dSJacob Faibussowitsch         PetscCall(MatICCFactor(B,perm,&factinfo));
2589566063dSJacob Faibussowitsch         PetscCall(MatEqual(sC,B,&equal));
259c4762a1bSJed Brown         if (!equal) {
260c4762a1bSJed Brown           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor");
261c4762a1bSJed Brown         }
2629566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&B));
263c4762a1bSJed Brown       */
264c4762a1bSJed Brown     }
265c4762a1bSJed Brown 
2669566063dSJacob Faibussowitsch     PetscCall(MatMult(sA,x,b));
2679566063dSJacob Faibussowitsch     PetscCall(MatSolve(sC,b,y));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown     /* Test MatSolves() */
270c4762a1bSJed Brown     if (bs == 1) {
271c4762a1bSJed Brown       Vecs xx,bb;
2729566063dSJacob Faibussowitsch       PetscCall(VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx));
2739566063dSJacob Faibussowitsch       PetscCall(VecsDuplicate(xx,&bb));
2749566063dSJacob Faibussowitsch       PetscCall(MatSolves(sC,bb,xx));
2759566063dSJacob Faibussowitsch       PetscCall(VecsDestroy(xx));
2769566063dSJacob Faibussowitsch       PetscCall(VecsDestroy(bb));
277c4762a1bSJed Brown     }
2789566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&sC));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown     /* Check the residual */
2819566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y,neg_one,x));
2829566063dSJacob Faibussowitsch     PetscCall(VecNorm(y,NORM_2,&norm2));
283c4762a1bSJed Brown     if (displ) {
2849566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  lvl: %" PetscInt_FMT ", residual: %g\n", lvl,(double)norm2));
285c4762a1bSJed Brown     }
286c4762a1bSJed Brown   }
287c4762a1bSJed Brown 
2889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
2899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
2919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
2939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
2949566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
295c4762a1bSJed Brown 
2969566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
297b122ec5aSJacob Faibussowitsch   return 0;
298c4762a1bSJed Brown }
299c4762a1bSJed Brown 
300c4762a1bSJed Brown /*TEST
301c4762a1bSJed Brown 
302c4762a1bSJed Brown    test:
303c4762a1bSJed Brown       args: -bs {{1 2 3 4 5 6 7 8}}
304c4762a1bSJed Brown 
305c4762a1bSJed Brown    test:
306c4762a1bSJed Brown       suffix: 3
307c4762a1bSJed Brown       args: -testaij
308c4762a1bSJed Brown       output_file: output/ex76_1.out
309c4762a1bSJed Brown 
310c4762a1bSJed Brown TEST*/
311