xref: /petsc/src/mat/tests/ex48.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscInt       m=45,rows[2],cols[2],bs=1,i,row,col,*idx,M;
12c4762a1bSJed Brown   PetscScalar    rval,vals1[4],vals2[4];
13c4762a1bSJed Brown   PetscRandom    rdm;
14c4762a1bSJed Brown   IS             is1,is2;
15c4762a1bSJed Brown   PetscReal      s1norm,s2norm,rnorm,tol = 1.e-4;
16c4762a1bSJed Brown   PetscBool      flg;
17c4762a1bSJed Brown   MatFactorInfo  info;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20c4762a1bSJed Brown   /* Test MatSetValues() and MatGetValues() */
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
23c4762a1bSJed Brown   M    = m*bs;
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,M,&xx));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(xx,&s1));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(xx,&s2));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(xx,&yy));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /* For each row add atleast 15 elements */
36c4762a1bSJed Brown   for (row=0; row<M; row++) {
37c4762a1bSJed Brown     for (i=0; i<25*bs; i++) {
38*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rdm,&rval));
39c4762a1bSJed Brown       col  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
40*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A,1,&row,1,&col,&rval,INSERT_VALUES));
41*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(B,1,&row,1,&col,&rval,INSERT_VALUES));
42c4762a1bSJed Brown     }
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Now set blocks of values */
46c4762a1bSJed Brown   for (i=0; i<20*bs; i++) {
47*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
48c4762a1bSJed Brown     cols[0]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
49c4762a1bSJed Brown     vals1[0] = rval;
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
51c4762a1bSJed Brown     cols[1]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
52c4762a1bSJed Brown     vals1[1] = rval;
53*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
54c4762a1bSJed Brown     rows[0]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
55c4762a1bSJed Brown     vals1[2] = rval;
56*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
57c4762a1bSJed Brown     rows[1]  = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
58c4762a1bSJed Brown     vals1[3] = rval;
59*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,2,rows,2,cols,vals1,INSERT_VALUES));
60*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(B,2,rows,2,cols,vals1,INSERT_VALUES));
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown 
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Test MatNorm() */
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A,NORM_FROBENIUS,&s1norm));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(B,NORM_FROBENIUS,&s2norm));
71c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm-s1norm)/s2norm;
72c4762a1bSJed Brown   if (rnorm>tol) {
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
74c4762a1bSJed Brown   }
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A,NORM_INFINITY,&s1norm));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(B,NORM_INFINITY,&s2norm));
77c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm-s1norm)/s2norm;
78c4762a1bSJed Brown   if (rnorm>tol) {
79*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
80c4762a1bSJed Brown   }
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A,NORM_1,&s1norm));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(B,NORM_1,&s2norm));
83c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm-s1norm)/s2norm;
84c4762a1bSJed Brown   if (rnorm>tol) {
85*5f80ce2aSJacob 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));
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* MatShift() */
89c4762a1bSJed Brown   rval = 10*s1norm;
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,rval));
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(B,rval));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* Test MatTranspose() */
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INPLACE_MATRIX,&A));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(B,MAT_INPLACE_MATRIX,&B));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* Now do MatGetValues()  */
98c4762a1bSJed Brown   for (i=0; i<30; i++) {
99*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
100c4762a1bSJed Brown     cols[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
101*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
102c4762a1bSJed Brown     cols[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
103*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
104c4762a1bSJed Brown     rows[0] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rdm,&rval));
106c4762a1bSJed Brown     rows[1] = PetscMin(M-1,(PetscInt)(PetscRealPart(rval)*M));
107*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetValues(A,2,rows,2,cols,vals1));
108*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetValues(B,2,rows,2,cols,vals2));
109*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycmp(vals1,vals2,4,&flg));
110c4762a1bSJed Brown     if (!flg) {
111*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatGetValues bs = %" PetscInt_FMT "\n",bs));
112c4762a1bSJed Brown     }
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* Test MatMult(), MatMultAdd() */
116c4762a1bSJed Brown   for (i=0; i<40; i++) {
117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
118*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(s2,0.0));
119*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,xx,s1));
120*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultAdd(A,xx,s2,s2));
121*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
122*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
123c4762a1bSJed Brown     rnorm = s2norm-s1norm;
124c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
125*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatMult not equal to MatMultAdd Norm1=%e Norm2=%e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
126c4762a1bSJed Brown     }
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Test MatMult() */
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(A,B,10,&flg));
131c4762a1bSJed Brown   if (!flg) {
132*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMult()\n"));
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* Test MatMultAdd() */
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAddEqual(A,B,10,&flg));
137c4762a1bSJed Brown   if (!flg) {
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultAdd()\n"));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Test MatMultTranspose() */
142*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeEqual(A,B,10,&flg));
143c4762a1bSJed Brown   if (!flg) {
144*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTranspose()\n"));
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTransposeAddEqual(A,B,10,&flg));
149c4762a1bSJed Brown   if (!flg) {
150*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMultTransposeAdd()\n"));
151c4762a1bSJed Brown   }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* Test MatMatMult() */
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&C));
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(C,rdm));
156*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D));
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,C,D,40,&flg));
158c4762a1bSJed Brown   if (!flg) {
159*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n"));
160c4762a1bSJed Brown   }
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,M,40,NULL,&D));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&D));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,C,D,40,&flg));
165c4762a1bSJed Brown   if (!flg) {
166*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult()\n"));
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* Do LUFactor() on both the matrices */
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(M,&idx));
171c4762a1bSJed Brown   for (i=0; i<M; i++) idx[i] = i;
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is1));
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,M,idx,PETSC_COPY_VALUES,&is2));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(is1));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(is2));
177c4762a1bSJed Brown 
178*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFactorInfoInitialize(&info));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   info.fill          = 2.0;
181c4762a1bSJed Brown   info.dtcol         = 0.0;
182c4762a1bSJed Brown   info.zeropivot     = 1.e-14;
183c4762a1bSJed Brown   info.pivotinblocks = 1.0;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   if (bs < 4) {
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(A,"petsc",MAT_FACTOR_LU,&Fact));
187*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactorSymbolic(Fact,A,is1,is2,&info));
188*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactorNumeric(Fact,A,&info));
189*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(yy,rdm));
190*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatForwardSolve(Fact,yy,xx));
191*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatBackwardSolve(Fact,xx,s1));
192*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Fact));
193*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(s1,-1.0));
194*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultAdd(A,s1,yy,yy));
195*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(yy,NORM_2,&rnorm));
196c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
197*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatForwardSolve/MatBackwardSolve - Norm1=%16.14e bs = %" PetscInt_FMT "\n",(double)rnorm,bs));
198c4762a1bSJed Brown     }
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown 
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(B,is1,is2,&info));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(A,is1,is2,&info));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* Test MatSolveAdd() */
205c4762a1bSJed Brown   for (i=0; i<10; i++) {
206*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
207*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(yy,rdm));
208*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(B,xx,yy,s2));
209*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(A,xx,yy,s1));
210*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
211*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
212c4762a1bSJed Brown     rnorm = s2norm-s1norm;
213c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
214*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
215c4762a1bSJed Brown     }
216c4762a1bSJed Brown   }
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* Test MatSolveAdd() when x = A'b +x */
219c4762a1bSJed Brown   for (i=0; i<10; i++) {
220*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
221*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(s1,rdm));
222*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(s2,s1));
223*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(B,xx,s2,s2));
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolveAdd(A,xx,s1,s1));
225*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
226*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
227c4762a1bSJed Brown     rnorm = s2norm-s1norm;
228c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
229*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveAdd(same) - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
230c4762a1bSJed Brown     }
231c4762a1bSJed Brown   }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* Test MatSolve() */
234c4762a1bSJed Brown   for (i=0; i<10; i++) {
235*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(xx,rdm));
236*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(B,xx,s2));
237*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(A,xx,s1));
238*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
239*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
240c4762a1bSJed Brown     rnorm = s2norm-s1norm;
241c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
242*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolve - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
243c4762a1bSJed Brown     }
244c4762a1bSJed Brown   }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /* Test MatSolveTranspose() */
247c4762a1bSJed Brown   if (bs < 8) {
248c4762a1bSJed Brown     for (i=0; i<10; i++) {
249*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(xx,rdm));
250*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSolveTranspose(B,xx,s2));
251*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSolveTranspose(A,xx,s1));
252*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_2,&s1norm));
253*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_2,&s2norm));
254c4762a1bSJed Brown       rnorm = s2norm-s1norm;
255c4762a1bSJed Brown       if (rnorm<-tol || rnorm>tol) {
256*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error:MatSolveTranspose - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n",(double)s1norm,(double)s2norm,bs));
257c4762a1bSJed Brown       }
258c4762a1bSJed Brown     }
259c4762a1bSJed Brown   }
260c4762a1bSJed Brown 
261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
262*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
263*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
264*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&xx));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s1));
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s2));
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&yy));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
270*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2));
271*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
272c4762a1bSJed Brown   ierr = PetscFinalize();
273c4762a1bSJed Brown   return ierr;
274c4762a1bSJed Brown }
275c4762a1bSJed Brown 
276c4762a1bSJed Brown /*TEST
277c4762a1bSJed Brown 
278c4762a1bSJed Brown    test:
279c4762a1bSJed Brown       args: -mat_block_size {{1 2 3 4 5 6 7 8}}
280c4762a1bSJed Brown 
281c4762a1bSJed Brown TEST*/
282