xref: /petsc/src/mat/tests/ex48.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatSeqBAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,B,C,D,Fact;
9c4762a1bSJed Brown   Vec            xx,s1,s2,yy;
10c4762a1bSJed Brown   PetscInt       m=45,rows[2],cols[2],bs=1,i,row,col,*idx,M;
11c4762a1bSJed Brown   PetscScalar    rval,vals1[4],vals2[4];
12c4762a1bSJed Brown   PetscRandom    rdm;
13c4762a1bSJed Brown   IS             is1,is2;
14c4762a1bSJed Brown   PetscReal      s1norm,s2norm,rnorm,tol = 1.e-4;
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown   MatFactorInfo  info;
17c4762a1bSJed Brown 
18*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
19c4762a1bSJed Brown   /* Test MatSetValues() and MatGetValues() */
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
22c4762a1bSJed Brown   M    = m*bs;
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,M,&xx));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(xx,&s1));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(xx,&s2));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(xx,&yy));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* For each row add atleast 15 elements */
35c4762a1bSJed Brown   for (row=0; row<M; row++) {
36c4762a1bSJed Brown     for (i=0; i<25*bs; i++) {
375f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rdm,&rval));
38c4762a1bSJed Brown       col  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
395f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A,1,&row,1,&col,&rval,INSERT_VALUES));
405f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(B,1,&row,1,&col,&rval,INSERT_VALUES));
41c4762a1bSJed Brown     }
42c4762a1bSJed Brown   }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Now set blocks of values */
45c4762a1bSJed Brown   for (i=0; i<20*bs; i++) {
465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
47c4762a1bSJed Brown     cols[0]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
48c4762a1bSJed Brown     vals1[0] = rval;
495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
50c4762a1bSJed Brown     cols[1]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
51c4762a1bSJed Brown     vals1[1] = rval;
525f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
53c4762a1bSJed Brown     rows[0]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
54c4762a1bSJed Brown     vals1[2] = rval;
555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
56c4762a1bSJed Brown     rows[1]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
57c4762a1bSJed Brown     vals1[3] = rval;
585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,2,rows,2,cols,vals1,INSERT_VALUES));
595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,2,rows,2,cols,vals1,INSERT_VALUES));
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Test MatNorm() */
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A,NORM_FROBENIUS,&s1norm));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(B,NORM_FROBENIUS,&s2norm));
70c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm-s1norm)/s2norm;
71c4762a1bSJed Brown   if (rnorm>tol) {
725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
73c4762a1bSJed Brown   }
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A,NORM_INFINITY,&s1norm));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(B,NORM_INFINITY,&s2norm));
76c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm-s1norm)/s2norm;
77c4762a1bSJed Brown   if (rnorm>tol) {
785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
79c4762a1bSJed Brown   }
805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A,NORM_1,&s1norm));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(B,NORM_1,&s2norm));
82c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm-s1norm)/s2norm;
83c4762a1bSJed Brown   if (rnorm>tol) {
845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* MatShift() */
88c4762a1bSJed Brown   rval = 10*s1norm;
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,rval));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(B,rval));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Test MatTranspose() */
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INPLACE_MATRIX,&A));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(B,MAT_INPLACE_MATRIX,&B));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Now do MatGetValues()  */
97c4762a1bSJed Brown   for (i=0; i<30; i++) {
985f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
99c4762a1bSJed Brown     cols[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
101c4762a1bSJed Brown     cols[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
103c4762a1bSJed Brown     rows[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
1045f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
105c4762a1bSJed Brown     rows[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
1065f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetValues(A,2,rows,2,cols,vals1));
1075f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetValues(B,2,rows,2,cols,vals2));
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycmp(vals1,vals2,4,&flg));
109c4762a1bSJed Brown     if (!flg) {
1105f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetValues bs = %" PetscInt_FMT "\n",bs));
111c4762a1bSJed Brown     }
112c4762a1bSJed Brown   }
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Test MatMult(), MatMultAdd() */
115c4762a1bSJed Brown   for (i=0; i<40; i++) {
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(s2,0.0));
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,xx,s1));
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultAdd(A,xx,s2,s2));
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
122c4762a1bSJed Brown     rnorm = s2norm-s1norm;
123c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
1245f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatMult not equal to MatMultAdd Norm1=%e Norm2=%e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
125c4762a1bSJed Brown     }
126c4762a1bSJed Brown   }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* Test MatMult() */
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(A,B,10,&flg));
130c4762a1bSJed Brown   if (!flg) {
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult()\n"));
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Test MatMultAdd() */
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAddEqual(A,B,10,&flg));
136c4762a1bSJed Brown   if (!flg) {
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultAdd()\n"));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* Test MatMultTranspose() */
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeEqual(A,B,10,&flg));
142c4762a1bSJed Brown   if (!flg) {
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTranspose()\n"));
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeAddEqual(A,B,10,&flg));
148c4762a1bSJed Brown   if (!flg) {
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTransposeAdd()\n"));
150c4762a1bSJed Brown   }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Test MatMatMult() */
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&C));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(C,rdm));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,C,D,40,&flg));
157c4762a1bSJed Brown   if (!flg) {
1585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n"));
159c4762a1bSJed Brown   }
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&D));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&D));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,C,D,40,&flg));
164c4762a1bSJed Brown   if (!flg) {
1655f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n"));
166c4762a1bSJed Brown   }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* Do LUFactor() on both the matrices */
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(M,&idx));
170c4762a1bSJed Brown   for (i=0; i<M; i++) idx[i] = i;
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is1));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is2));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(is1));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(is2));
176c4762a1bSJed Brown 
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFactorInfoInitialize(&info));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   info.fill          = 2.0;
180c4762a1bSJed Brown   info.dtcol         = 0.0;
181c4762a1bSJed Brown   info.zeropivot     = 1.e-14;
182c4762a1bSJed Brown   info.pivotinblocks = 1.0;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   if (bs < 4) {
1855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(A,"petsc",MAT_FACTOR_LU,&Fact));
1865f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactorSymbolic(Fact,A,is1,is2,&info));
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactorNumeric(Fact,A,&info));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(yy,rdm));
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(MatForwardSolve(Fact,yy,xx));
1905f80ce2aSJacob Faibussowitsch     CHKERRQ(MatBackwardSolve(Fact,xx,s1));
1915f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Fact));
1925f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(s1,-1.0));
1935f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultAdd(A,s1,yy,yy));
1945f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(yy,NORM_2,&rnorm));
195c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
1965f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatForwardSolve/MatBackwardSolve - Norm1=%16.14e bs = %" PetscInt_FMT "\n",(double)rnorm,bs));
197c4762a1bSJed Brown     }
198c4762a1bSJed Brown   }
199c4762a1bSJed Brown 
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(B,is1,is2,&info));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(A,is1,is2,&info));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* Test MatSolveAdd() */
204c4762a1bSJed Brown   for (i=0; i<10; i++) {
2055f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
2065f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(yy,rdm));
2075f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(B,xx,yy,s2));
2085f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(A,xx,yy,s1));
2095f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
2105f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
211c4762a1bSJed Brown     rnorm = s2norm-s1norm;
212c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
2135f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
214c4762a1bSJed Brown     }
215c4762a1bSJed Brown   }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* Test MatSolveAdd() when x = A'b +x */
218c4762a1bSJed Brown   for (i=0; i<10; i++) {
2195f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
2205f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(s1,rdm));
2215f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(s2,s1));
2225f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(B,xx,s2,s2));
2235f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(A,xx,s1,s1));
2245f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
2255f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
226c4762a1bSJed Brown     rnorm = s2norm-s1norm;
227c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
2285f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd(same) - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
229c4762a1bSJed Brown     }
230c4762a1bSJed Brown   }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* Test MatSolve() */
233c4762a1bSJed Brown   for (i=0; i<10; i++) {
2345f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
2355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(B,xx,s2));
2365f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(A,xx,s1));
2375f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
2385f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
239c4762a1bSJed Brown     rnorm = s2norm-s1norm;
240c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
2415f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
242c4762a1bSJed Brown     }
243c4762a1bSJed Brown   }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /* Test MatSolveTranspose() */
246c4762a1bSJed Brown   if (bs < 8) {
247c4762a1bSJed Brown     for (i=0; i<10; i++) {
2485f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rdm));
2495f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSolveTranspose(B,xx,s2));
2505f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSolveTranspose(A,xx,s1));
2515f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
2525f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
253c4762a1bSJed Brown       rnorm = s2norm-s1norm;
254c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
2555f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveTranspose - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
256c4762a1bSJed Brown       }
257c4762a1bSJed Brown     }
258c4762a1bSJed Brown   }
259c4762a1bSJed Brown 
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
2625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&xx));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s1));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s2));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&yy));
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2));
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
271*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
272*b122ec5aSJacob Faibussowitsch   return 0;
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
275c4762a1bSJed Brown /*TEST
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    test:
278c4762a1bSJed Brown       args: -mat_block_size {{1 2 3 4 5 6 7 8}}
279c4762a1bSJed Brown 
280c4762a1bSJed Brown TEST*/
281