xref: /petsc/src/mat/tests/ex53.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPIBAIJ format.\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown #include <petscmat.h>
6*c4762a1bSJed Brown #define IMAX 15
7*c4762a1bSJed Brown int main(int argc,char **args)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   Mat               A,B,C,At,Bt;
10*c4762a1bSJed Brown   PetscViewer       fd;
11*c4762a1bSJed Brown   char              file[PETSC_MAX_PATH_LEN];
12*c4762a1bSJed Brown   PetscRandom       rand;
13*c4762a1bSJed Brown   Vec               xx,yy,s1,s2;
14*c4762a1bSJed Brown   PetscReal         s1norm,s2norm,rnorm,tol=1.e-10;
15*c4762a1bSJed Brown   PetscInt          rstart,rend,rows[2],cols[2],m,n,i,j,M,N,ct,row,ncols1,ncols2,bs;
16*c4762a1bSJed Brown   PetscMPIInt       rank,size;
17*c4762a1bSJed Brown   PetscErrorCode    ierr = 0;
18*c4762a1bSJed Brown   const PetscInt    *cols1,*cols2;
19*c4762a1bSJed Brown   PetscScalar       vals1[4],vals2[4],v;
20*c4762a1bSJed Brown   const PetscScalar *v1,*v2;
21*c4762a1bSJed Brown   PetscBool         flg;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
24*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown   /* Check out if MatLoad() works */
28*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
29*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Input file not specified");
30*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatSetType(A,MATBAIJ);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&xx);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = VecSetSizes(xx,m,PETSC_DECIDE);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = VecSetFromOptions(xx);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = VecDuplicate(xx,&s1);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = VecDuplicate(xx,&s2);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = VecDuplicate(xx,&yy);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown   /* Test MatNorm() */
50*c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_FROBENIUS,&s1norm);CHKERRQ(ierr);
51*c4762a1bSJed Brown   ierr  = MatNorm(B,NORM_FROBENIUS,&s2norm);CHKERRQ(ierr);
52*c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm;
53*c4762a1bSJed Brown   if (rnorm>tol) {
54*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
55*c4762a1bSJed Brown   }
56*c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_INFINITY,&s1norm);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr  = MatNorm(B,NORM_INFINITY,&s2norm);CHKERRQ(ierr);
58*c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm;
59*c4762a1bSJed Brown   if (rnorm>tol) {
60*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
61*c4762a1bSJed Brown   }
62*c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_1,&s1norm);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr  = MatNorm(B,NORM_1,&s2norm);CHKERRQ(ierr);
64*c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm;
65*c4762a1bSJed Brown   if (rnorm>tol) {
66*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
67*c4762a1bSJed Brown   }
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown   /* Test MatMult() */
70*c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
71*c4762a1bSJed Brown     ierr = VecSetRandom(xx,rand);CHKERRQ(ierr);
72*c4762a1bSJed Brown     ierr = MatMult(A,xx,s1);CHKERRQ(ierr);
73*c4762a1bSJed Brown     ierr = MatMult(B,xx,s2);CHKERRQ(ierr);
74*c4762a1bSJed Brown     ierr = VecAXPY(s2,-1.0,s1);CHKERRQ(ierr);
75*c4762a1bSJed Brown     ierr = VecNorm(s2,NORM_2,&rnorm);CHKERRQ(ierr);
76*c4762a1bSJed Brown     if (rnorm>tol) {
77*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMult - Norm2=%16.14e bs = %D\n",rank,rnorm,bs);CHKERRQ(ierr);
78*c4762a1bSJed Brown     }
79*c4762a1bSJed Brown   }
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   /* test MatMultAdd() */
82*c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
83*c4762a1bSJed Brown     ierr = VecSetRandom(xx,rand);CHKERRQ(ierr);
84*c4762a1bSJed Brown     ierr = VecSetRandom(yy,rand);CHKERRQ(ierr);
85*c4762a1bSJed Brown     ierr = MatMultAdd(A,xx,yy,s1);CHKERRQ(ierr);
86*c4762a1bSJed Brown     ierr = MatMultAdd(B,xx,yy,s2);CHKERRQ(ierr);
87*c4762a1bSJed Brown     ierr = VecAXPY(s2,-1.0,s1);CHKERRQ(ierr);
88*c4762a1bSJed Brown     ierr = VecNorm(s2,NORM_2,&rnorm);CHKERRQ(ierr);
89*c4762a1bSJed Brown     if (rnorm>tol) {
90*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultAdd - Norm2=%16.14e bs = %D\n",rank,rnorm,bs);CHKERRQ(ierr);
91*c4762a1bSJed Brown     }
92*c4762a1bSJed Brown   }
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown   /* Test MatMultTranspose() */
95*c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
96*c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
97*c4762a1bSJed Brown     ierr  = MatMultTranspose(A,xx,s1);CHKERRQ(ierr);
98*c4762a1bSJed Brown     ierr  = MatMultTranspose(B,xx,s2);CHKERRQ(ierr);
99*c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
100*c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
101*c4762a1bSJed Brown     rnorm = s2norm-s1norm;
102*c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
103*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultTranspose - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
104*c4762a1bSJed Brown     }
105*c4762a1bSJed Brown   }
106*c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
107*c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
108*c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
109*c4762a1bSJed Brown     ierr  = VecSetRandom(yy,rand);CHKERRQ(ierr);
110*c4762a1bSJed Brown     ierr  = MatMultTransposeAdd(A,xx,yy,s1);CHKERRQ(ierr);
111*c4762a1bSJed Brown     ierr  = MatMultTransposeAdd(B,xx,yy,s2);CHKERRQ(ierr);
112*c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
113*c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
114*c4762a1bSJed Brown     rnorm = s2norm-s1norm;
115*c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
116*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultTransposeAdd - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
117*c4762a1bSJed Brown     }
118*c4762a1bSJed Brown   }
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   /* Check MatGetValues() */
121*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
126*c4762a1bSJed Brown     /* Create random row numbers ad col numbers */
127*c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
128*c4762a1bSJed Brown     cols[0] = (int)(PetscRealPart(v)*N);
129*c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
130*c4762a1bSJed Brown     cols[1] = (int)(PetscRealPart(v)*N);
131*c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
132*c4762a1bSJed Brown     rows[0] = rstart + (int)(PetscRealPart(v)*m);
133*c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
134*c4762a1bSJed Brown     rows[1] = rstart + (int)(PetscRealPart(v)*m);
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown     ierr = MatGetValues(A,2,rows,2,cols,vals1);CHKERRQ(ierr);
137*c4762a1bSJed Brown     ierr = MatGetValues(B,2,rows,2,cols,vals2);CHKERRQ(ierr);
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown     for (j=0; j<4; j++) {
141*c4762a1bSJed Brown       if (vals1[j] != vals2[j]) {
142*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"[%d]: Error: MatGetValues rstart = %2d  row = %2d col = %2d val1 = %e val2 = %e bs = %D\n",rank,rstart,rows[j/2],cols[j%2],PetscRealPart(vals1[j]),PetscRealPart(vals2[j]),bs);CHKERRQ(ierr);
143*c4762a1bSJed Brown       }
144*c4762a1bSJed Brown     }
145*c4762a1bSJed Brown   }
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown   /* Test MatGetRow()/ MatRestoreRow() */
148*c4762a1bSJed Brown   for (ct=0; ct<100; ct++) {
149*c4762a1bSJed Brown     ierr = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
150*c4762a1bSJed Brown     row  = rstart + (int)(PetscRealPart(v)*m);
151*c4762a1bSJed Brown     ierr = MatGetRow(A,row,&ncols1,&cols1,&v1);CHKERRQ(ierr);
152*c4762a1bSJed Brown     ierr = MatGetRow(B,row,&ncols2,&cols2,&v2);CHKERRQ(ierr);
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown     for (i=0,j=0; i<ncols1 && j<ncols2; j++) {
155*c4762a1bSJed Brown       while (cols2[j] != cols1[i]) i++;
156*c4762a1bSJed Brown       if (v1[i] != v2[j]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect.");
157*c4762a1bSJed Brown     }
158*c4762a1bSJed Brown     if (j<ncols2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect");
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown     ierr = MatRestoreRow(A,row,&ncols1,&cols1,&v1);CHKERRQ(ierr);
161*c4762a1bSJed Brown     ierr = MatRestoreRow(B,row,&ncols2,&cols2,&v2);CHKERRQ(ierr);
162*c4762a1bSJed Brown   }
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown   /* Test MatConvert() */
165*c4762a1bSJed Brown   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr);
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown   /* See if MatMult Says both are same */
168*c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
169*c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
170*c4762a1bSJed Brown     ierr  = MatMult(A,xx,s1);CHKERRQ(ierr);
171*c4762a1bSJed Brown     ierr  = MatMult(C,xx,s2);CHKERRQ(ierr);
172*c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
173*c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
174*c4762a1bSJed Brown     rnorm = s2norm-s1norm;
175*c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
176*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error in MatConvert: MatMult - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
177*c4762a1bSJed Brown     }
178*c4762a1bSJed Brown   }
179*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown   /* Test MatTranspose() */
182*c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = MatTranspose(B,MAT_INITIAL_MATRIX,&Bt);CHKERRQ(ierr);
184*c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
185*c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
186*c4762a1bSJed Brown     ierr  = MatMult(At,xx,s1);CHKERRQ(ierr);
187*c4762a1bSJed Brown     ierr  = MatMult(Bt,xx,s2);CHKERRQ(ierr);
188*c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
189*c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
190*c4762a1bSJed Brown     rnorm = s2norm-s1norm;
191*c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
192*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error in MatConvert:MatMult - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
193*c4762a1bSJed Brown     }
194*c4762a1bSJed Brown   }
195*c4762a1bSJed Brown   ierr = MatDestroy(&At);CHKERRQ(ierr);
196*c4762a1bSJed Brown   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
199*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr = VecDestroy(&xx);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr = VecDestroy(&yy);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = VecDestroy(&s1);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = VecDestroy(&s2);CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
205*c4762a1bSJed Brown   ierr = PetscFinalize();
206*c4762a1bSJed Brown   return ierr;
207*c4762a1bSJed Brown }
208*c4762a1bSJed Brown 
209*c4762a1bSJed Brown 
210*c4762a1bSJed Brown /*TEST
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown    build:
213*c4762a1bSJed Brown       requires: !complex
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown    test:
216*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
217*c4762a1bSJed Brown       nsize: 3
218*c4762a1bSJed Brown       args: -matload_block_size 1 -f ${DATAFILESPATH}/matrices/small
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown    test:
221*c4762a1bSJed Brown       suffix: 2
222*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
223*c4762a1bSJed Brown       nsize: 3
224*c4762a1bSJed Brown       args: -matload_block_size 2 -f ${DATAFILESPATH}/matrices/small
225*c4762a1bSJed Brown       output_file: output/ex53_1.out
226*c4762a1bSJed Brown 
227*c4762a1bSJed Brown    test:
228*c4762a1bSJed Brown       suffix: 3
229*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
230*c4762a1bSJed Brown       nsize: 3
231*c4762a1bSJed Brown       args: -matload_block_size 4 -f ${DATAFILESPATH}/matrices/small
232*c4762a1bSJed Brown       output_file: output/ex53_1.out
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown    test:
235*c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
236*c4762a1bSJed Brown       suffix: 4
237*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
238*c4762a1bSJed Brown       nsize: 3
239*c4762a1bSJed Brown       args: -matload_block_size 5 -f ${DATAFILESPATH}/matrices/small
240*c4762a1bSJed Brown       output_file: output/ex53_1.out
241*c4762a1bSJed Brown 
242*c4762a1bSJed Brown    test:
243*c4762a1bSJed Brown       suffix: 5
244*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
245*c4762a1bSJed Brown       nsize: 3
246*c4762a1bSJed Brown       args: -matload_block_size 6 -f ${DATAFILESPATH}/matrices/small
247*c4762a1bSJed Brown       output_file: output/ex53_1.out
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown    test:
250*c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
251*c4762a1bSJed Brown       suffix: 6
252*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
253*c4762a1bSJed Brown       nsize: 3
254*c4762a1bSJed Brown       args: -matload_block_size 7 -f ${DATAFILESPATH}/matrices/small
255*c4762a1bSJed Brown       output_file: output/ex53_1.out
256*c4762a1bSJed Brown 
257*c4762a1bSJed Brown    test:
258*c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
259*c4762a1bSJed Brown       suffix: 7
260*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
261*c4762a1bSJed Brown       nsize: 3
262*c4762a1bSJed Brown       args: -matload_block_size 8 -f ${DATAFILESPATH}/matrices/small
263*c4762a1bSJed Brown       output_file: output/ex53_1.out
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown    test:
266*c4762a1bSJed Brown       suffix: 8
267*c4762a1bSJed Brown       requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
268*c4762a1bSJed Brown       nsize: 4
269*c4762a1bSJed Brown       args: -matload_block_size 3 -f ${DATAFILESPATH}/matrices/small
270*c4762a1bSJed Brown       output_file: output/ex53_1.out
271*c4762a1bSJed Brown 
272*c4762a1bSJed Brown TEST*/
273