1c4762a1bSJed Brown static char help[] = "Tests cholesky, icc factorization and solve on sequential aij, baij and sbaij matrices. \n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Vec x, y, b;
8c4762a1bSJed Brown Mat A; /* linear system matrix */
9c4762a1bSJed Brown Mat sA, sC; /* symmetric part of the matrices */
10c4762a1bSJed Brown PetscInt n, mbs = 16, bs = 1, nz = 3, prob = 1, i, j, col[3], block, row, Ii, J, n1, lvl;
11c4762a1bSJed Brown PetscMPIInt size;
12c4762a1bSJed Brown PetscReal norm2;
13c4762a1bSJed Brown PetscScalar neg_one = -1.0, four = 4.0, value[3];
14c4762a1bSJed Brown IS perm, cperm;
15c4762a1bSJed Brown PetscRandom rdm;
16c4762a1bSJed Brown PetscBool reorder = PETSC_FALSE, displ = PETSC_FALSE;
17c4762a1bSJed Brown MatFactorInfo factinfo;
18c4762a1bSJed Brown PetscBool equal;
19c4762a1bSJed Brown PetscBool TestAIJ = PETSC_FALSE, TestBAIJ = PETSC_TRUE;
20c4762a1bSJed Brown PetscInt TestShift = 0;
21c4762a1bSJed Brown
22327415f7SBarry Smith PetscFunctionBeginUser;
23c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
25be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL));
289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-reorder", &reorder, NULL));
299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-testaij", &TestAIJ, NULL));
309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-testShift", &TestShift, NULL));
319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-displ", &displ, NULL));
32c4762a1bSJed Brown
33c4762a1bSJed Brown n = mbs * bs;
34c4762a1bSJed Brown if (TestAIJ) { /* A is in aij format */
359566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, n, n, nz, NULL, &A));
36c4762a1bSJed Brown TestBAIJ = PETSC_FALSE;
37c4762a1bSJed Brown } else { /* A is in baij format */
389566063dSJacob 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) {
449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_problem", &prob, NULL));
45c4762a1bSJed Brown if (prob == 1) { /* tridiagonal matrix */
469371c9d4SSatish Balay value[0] = -1.0;
479371c9d4SSatish Balay value[1] = 2.0;
489371c9d4SSatish Balay value[2] = -1.0;
49c4762a1bSJed Brown for (i = 1; i < n - 1; i++) {
509371c9d4SSatish Balay col[0] = i - 1;
519371c9d4SSatish Balay col[1] = i;
529371c9d4SSatish Balay col[2] = i + 1;
539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
54c4762a1bSJed Brown }
559371c9d4SSatish Balay i = n - 1;
569371c9d4SSatish Balay col[0] = 0;
579371c9d4SSatish Balay col[1] = n - 2;
589371c9d4SSatish Balay col[2] = n - 1;
59c4762a1bSJed Brown
609371c9d4SSatish Balay value[0] = 0.1;
619371c9d4SSatish Balay value[1] = -1;
629371c9d4SSatish Balay value[2] = 2;
639566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
64c4762a1bSJed Brown
659371c9d4SSatish Balay i = 0;
669371c9d4SSatish Balay col[0] = 0;
679371c9d4SSatish Balay col[1] = 1;
689371c9d4SSatish Balay col[2] = n - 1;
69c4762a1bSJed Brown
709371c9d4SSatish Balay value[0] = 2.0;
719371c9d4SSatish Balay value[1] = -1.0;
729371c9d4SSatish Balay value[2] = 0.1;
739566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
74c4762a1bSJed Brown } else if (prob == 2) { /* matrix for the five point stencil */
75c4762a1bSJed Brown n1 = (PetscInt)(PetscSqrtReal((PetscReal)n) + 0.001);
76bc30f867SBarry Smith PetscCheck(n1 * n1 == n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "sqrt(n) must be a positive integer!");
77c4762a1bSJed Brown for (i = 0; i < n1; i++) {
78c4762a1bSJed Brown for (j = 0; j < n1; j++) {
79c4762a1bSJed Brown Ii = j + n1 * i;
80c4762a1bSJed Brown if (i > 0) {
81c4762a1bSJed Brown J = Ii - n1;
829566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
83c4762a1bSJed Brown }
84c4762a1bSJed Brown if (i < n1 - 1) {
85c4762a1bSJed Brown J = Ii + n1;
869566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
87c4762a1bSJed Brown }
88c4762a1bSJed Brown if (j > 0) {
89c4762a1bSJed Brown J = Ii - 1;
909566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
91c4762a1bSJed Brown }
92c4762a1bSJed Brown if (j < n1 - 1) {
93c4762a1bSJed Brown J = Ii + 1;
949566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
95c4762a1bSJed Brown }
969566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &four, INSERT_VALUES));
97c4762a1bSJed Brown }
98c4762a1bSJed Brown }
99c4762a1bSJed Brown }
100c4762a1bSJed Brown } else { /* bs > 1 */
101c4762a1bSJed Brown for (block = 0; block < n / bs; block++) {
102c4762a1bSJed Brown /* diagonal blocks */
1039371c9d4SSatish Balay value[0] = -1.0;
1049371c9d4SSatish Balay value[1] = 4.0;
1059371c9d4SSatish Balay value[2] = -1.0;
106c4762a1bSJed Brown for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) {
1079371c9d4SSatish Balay col[0] = i - 1;
1089371c9d4SSatish Balay col[1] = i;
1099371c9d4SSatish Balay col[2] = i + 1;
1109566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
111c4762a1bSJed Brown }
1129371c9d4SSatish Balay i = bs - 1 + block * bs;
1139371c9d4SSatish Balay col[0] = bs - 2 + block * bs;
1149371c9d4SSatish Balay col[1] = bs - 1 + block * bs;
115c4762a1bSJed Brown
1169371c9d4SSatish Balay value[0] = -1.0;
1179371c9d4SSatish Balay value[1] = 4.0;
1189566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
119c4762a1bSJed Brown
1209371c9d4SSatish Balay i = 0 + block * bs;
1219371c9d4SSatish Balay col[0] = 0 + block * bs;
1229371c9d4SSatish Balay col[1] = 1 + block * bs;
123c4762a1bSJed Brown
1249371c9d4SSatish Balay value[0] = 4.0;
1259371c9d4SSatish Balay value[1] = -1.0;
1269566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
127c4762a1bSJed Brown }
128c4762a1bSJed Brown /* off-diagonal blocks */
129c4762a1bSJed Brown value[0] = -1.0;
130c4762a1bSJed Brown for (i = 0; i < (n / bs - 1) * bs; i++) {
131c4762a1bSJed Brown col[0] = i + bs;
1329566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, col, value, INSERT_VALUES));
1339371c9d4SSatish Balay col[0] = i;
1349371c9d4SSatish Balay row = i + bs;
1359566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES));
136c4762a1bSJed Brown }
137c4762a1bSJed Brown }
138c4762a1bSJed Brown
139c4762a1bSJed Brown if (TestShift) {
140c4762a1bSJed Brown /* set diagonals in the 0-th block as 0 for testing shift numerical factor */
141c4762a1bSJed Brown for (i = 0; i < bs; i++) {
1429371c9d4SSatish Balay row = i;
1439371c9d4SSatish Balay col[0] = i;
1449371c9d4SSatish Balay value[0] = 0.0;
1459566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES));
146c4762a1bSJed Brown }
147c4762a1bSJed Brown }
148c4762a1bSJed Brown
1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
151c4762a1bSJed Brown
152c4762a1bSJed Brown /* Test MatConvert */
1539566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
1549566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &sA));
1559566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, sA, 20, &equal));
15628b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_USER, "A != sA");
157c4762a1bSJed Brown
158c4762a1bSJed Brown /* Test MatGetOwnershipRange() */
1599566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &Ii, &J));
1609566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(sA, &i, &j));
161bc30f867SBarry Smith PetscCheck(i == Ii && j == J, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatGetOwnershipRange() in MatSBAIJ format");
162c4762a1bSJed Brown
163c4762a1bSJed Brown /* Vectors */
1649566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
1659566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm));
1669566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x));
1679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &b));
1689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y));
1699566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm));
170c4762a1bSJed Brown
171c4762a1bSJed Brown /* Test MatReordering() - not work on sbaij matrix */
172c4762a1bSJed Brown if (reorder) {
1739566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A, MATORDERINGRCM, &perm, &cperm));
174c4762a1bSJed Brown } else {
1759566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A, MATORDERINGNATURAL, &perm, &cperm));
176c4762a1bSJed Brown }
1779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cperm));
178c4762a1bSJed Brown
179c4762a1bSJed Brown /* initialize factinfo */
1809566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&factinfo));
181c4762a1bSJed Brown if (TestShift == 1) {
182c4762a1bSJed Brown factinfo.shifttype = (PetscReal)MAT_SHIFT_NONZERO;
183c4762a1bSJed Brown factinfo.shiftamount = 0.1;
184c4762a1bSJed Brown } else if (TestShift == 2) {
185c4762a1bSJed Brown factinfo.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
186c4762a1bSJed Brown }
187c4762a1bSJed Brown
188c4762a1bSJed Brown /* Test MatCholeskyFactor(), MatICCFactor() */
189c4762a1bSJed Brown /*------------------------------------------*/
190c4762a1bSJed Brown /* Test aij matrix A */
191c4762a1bSJed Brown if (TestAIJ) {
19248a46eb9SPierre Jolivet if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "AIJ: \n"));
193c4762a1bSJed Brown i = 0;
194c4762a1bSJed Brown for (lvl = -1; lvl < 10; lvl++) {
195c4762a1bSJed Brown if (lvl == -1) { /* Cholesky factor */
196c4762a1bSJed Brown factinfo.fill = 5.0;
197c4762a1bSJed Brown
1989566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sC));
1999566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(sC, A, perm, &factinfo));
200c4762a1bSJed Brown } else { /* incomplete Cholesky factor */
201c4762a1bSJed Brown factinfo.fill = 5.0;
202c4762a1bSJed Brown factinfo.levels = lvl;
203c4762a1bSJed Brown
2049566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_ICC, &sC));
2059566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(sC, A, perm, &factinfo));
206c4762a1bSJed Brown }
2079566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(sC, A, &factinfo));
208c4762a1bSJed Brown
2099566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, b));
2109566063dSJacob Faibussowitsch PetscCall(MatSolve(sC, b, y));
2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sC));
212c4762a1bSJed Brown
213c4762a1bSJed Brown /* Check the residual */
2149566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, neg_one, x));
2159566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2));
216c4762a1bSJed Brown
21748a46eb9SPierre Jolivet if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " lvl: %" PetscInt_FMT ", residual: %g\n", lvl, (double)norm2));
218c4762a1bSJed Brown }
219c4762a1bSJed Brown }
220c4762a1bSJed Brown
221c4762a1bSJed Brown /* Test baij matrix A */
222c4762a1bSJed Brown if (TestBAIJ) {
22348a46eb9SPierre Jolivet if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "BAIJ: \n"));
224c4762a1bSJed Brown i = 0;
225c4762a1bSJed Brown for (lvl = -1; lvl < 10; lvl++) {
226c4762a1bSJed Brown if (lvl == -1) { /* Cholesky factor */
227c4762a1bSJed Brown factinfo.fill = 5.0;
228c4762a1bSJed Brown
2299566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sC));
2309566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(sC, A, perm, &factinfo));
231c4762a1bSJed Brown } else { /* incomplete Cholesky factor */
232c4762a1bSJed Brown factinfo.fill = 5.0;
233c4762a1bSJed Brown factinfo.levels = lvl;
234c4762a1bSJed Brown
2359566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_ICC, &sC));
2369566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(sC, A, perm, &factinfo));
237c4762a1bSJed Brown }
2389566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(sC, A, &factinfo));
239c4762a1bSJed Brown
2409566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, b));
2419566063dSJacob Faibussowitsch PetscCall(MatSolve(sC, b, y));
2429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sC));
243c4762a1bSJed Brown
244c4762a1bSJed Brown /* Check the residual */
2459566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, neg_one, x));
2469566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2));
24748a46eb9SPierre Jolivet if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " lvl: %" PetscInt_FMT ", residual: %g\n", lvl, (double)norm2));
248c4762a1bSJed Brown }
249c4762a1bSJed Brown }
250c4762a1bSJed Brown
251c4762a1bSJed Brown /* Test sbaij matrix sA */
25248a46eb9SPierre Jolivet if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SBAIJ: \n"));
253c4762a1bSJed Brown i = 0;
254c4762a1bSJed Brown for (lvl = -1; lvl < 10; lvl++) {
255c4762a1bSJed Brown if (lvl == -1) { /* Cholesky factor */
256c4762a1bSJed Brown factinfo.fill = 5.0;
257c4762a1bSJed Brown
2589566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sA, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sC));
2599566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(sC, sA, perm, &factinfo));
260c4762a1bSJed Brown } else { /* incomplete Cholesky factor */
261c4762a1bSJed Brown factinfo.fill = 5.0;
262c4762a1bSJed Brown factinfo.levels = lvl;
263c4762a1bSJed Brown
2649566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sA, MATSOLVERPETSC, MAT_FACTOR_ICC, &sC));
2659566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(sC, sA, perm, &factinfo));
266c4762a1bSJed Brown }
2679566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(sC, sA, &factinfo));
268c4762a1bSJed Brown
269c4762a1bSJed Brown if (lvl == 0 && bs == 1) { /* Test inplace ICC(0) for sbaij sA - does not work for new datastructure */
270c4762a1bSJed Brown /*
271c4762a1bSJed Brown Mat B;
2729566063dSJacob Faibussowitsch PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&B));
2739566063dSJacob Faibussowitsch PetscCall(MatICCFactor(B,perm,&factinfo));
2749566063dSJacob Faibussowitsch PetscCall(MatEqual(sC,B,&equal));
275*966bd95aSPierre Jolivet PetscCheck(equal, PETSC_COMM_SELF,PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor");
2769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
277c4762a1bSJed Brown */
278c4762a1bSJed Brown }
279c4762a1bSJed Brown
2809566063dSJacob Faibussowitsch PetscCall(MatMult(sA, x, b));
2819566063dSJacob Faibussowitsch PetscCall(MatSolve(sC, b, y));
282c4762a1bSJed Brown
283c4762a1bSJed Brown /* Test MatSolves() */
284c4762a1bSJed Brown if (bs == 1) {
285c4762a1bSJed Brown Vecs xx, bb;
2869566063dSJacob Faibussowitsch PetscCall(VecsCreateSeq(PETSC_COMM_SELF, n, 4, &xx));
2879566063dSJacob Faibussowitsch PetscCall(VecsDuplicate(xx, &bb));
2889566063dSJacob Faibussowitsch PetscCall(MatSolves(sC, bb, xx));
2899566063dSJacob Faibussowitsch PetscCall(VecsDestroy(xx));
2909566063dSJacob Faibussowitsch PetscCall(VecsDestroy(bb));
291c4762a1bSJed Brown }
2929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sC));
293c4762a1bSJed Brown
294c4762a1bSJed Brown /* Check the residual */
2959566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, neg_one, x));
2969566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2));
29748a46eb9SPierre Jolivet if (displ) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " lvl: %" PetscInt_FMT ", residual: %g\n", lvl, (double)norm2));
298c4762a1bSJed Brown }
299c4762a1bSJed Brown
3009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm));
3019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
3029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA));
3039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
3049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
3059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
3069566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm));
307c4762a1bSJed Brown
3089566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
309b122ec5aSJacob Faibussowitsch return 0;
310c4762a1bSJed Brown }
311c4762a1bSJed Brown
312c4762a1bSJed Brown /*TEST
313c4762a1bSJed Brown
314c4762a1bSJed Brown test:
315c4762a1bSJed Brown args: -bs {{1 2 3 4 5 6 7 8}}
3163886731fSPierre Jolivet output_file: output/empty.out
317c4762a1bSJed Brown
318c4762a1bSJed Brown test:
319c4762a1bSJed Brown suffix: 3
320c4762a1bSJed Brown args: -testaij
3213886731fSPierre Jolivet output_file: output/empty.out
322c4762a1bSJed Brown
323c4762a1bSJed Brown TEST*/
324