xref: /petsc/src/mat/tests/ex114.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatGetRowMax(), MatGetRowMin(), MatGetRowMaxAbs()\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #define M 5
7c4762a1bSJed Brown #define N 6
8c4762a1bSJed Brown 
9c4762a1bSJed Brown int main(int argc,char **args)
10c4762a1bSJed Brown {
11c4762a1bSJed Brown   Mat            A;
12fa213d2fSHong Zhang   Vec            min,max,maxabs,e;
134e879edeSHong Zhang   PetscInt       m,n,j,imin[M],imax[M],imaxabs[M],indices[N],row,testcase=0;
14c4762a1bSJed Brown   PetscScalar    values[N];
151a254869SHong Zhang   PetscMPIInt    size,rank;
16fa213d2fSHong Zhang   PetscReal      enorm;
17c4762a1bSJed Brown 
18*327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-testcase",&testcase,NULL));
23c4762a1bSJed Brown 
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
251a254869SHong Zhang   if (testcase == 1) { /* proc[0] holds entire A and other processes have no entry */
26dd400576SPatrick Sanan     if (rank == 0) {
279566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A,M,N,PETSC_DECIDE,PETSC_DECIDE));
281a254869SHong Zhang     } else {
299566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A,0,0,PETSC_DECIDE,PETSC_DECIDE));
301a254869SHong Zhang     }
311a254869SHong Zhang     testcase = 0;
321a254869SHong Zhang   } else {
339566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
341a254869SHong Zhang   }
359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
369566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
371a254869SHong Zhang 
38dd400576SPatrick Sanan   if (rank == 0) { /* proc[0] sets matrix A */
391a254869SHong Zhang     for (j=0; j<N; j++) indices[j] = j;
401a254869SHong Zhang     switch (testcase) {
411a254869SHong Zhang     case 1: /* see testcast 0 */
421a254869SHong Zhang       break;
431a254869SHong Zhang     case 2:
44c4762a1bSJed Brown       row = 0;
45f07e67edSHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2] = -2.0; values[3] = -4.0; values[4] = 1.0; values[5] = 1.0;
469566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES));
471a254869SHong Zhang       row = 2;
481a254869SHong Zhang       indices[0] = 0;    indices[1] = 3;    indices[2] = 5;
491a254869SHong Zhang       values[0]  = -2.0; values[1]  = -2.0; values[2]  = -2.0;
509566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES));
511a254869SHong Zhang       row = 3;
521a254869SHong Zhang       indices[0] = 0;    indices[1] = 1;    indices[2] = 4;
531a254869SHong Zhang       values[0]  = -2.0; values[1]  = -2.0; values[2]  = -2.0;
549566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES));
55c4762a1bSJed Brown       row = 4;
561a254869SHong Zhang       indices[0] = 0;    indices[1] = 1;    indices[2] = 2;
571a254869SHong Zhang       values[0]  = -2.0; values[1]  = -2.0; values[2]  = -2.0;
589566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES));
591a254869SHong Zhang       break;
601a254869SHong Zhang     case 3:
611a254869SHong Zhang       row = 0;
621a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2] = -2.0;
639566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,3,indices+1,values,INSERT_VALUES));
641a254869SHong Zhang       row = 1;
651a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2] = -2.0;
669566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES));
671a254869SHong Zhang       row = 2;
681a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2]  = -2.0;
699566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES));
701a254869SHong Zhang       row = 3;
711a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2]  = -2.0; values[3] = -1.0;
729566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES));
731a254869SHong Zhang       row = 4;
741a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2]  = -2.0; values[3] = -1.0;
759566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES));
761a254869SHong Zhang       break;
771a254869SHong Zhang 
781a254869SHong Zhang     default:
791a254869SHong Zhang       row  = 0;
801a254869SHong Zhang       values[0]  = -1.0; values[1] = 0.0; values[2] = 1.0; values[3] = 3.0; values[4] = 4.0; values[5] = -5.0;
819566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES));
821a254869SHong Zhang       row  = 1;
839566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES));
841a254869SHong Zhang       row  = 3;
859566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,indices+4,values+4,INSERT_VALUES));
86c4762a1bSJed Brown       row  = 4;
879566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,2,indices+4,values+4,INSERT_VALUES));
881a254869SHong Zhang     }
891a254869SHong Zhang   }
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m,&n));
959566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&min));
969566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(min,m,PETSC_DECIDE));
979566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(min));
989566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(min,&max));
999566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(min,&maxabs));
1009566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(min,&e));
101c4762a1bSJed Brown 
102f07e67edSHong Zhang   /* MatGetRowMax() */
1039566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMax\n"));
1049566063dSJacob Faibussowitsch   PetscCall(MatGetRowMax(A,max,NULL));
1059566063dSJacob Faibussowitsch   PetscCall(MatGetRowMax(A,max,imax));
1069566063dSJacob Faibussowitsch   PetscCall(VecView(max,PETSC_VIEWER_STDOUT_WORLD));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(max,&n));
1089566063dSJacob Faibussowitsch   PetscCall(PetscIntView(n,imax,PETSC_VIEWER_STDOUT_WORLD));
109fa213d2fSHong Zhang 
110f07e67edSHong Zhang   /* MatGetRowMin() */
1119566063dSJacob Faibussowitsch   PetscCall(MatScale(A,-1.0));
1129566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
1139566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMin\n"));
1149566063dSJacob Faibussowitsch   PetscCall(MatGetRowMin(A,min,NULL));
1159566063dSJacob Faibussowitsch   PetscCall(MatGetRowMin(A,min,imin));
116fa213d2fSHong Zhang 
1179566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(e,1.0,max,min)); /* e = max + min */
1189566063dSJacob Faibussowitsch   PetscCall(VecNorm(e,NORM_INFINITY,&enorm));
119e00437b9SBarry Smith   PetscCheck(enorm <= PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON ");
120f07e67edSHong Zhang   for (j = 0; j < n; j++) {
121e00437b9SBarry Smith     PetscCheck(imin[j] == imax[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%" PetscInt_FMT "] %" PetscInt_FMT " != imax %" PetscInt_FMT,j,imin[j],imax[j]);
122f07e67edSHong Zhang   }
123fa213d2fSHong Zhang 
124f07e67edSHong Zhang   /* MatGetRowMaxAbs() */
1259566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs\n"));
1269566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(A,maxabs,NULL));
1279566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(A,maxabs,imaxabs));
1289566063dSJacob Faibussowitsch   PetscCall(VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD));
1299566063dSJacob Faibussowitsch   PetscCall(PetscIntView(n,imaxabs,PETSC_VIEWER_STDOUT_WORLD));
1301a254869SHong Zhang 
131f07e67edSHong Zhang   /* MatGetRowMinAbs() */
1329566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMinAbs\n"));
1339566063dSJacob Faibussowitsch   PetscCall(MatGetRowMinAbs(A,min,NULL));
1349566063dSJacob Faibussowitsch   PetscCall(MatGetRowMinAbs(A,min,imin));
1359566063dSJacob Faibussowitsch   PetscCall(VecView(min,PETSC_VIEWER_STDOUT_WORLD));
1369566063dSJacob Faibussowitsch   PetscCall(PetscIntView(n,imin,PETSC_VIEWER_STDOUT_WORLD));
137f07e67edSHong Zhang 
138f07e67edSHong Zhang   if (size == 1) {
139f07e67edSHong Zhang     /* Test MatGetRowMax, MatGetRowMin and MatGetRowMaxAbs for SeqDense and MPIBAIJ matrix */
1404e879edeSHong Zhang     Mat Adense;
1414e879edeSHong Zhang     Vec max_d,maxabs_d;
1429566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(min,&max_d));
1439566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(min,&maxabs_d));
144f07e67edSHong Zhang 
1459566063dSJacob Faibussowitsch     PetscCall(MatScale(A,-1.0));
1469566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense));
1479566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMax for seqdense matrix\n"));
1489566063dSJacob Faibussowitsch     PetscCall(MatGetRowMax(Adense,max_d,imax));
1494e879edeSHong Zhang 
1509566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(e,-1.0,max,max_d)); /* e = -max + max_d */
1519566063dSJacob Faibussowitsch     PetscCall(VecNorm(e,NORM_INFINITY,&enorm));
152e00437b9SBarry Smith     PetscCheck(enorm <= PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-max + max_d) %g > PETSC_MACHINE_EPSILON",(double)enorm);
1534e879edeSHong Zhang 
1549566063dSJacob Faibussowitsch     PetscCall(MatScale(Adense,-1.0));
1559566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMin for seqdense matrix\n"));
1569566063dSJacob Faibussowitsch     PetscCall(MatGetRowMin(Adense,min,imin));
157f07e67edSHong Zhang 
1589566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(e,1.0,max,min)); /* e = max + min */
1599566063dSJacob Faibussowitsch     PetscCall(VecNorm(e,NORM_INFINITY,&enorm));
160e00437b9SBarry Smith     PetscCheck(enorm <= PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON ");
161f07e67edSHong Zhang     for (j = 0; j < n; j++) {
162f07e67edSHong Zhang       if (imin[j] != imax[j]) {
16398921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%" PetscInt_FMT "] %" PetscInt_FMT " != imax %" PetscInt_FMT " for seqdense matrix",j,imin[j],imax[j]);
164f07e67edSHong Zhang       }
165f07e67edSHong Zhang     }
166f07e67edSHong Zhang 
1679566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMaxAbs for seqdense matrix\n"));
1689566063dSJacob Faibussowitsch     PetscCall(MatGetRowMaxAbs(Adense,maxabs_d,imaxabs));
1699566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(e,-1.0,maxabs,maxabs_d)); /* e = -maxabs + maxabs_d */
1709566063dSJacob Faibussowitsch     PetscCall(VecNorm(e,NORM_INFINITY,&enorm));
171e00437b9SBarry Smith     PetscCheck(enorm <= PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm);
172c4762a1bSJed Brown 
1739566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Adense));
1749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&max_d));
1759566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&maxabs_d));
17643359b5eSHong Zhang   }
17743359b5eSHong Zhang 
17843359b5eSHong Zhang   { /* BAIJ matrix */
1794e879edeSHong Zhang     Mat               B;
18043359b5eSHong Zhang     Vec               maxabsB,maxabsB2;
18143359b5eSHong Zhang     PetscInt          bs=2,*imaxabsB,*imaxabsB2,rstart,rend,cstart,cend,ncols,col,Brows[2],Bcols[2];
18243359b5eSHong Zhang     const PetscInt    *cols;
18343359b5eSHong Zhang     const PetscScalar *vals,*vals2;
18443359b5eSHong Zhang     PetscScalar       Bvals[4];
18543359b5eSHong Zhang 
1869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(M,&imaxabsB,bs*M,&imaxabsB2));
18743359b5eSHong Zhang 
18843359b5eSHong Zhang     /* bs = 1 */
1899566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATMPIBAIJ,MAT_INITIAL_MATRIX,&B));
1909566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(min,&maxabsB));
1919566063dSJacob Faibussowitsch     PetscCall(MatGetRowMaxAbs(B,maxabsB,NULL));
1929566063dSJacob Faibussowitsch     PetscCall(MatGetRowMaxAbs(B,maxabsB,imaxabsB));
1939566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs for BAIJ matrix\n"));
1949566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(e,-1.0,maxabs,maxabsB)); /* e = -maxabs + maxabsB */
1959566063dSJacob Faibussowitsch     PetscCall(VecNorm(e,NORM_INFINITY,&enorm));
196e00437b9SBarry Smith     PetscCheck(enorm <= PETSC_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",(double)enorm);
1974e879edeSHong Zhang 
1984e879edeSHong Zhang     for (j = 0; j < n; j++) {
199e00437b9SBarry Smith       PetscCheck(imaxabs[j] == imaxabsB[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"imaxabs[%" PetscInt_FMT "] %" PetscInt_FMT " != imaxabsB %" PetscInt_FMT,j,imin[j],imax[j]);
200c4762a1bSJed Brown     }
2019566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
20243359b5eSHong Zhang 
20343359b5eSHong Zhang     /* Test bs = 2: Create B with bs*bs block structure of A */
2049566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD,&maxabsB2));
2059566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(maxabsB2,bs*m,PETSC_DECIDE));
2069566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(maxabsB2));
20743359b5eSHong Zhang 
2089566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
2099566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangeColumn(A,&cstart,&cend));
2109566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
2119566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,bs*(rend-rstart),bs*(cend-cstart),PETSC_DECIDE,PETSC_DECIDE));
2129566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(B));
2139566063dSJacob Faibussowitsch     PetscCall(MatSetUp(B));
21443359b5eSHong Zhang 
21543359b5eSHong Zhang     for (row=rstart; row<rend; row++) {
2169566063dSJacob Faibussowitsch       PetscCall(MatGetRow(A,row,&ncols,&cols,&vals));
21743359b5eSHong Zhang       for (col=0; col<ncols; col++) {
21843359b5eSHong Zhang         for (j=0; j<bs; j++) {
21943359b5eSHong Zhang           Brows[j] = bs*row + j;
22043359b5eSHong Zhang           Bcols[j] = bs*cols[col]+j;
22143359b5eSHong Zhang         }
22243359b5eSHong Zhang         for (j=0; j<bs*bs; j++) Bvals[j] = vals[col];
2239566063dSJacob Faibussowitsch         PetscCall(MatSetValues(B,bs,Brows,bs,Bcols,Bvals,INSERT_VALUES));
22443359b5eSHong Zhang       }
2259566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(A,row,&ncols,&cols,&vals));
22643359b5eSHong Zhang     }
2279566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2289566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
22943359b5eSHong Zhang 
2309566063dSJacob Faibussowitsch     PetscCall(MatGetRowMaxAbs(B,maxabsB2,imaxabsB2));
23143359b5eSHong Zhang 
23243359b5eSHong Zhang     /* Check maxabsB2 and imaxabsB2 */
2339566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(maxabsB,&vals));
2349566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(maxabsB2,&vals2));
23543359b5eSHong Zhang     for (row=0; row<m; row++) {
23643359b5eSHong Zhang       if (PetscAbsScalar(vals[row] - vals2[bs*row]) > PETSC_MACHINE_EPSILON)
23798921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"row %" PetscInt_FMT " maxabsB != maxabsB2",row);
23843359b5eSHong Zhang     }
2399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(maxabsB,&vals));
2409566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(maxabsB2,&vals2));
24143359b5eSHong Zhang 
24243359b5eSHong Zhang     for (col=0; col<n; col++) {
24343359b5eSHong Zhang       if (imaxabsB[col] != imaxabsB2[bs*col]/bs)
24498921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"col %" PetscInt_FMT " imaxabsB != imaxabsB2",col);
24543359b5eSHong Zhang     }
2469566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&maxabsB));
2479566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
2489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&maxabsB2));
2499566063dSJacob Faibussowitsch     PetscCall(PetscFree2(imaxabsB,imaxabsB2));
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown 
2529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&min));
2539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&max));
2549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&maxabs));
2559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&e));
2569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2579566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
258b122ec5aSJacob Faibussowitsch   return 0;
259c4762a1bSJed Brown }
260c4762a1bSJed Brown 
261c4762a1bSJed Brown /*TEST
262c4762a1bSJed Brown 
263c4762a1bSJed Brown    test:
264c4762a1bSJed Brown       output_file: output/ex114.out
265c4762a1bSJed Brown 
266c4762a1bSJed Brown    test:
267c4762a1bSJed Brown       suffix: 2
2684e879edeSHong Zhang       args: -testcase 1
2694e879edeSHong Zhang       output_file: output/ex114.out
270c4762a1bSJed Brown 
271c4762a1bSJed Brown    test:
272c4762a1bSJed Brown       suffix: 3
2734e879edeSHong Zhang       args: -testcase 2
2744e879edeSHong Zhang       output_file: output/ex114_3.out
2754e879edeSHong Zhang 
2764e879edeSHong Zhang    test:
2774e879edeSHong Zhang       suffix: 4
2784e879edeSHong Zhang       args: -testcase 3
2794e879edeSHong Zhang       output_file: output/ex114_4.out
2804e879edeSHong Zhang 
2814e879edeSHong Zhang    test:
2824e879edeSHong Zhang       suffix: 5
2834e879edeSHong Zhang       nsize: 3
2844e879edeSHong Zhang       args: -testcase 0
2854e879edeSHong Zhang       output_file: output/ex114_5.out
2864e879edeSHong Zhang 
2874e879edeSHong Zhang    test:
2884e879edeSHong Zhang       suffix: 6
2894e879edeSHong Zhang       nsize: 3
2904e879edeSHong Zhang       args: -testcase 1
2914e879edeSHong Zhang       output_file: output/ex114_6.out
2924e879edeSHong Zhang 
2934e879edeSHong Zhang    test:
2944e879edeSHong Zhang       suffix: 7
2954e879edeSHong Zhang       nsize: 3
2964e879edeSHong Zhang       args: -testcase 2
2974e879edeSHong Zhang       output_file: output/ex114_7.out
2984e879edeSHong Zhang 
2994e879edeSHong Zhang    test:
3004e879edeSHong Zhang       suffix: 8
3014e879edeSHong Zhang       nsize: 3
3024e879edeSHong Zhang       args: -testcase 3
3034e879edeSHong Zhang       output_file: output/ex114_8.out
304c4762a1bSJed Brown 
305c4762a1bSJed Brown TEST*/
306