xref: /petsc/src/mat/tests/ex75.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPISBAIJ format.\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,u,s1,s2;
9c4762a1bSJed Brown   Mat            A,sA,sB;
10c4762a1bSJed Brown   PetscRandom    rctx;
11c4762a1bSJed Brown   PetscReal      r1,r2,rnorm,tol = PETSC_SQRT_MACHINE_EPSILON;
12c4762a1bSJed Brown   PetscScalar    one=1.0, neg_one=-1.0, value[3], four=4.0,alpha=0.1;
1343359b5eSHong Zhang   PetscInt       n,col[3],n1,block,row,i,j,i2,j2,Ii,J,rstart,rend,bs=1,mbs=16,d_nz=3,o_nz=3,prob=1;
14c4762a1bSJed Brown   PetscMPIInt    size,rank;
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown 
17*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
18*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL));
19*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
20*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-prob",&prob,NULL));
21c4762a1bSJed Brown 
22*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
24c4762a1bSJed Brown 
2543359b5eSHong Zhang   /* Create a BAIJ matrix A */
26c4762a1bSJed Brown   n = mbs*bs;
27*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
28*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
29*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATBAIJ));
30*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
31*9566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A,bs,d_nz,NULL,o_nz,NULL));
32*9566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A,bs,d_nz,NULL));
33*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   if (bs == 1) {
36c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
3743359b5eSHong Zhang       value[0] = -1.0; value[1] = 0.0; value[2] = -1.0;
38c4762a1bSJed Brown       for (i=1; i<n-1; i++) {
39c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
40*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
41c4762a1bSJed Brown       }
42c4762a1bSJed Brown       i       = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
4343359b5eSHong Zhang       value[0]= 0.1; value[1]=-1.0; value[2]=0.0;
44*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown       i        = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
4743359b5eSHong Zhang       value[0] = 0.0; value[1] = -1.0; value[2]=0.1;
48*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
49c4762a1bSJed Brown     } else if (prob ==2) { /* matrix for the five point stencil */
50c4762a1bSJed Brown       n1 = (int) PetscSqrtReal((PetscReal)n);
51c4762a1bSJed Brown       for (i=0; i<n1; i++) {
52c4762a1bSJed Brown         for (j=0; j<n1; j++) {
53c4762a1bSJed Brown           Ii = j + n1*i;
54*9566063dSJacob Faibussowitsch           if (i>0)    {J = Ii - n1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
55*9566063dSJacob Faibussowitsch           if (i<n1-1) {J = Ii + n1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
56*9566063dSJacob Faibussowitsch           if (j>0)    {J = Ii - 1;  PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
57*9566063dSJacob Faibussowitsch           if (j<n1-1) {J = Ii + 1;  PetscCall(MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES));}
58*9566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES));
59c4762a1bSJed Brown         }
60c4762a1bSJed Brown       }
61c4762a1bSJed Brown     }
62c4762a1bSJed Brown     /* end of if (bs == 1) */
63c4762a1bSJed Brown   } else {  /* bs > 1 */
64c4762a1bSJed Brown     for (block=0; block<n/bs; block++) {
65c4762a1bSJed Brown       /* diagonal blocks */
66c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
67c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
68c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
69*9566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES));
70c4762a1bSJed Brown       }
71c4762a1bSJed Brown       i       = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
72c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
73*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown       i       = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
76c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
77*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
78c4762a1bSJed Brown     }
79c4762a1bSJed Brown     /* off-diagonal blocks */
80c4762a1bSJed Brown     value[0]=-1.0;
81c4762a1bSJed Brown     for (i=0; i<(n/bs-1)*bs; i++) {
82c4762a1bSJed Brown       col[0]=i+bs;
83*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&i,1,col,value,INSERT_VALUES));
84c4762a1bSJed Brown       col[0]=i; row=i+bs;
85*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,1,col,value,INSERT_VALUES));
86c4762a1bSJed Brown     }
87c4762a1bSJed Brown   }
88*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
89*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
90*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
9143359b5eSHong Zhang 
9243359b5eSHong Zhang   /* Get SBAIJ matrix sA from A */
93*9566063dSJacob Faibussowitsch   PetscCall(MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Test MatGetSize(), MatGetLocalSize() */
96*9566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sA, &i,&j));
97*9566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &i2,&j2));
98c4762a1bSJed Brown   i   -= i2; j -= j2;
99c4762a1bSJed Brown   if (i || j) {
100*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetSize()\n",rank));
101*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown 
104*9566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(sA, &i,&j));
105*9566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &i2,&j2));
106c4762a1bSJed Brown   i2  -= i; j2 -= j;
107c4762a1bSJed Brown   if (i2 || j2) {
108*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetLocalSize()\n",rank));
109*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* vectors */
113c4762a1bSJed Brown   /*--------------------*/
114c4762a1bSJed Brown   /* i is obtained from MatGetLocalSize() */
115*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
116*9566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,i,PETSC_DECIDE));
117*9566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
118*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&y));
119*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&u));
120*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&s1));
121*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&s2));
122c4762a1bSJed Brown 
123*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
124*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
125*9566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x,rctx));
126*9566063dSJacob Faibussowitsch   PetscCall(VecSet(u,one));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* Test MatNorm() */
129*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_FROBENIUS,&r1));
130*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(sA,NORM_FROBENIUS,&r2));
131c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
132dd400576SPatrick Sanan   if (rnorm > tol && rank == 0) {
133*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs));
134c4762a1bSJed Brown   }
135*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_INFINITY,&r1));
136*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(sA,NORM_INFINITY,&r2));
137c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
138dd400576SPatrick Sanan   if (rnorm > tol && rank == 0) {
139*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_INFINITY(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs));
140c4762a1bSJed Brown   }
141*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_1,&r1));
142*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(sA,NORM_1,&r2));
143c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
144dd400576SPatrick Sanan   if (rnorm > tol && rank == 0) {
145*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_1(), Anorm=%16.14e, sAnorm=%16.14e bs=%" PetscInt_FMT "\n",(double)r1,(double)r2,bs));
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
149*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(sA,&rstart,&rend));
150*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&i2,&j2));
151c4762a1bSJed Brown   i2  -= rstart; j2 -= rend;
152c4762a1bSJed Brown   if (i2 || j2) {
153*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MaGetOwnershipRange()\n",rank));
154*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* Test MatDiagonalScale() */
158*9566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A,x,x));
159*9566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(sA,x,x));
160*9566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A,sA,10,&flg));
16128b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale");
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Test MatGetDiagonal(), MatScale() */
164*9566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A,s1));
165*9566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(sA,s2));
166*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(s1,NORM_1,&r1));
167*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2,NORM_1,&r2));
168c4762a1bSJed Brown   r1  -= r2;
169c4762a1bSJed Brown   if (r1<-tol || r1>tol) {
170*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n",rank,(double)r1));
171*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown 
174*9566063dSJacob Faibussowitsch   PetscCall(MatScale(A,alpha));
175*9566063dSJacob Faibussowitsch   PetscCall(MatScale(sA,alpha));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* Test MatGetRowMaxAbs() */
178*9566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(A,s1,NULL));
179*9566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(sA,s2,NULL));
180c4762a1bSJed Brown 
181*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(s1,NORM_1,&r1));
182*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2,NORM_1,&r2));
183c4762a1bSJed Brown   r1  -= r2;
184c4762a1bSJed Brown   if (r1<-tol || r1>tol) {
185*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRowMaxAbs() \n"));
186c4762a1bSJed Brown   }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* Test MatMult(), MatMultAdd() */
189*9566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A,sA,10,&flg));
190c4762a1bSJed Brown   if (!flg) {
191*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale()\n",rank));
192*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
193c4762a1bSJed Brown   }
194c4762a1bSJed Brown 
195*9566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(A,sA,10,&flg));
196c4762a1bSJed Brown   if (!flg) {
197*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd()\n",rank));
198*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* Test MatMultTranspose(), MatMultTransposeAdd() */
202c4762a1bSJed Brown   for (i=0; i<10; i++) {
203*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x,rctx));
204*9566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,x,s1));
205*9566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(sA,x,s2));
206*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1,NORM_1,&r1));
207*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2,NORM_1,&r2));
208c4762a1bSJed Brown     r1  -= r2;
209c4762a1bSJed Brown     if (r1<-tol || r1>tol) {
210*9566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%g\n",rank,(double)r1));
211*9566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
212c4762a1bSJed Brown     }
213c4762a1bSJed Brown   }
214c4762a1bSJed Brown   for (i=0; i<10; i++) {
215*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x,rctx));
216*9566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y,rctx));
217*9566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(A,x,y,s1));
218*9566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(sA,x,y,s2));
219*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1,NORM_1,&r1));
220*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2,NORM_1,&r2));
221c4762a1bSJed Brown     r1  -= r2;
222c4762a1bSJed Brown     if (r1<-tol || r1>tol) {
223*9566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g \n",rank,(double)r1));
224*9566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
225c4762a1bSJed Brown     }
226c4762a1bSJed Brown   }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* Test MatDuplicate() */
229*9566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(sA,MAT_COPY_VALUES,&sB));
230*9566063dSJacob Faibussowitsch   PetscCall(MatEqual(sA,sB,&flg));
231c4762a1bSJed Brown   if (!flg) {
232*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error in MatDuplicate(), sA != sB \n"));
233c4762a1bSJed Brown   }
234*9566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(sA,sB,5,&flg));
235c4762a1bSJed Brown   if (!flg) {
236*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMult()\n",rank));
237*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
238c4762a1bSJed Brown   }
239*9566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(sA,sB,5,&flg));
240c4762a1bSJed Brown   if (!flg) {
241*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMultAdd(()\n",rank));
242*9566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
243c4762a1bSJed Brown   }
244*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sB));
245*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
246*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
247*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
248*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
249*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
250*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
251*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
252*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
253*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
254b122ec5aSJacob Faibussowitsch   return 0;
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown /*TEST
258c4762a1bSJed Brown 
259c4762a1bSJed Brown    test:
260c4762a1bSJed Brown       nsize: {{1 3}}
26143359b5eSHong Zhang       args: -bs {{1 2 3  5  7 8}} -mat_ignore_lower_triangular -prob {{1 2}}
262c4762a1bSJed Brown 
263c4762a1bSJed Brown TEST*/
264