xref: /petsc/src/mat/tests/ex75.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPISBAIJ format.\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Vec            x,y,u,s1,s2;
9*c4762a1bSJed Brown   Mat            A,sA,sB;
10*c4762a1bSJed Brown   PetscRandom    rctx;
11*c4762a1bSJed Brown   PetscReal      r1,r2,rnorm,tol = PETSC_SQRT_MACHINE_EPSILON;
12*c4762a1bSJed Brown   PetscScalar    one=1.0, neg_one=-1.0, value[3], four=4.0,alpha=0.1;
13*c4762a1bSJed Brown   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=2;
14*c4762a1bSJed Brown   PetscErrorCode ierr;
15*c4762a1bSJed Brown   PetscMPIInt    size,rank;
16*c4762a1bSJed Brown   PetscBool      flg;
17*c4762a1bSJed Brown   MatType        type;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr);
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown   n = mbs*bs;
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   /* Assemble MPISBAIJ matrix sA */
29*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&sA);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatSetSizes(sA,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatSetType(sA,MATSBAIJ);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatSetFromOptions(sA);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = MatGetType(sA,&type);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = MatMPISBAIJSetPreallocation(sA,bs,d_nz,NULL,o_nz,NULL);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MatSeqSBAIJSetPreallocation(sA,bs,d_nz,NULL);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = MatSetOption(sA,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   if (bs == 1) {
39*c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
40*c4762a1bSJed Brown       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
41*c4762a1bSJed Brown       for (i=1; i<n-1; i++) {
42*c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
43*c4762a1bSJed Brown         ierr   = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
44*c4762a1bSJed Brown       }
45*c4762a1bSJed Brown       i       = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
46*c4762a1bSJed Brown       value[0]= 0.1; value[1]=-1; value[2]=2;
47*c4762a1bSJed Brown       ierr    = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown       i        = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
50*c4762a1bSJed Brown       value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
51*c4762a1bSJed Brown       ierr     = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
52*c4762a1bSJed Brown     } else if (prob ==2) { /* matrix for the five point stencil */
53*c4762a1bSJed Brown       n1 =  (int) PetscSqrtReal((PetscReal)n);
54*c4762a1bSJed Brown       if (n1*n1 != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"n must be a perfect square of n1");
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown       for (i=0; i<n1; i++) {
57*c4762a1bSJed Brown         for (j=0; j<n1; j++) {
58*c4762a1bSJed Brown           Ii = j + n1*i;
59*c4762a1bSJed Brown           if (i>0)    {J = Ii - n1; ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
60*c4762a1bSJed Brown           if (i<n1-1) {J = Ii + n1; ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
61*c4762a1bSJed Brown           if (j>0)    {J = Ii - 1;  ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
62*c4762a1bSJed Brown           if (j<n1-1) {J = Ii + 1;  ierr = MatSetValues(sA,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
63*c4762a1bSJed Brown           ierr = MatSetValues(sA,1,&Ii,1,&Ii,&four,INSERT_VALUES);CHKERRQ(ierr);
64*c4762a1bSJed Brown         }
65*c4762a1bSJed Brown       }
66*c4762a1bSJed Brown     }
67*c4762a1bSJed Brown     /* end of if (bs == 1) */
68*c4762a1bSJed Brown   } else {  /* bs > 1 */
69*c4762a1bSJed Brown     for (block=0; block<n/bs; block++) {
70*c4762a1bSJed Brown       /* diagonal blocks */
71*c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
72*c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
73*c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
74*c4762a1bSJed Brown         ierr   = MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
75*c4762a1bSJed Brown       }
76*c4762a1bSJed Brown       i       = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
77*c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
78*c4762a1bSJed Brown       ierr    = MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
79*c4762a1bSJed Brown 
80*c4762a1bSJed Brown       i       = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
81*c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
82*c4762a1bSJed Brown       ierr    = MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
83*c4762a1bSJed Brown     }
84*c4762a1bSJed Brown     /* off-diagonal blocks */
85*c4762a1bSJed Brown     value[0]=-1.0;
86*c4762a1bSJed Brown     for (i=0; i<(n/bs-1)*bs; i++) {
87*c4762a1bSJed Brown       col[0]=i+bs;
88*c4762a1bSJed Brown       ierr  = MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr);
89*c4762a1bSJed Brown       col[0]=i; row=i+bs;
90*c4762a1bSJed Brown       ierr  = MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr);
91*c4762a1bSJed Brown     }
92*c4762a1bSJed Brown   }
93*c4762a1bSJed Brown   ierr = MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   /* Test MatView() */
97*c4762a1bSJed Brown   ierr = MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,n,n,d_nz,NULL,o_nz,NULL,&A);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   if (bs == 1) {
101*c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
102*c4762a1bSJed Brown       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
103*c4762a1bSJed Brown       for (i=1; i<n-1; i++) {
104*c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
105*c4762a1bSJed Brown         ierr   = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
106*c4762a1bSJed Brown       }
107*c4762a1bSJed Brown       i       = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
108*c4762a1bSJed Brown       value[0]= 0.1; value[1]=-1; value[2]=2;
109*c4762a1bSJed Brown       ierr    = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown       i        = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
112*c4762a1bSJed Brown       value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
113*c4762a1bSJed Brown       ierr     = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
114*c4762a1bSJed Brown     } else if (prob ==2) { /* matrix for the five point stencil */
115*c4762a1bSJed Brown       n1 = (int) PetscSqrtReal((PetscReal)n);
116*c4762a1bSJed Brown       for (i=0; i<n1; i++) {
117*c4762a1bSJed Brown         for (j=0; j<n1; j++) {
118*c4762a1bSJed Brown           Ii = j + n1*i;
119*c4762a1bSJed Brown           if (i>0)    {J = Ii - n1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
120*c4762a1bSJed Brown           if (i<n1-1) {J = Ii + n1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
121*c4762a1bSJed Brown           if (j>0)    {J = Ii - 1;  ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
122*c4762a1bSJed Brown           if (j<n1-1) {J = Ii + 1;  ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr);}
123*c4762a1bSJed Brown           ierr = MatSetValues(A,1,&Ii,1,&Ii,&four,INSERT_VALUES);CHKERRQ(ierr);
124*c4762a1bSJed Brown         }
125*c4762a1bSJed Brown       }
126*c4762a1bSJed Brown     }
127*c4762a1bSJed Brown     /* end of if (bs == 1) */
128*c4762a1bSJed Brown   } else {  /* bs > 1 */
129*c4762a1bSJed Brown     for (block=0; block<n/bs; block++) {
130*c4762a1bSJed Brown       /* diagonal blocks */
131*c4762a1bSJed Brown       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
132*c4762a1bSJed Brown       for (i=1+block*bs; i<bs-1+block*bs; i++) {
133*c4762a1bSJed Brown         col[0] = i-1; col[1] = i; col[2] = i+1;
134*c4762a1bSJed Brown         ierr   = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
135*c4762a1bSJed Brown       }
136*c4762a1bSJed Brown       i       = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
137*c4762a1bSJed Brown       value[0]=-1.0; value[1]=4.0;
138*c4762a1bSJed Brown       ierr    = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown       i       = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
141*c4762a1bSJed Brown       value[0]=4.0; value[1] = -1.0;
142*c4762a1bSJed Brown       ierr    = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
143*c4762a1bSJed Brown     }
144*c4762a1bSJed Brown     /* off-diagonal blocks */
145*c4762a1bSJed Brown     value[0]=-1.0;
146*c4762a1bSJed Brown     for (i=0; i<(n/bs-1)*bs; i++) {
147*c4762a1bSJed Brown       col[0]=i+bs;
148*c4762a1bSJed Brown       ierr  = MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);CHKERRQ(ierr);
149*c4762a1bSJed Brown       col[0]=i; row=i+bs;
150*c4762a1bSJed Brown       ierr  = MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);CHKERRQ(ierr);
151*c4762a1bSJed Brown     }
152*c4762a1bSJed Brown   }
153*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown   /* Test MatGetSize(), MatGetLocalSize() */
157*c4762a1bSJed Brown   ierr = MatGetSize(sA, &i,&j); CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = MatGetSize(A, &i2,&j2);CHKERRQ(ierr);
159*c4762a1bSJed Brown   i   -= i2; j -= j2;
160*c4762a1bSJed Brown   if (i || j) {
161*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetSize()\n",rank);CHKERRQ(ierr);
162*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
163*c4762a1bSJed Brown   }
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown   ierr = MatGetLocalSize(sA, &i,&j);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = MatGetLocalSize(A, &i2,&j2);CHKERRQ(ierr);
167*c4762a1bSJed Brown   i2  -= i; j2 -= j;
168*c4762a1bSJed Brown   if (i2 || j2) {
169*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetLocalSize()\n",rank);CHKERRQ(ierr);
170*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
171*c4762a1bSJed Brown   }
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown   /* vectors */
174*c4762a1bSJed Brown   /*--------------------*/
175*c4762a1bSJed Brown   /* i is obtained from MatGetLocalSize() */
176*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = VecSetSizes(x,i,PETSC_DECIDE);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
180*c4762a1bSJed Brown   ierr = VecDuplicate(x,&u);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ierr = VecDuplicate(x,&s1);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = VecDuplicate(x,&s2);CHKERRQ(ierr);
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = VecSet(u,one);CHKERRQ(ierr);
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown   /* Test MatNorm() */
190*c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_FROBENIUS,&r1);CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr  = MatNorm(sA,NORM_FROBENIUS,&r2);CHKERRQ(ierr);
192*c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
193*c4762a1bSJed Brown   if (rnorm > tol && !rank) {
194*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm_FROBENIUS(), Anorm=%16.14e, sAnorm=%16.14e bs=%D\n",r1,r2,bs);CHKERRQ(ierr);
195*c4762a1bSJed Brown   }
196*c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_INFINITY,&r1);CHKERRQ(ierr);
197*c4762a1bSJed Brown   ierr  = MatNorm(sA,NORM_INFINITY,&r2);CHKERRQ(ierr);
198*c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
199*c4762a1bSJed Brown   if (rnorm > tol && !rank) {
200*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_INFINITY(), Anorm=%16.14e, sAnorm=%16.14e bs=%D\n",r1,r2,bs);CHKERRQ(ierr);
201*c4762a1bSJed Brown   }
202*c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_1,&r1);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr  = MatNorm(sA,NORM_1,&r2);CHKERRQ(ierr);
204*c4762a1bSJed Brown   rnorm = PetscAbsReal(r1-r2)/r2;
205*c4762a1bSJed Brown   if (rnorm > tol && !rank) {
206*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: MatNorm_1(), Anorm=%16.14e, sAnorm=%16.14e bs=%D\n",r1,r2,bs);CHKERRQ(ierr);
207*c4762a1bSJed Brown   }
208*c4762a1bSJed Brown 
209*c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
210*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(sA,&rstart,&rend);CHKERRQ(ierr);
211*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&i2,&j2);CHKERRQ(ierr);
212*c4762a1bSJed Brown   i2  -= rstart; j2 -= rend;
213*c4762a1bSJed Brown   if (i2 || j2) {
214*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MaGetOwnershipRange()\n",rank);CHKERRQ(ierr);
215*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
216*c4762a1bSJed Brown   }
217*c4762a1bSJed Brown 
218*c4762a1bSJed Brown   /* Test MatDiagonalScale() */
219*c4762a1bSJed Brown   ierr = MatDiagonalScale(A,x,x);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = MatDiagonalScale(sA,x,x);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = MatMultEqual(A,sA,10,&flg);CHKERRQ(ierr);
222*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatDiagonalScale");
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   /* Test MatGetDiagonal(), MatScale() */
225*c4762a1bSJed Brown   ierr = MatGetDiagonal(A,s1);CHKERRQ(ierr);
226*c4762a1bSJed Brown   ierr = MatGetDiagonal(sA,s2);CHKERRQ(ierr);
227*c4762a1bSJed Brown   ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
229*c4762a1bSJed Brown   r1  -= r2;
230*c4762a1bSJed Brown   if (r1<-tol || r1>tol) {
231*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n",rank,(double)r1);CHKERRQ(ierr);
232*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
233*c4762a1bSJed Brown   }
234*c4762a1bSJed Brown 
235*c4762a1bSJed Brown   ierr = MatScale(A,alpha);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = MatScale(sA,alpha);CHKERRQ(ierr);
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown   /* Test MatGetRowMaxAbs() */
239*c4762a1bSJed Brown   ierr = MatGetRowMaxAbs(A,s1,NULL);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = MatGetRowMaxAbs(sA,s2,NULL);CHKERRQ(ierr);
241*c4762a1bSJed Brown 
242*c4762a1bSJed Brown   ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
243*c4762a1bSJed Brown   ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
244*c4762a1bSJed Brown   r1  -= r2;
245*c4762a1bSJed Brown   if (r1<-tol || r1>tol) {
246*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRowMaxAbs() \n");CHKERRQ(ierr);
247*c4762a1bSJed Brown   }
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown   /* Test MatMult(), MatMultAdd() */
250*c4762a1bSJed Brown   ierr = MatMultEqual(A,sA,10,&flg);CHKERRQ(ierr);
251*c4762a1bSJed Brown   if (!flg) {
252*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale()\n",rank);CHKERRQ(ierr);
253*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
254*c4762a1bSJed Brown   }
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown   ierr = MatMultAddEqual(A,sA,10,&flg);CHKERRQ(ierr);
257*c4762a1bSJed Brown   if (!flg) {
258*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd()\n",rank);CHKERRQ(ierr);
259*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
260*c4762a1bSJed Brown   }
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   /* Test MatMultTranspose(), MatMultTransposeAdd() */
263*c4762a1bSJed Brown   for (i=0; i<10; i++) {
264*c4762a1bSJed Brown     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
265*c4762a1bSJed Brown     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
266*c4762a1bSJed Brown     ierr = MatMultTranspose(sA,x,s2);CHKERRQ(ierr);
267*c4762a1bSJed Brown     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
268*c4762a1bSJed Brown     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
269*c4762a1bSJed Brown     r1  -= r2;
270*c4762a1bSJed Brown     if (r1<-tol || r1>tol) {
271*c4762a1bSJed Brown       ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%g\n",rank,(double)r1);CHKERRQ(ierr);
272*c4762a1bSJed Brown       ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
273*c4762a1bSJed Brown     }
274*c4762a1bSJed Brown   }
275*c4762a1bSJed Brown   for (i=0; i<10; i++) {
276*c4762a1bSJed Brown     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
277*c4762a1bSJed Brown     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
278*c4762a1bSJed Brown     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
279*c4762a1bSJed Brown     ierr = MatMultTransposeAdd(sA,x,y,s2);CHKERRQ(ierr);
280*c4762a1bSJed Brown     ierr = VecNorm(s1,NORM_1,&r1);CHKERRQ(ierr);
281*c4762a1bSJed Brown     ierr = VecNorm(s2,NORM_1,&r2);CHKERRQ(ierr);
282*c4762a1bSJed Brown     r1  -= r2;
283*c4762a1bSJed Brown     if (r1<-tol || r1>tol) {
284*c4762a1bSJed Brown       ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g \n",rank,(double)r1);CHKERRQ(ierr);
285*c4762a1bSJed Brown       ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
286*c4762a1bSJed Brown     }
287*c4762a1bSJed Brown   }
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown   /* Test MatDuplicate() */
290*c4762a1bSJed Brown   ierr = MatDuplicate(sA,MAT_COPY_VALUES,&sB);CHKERRQ(ierr);
291*c4762a1bSJed Brown   ierr = MatEqual(sA,sB,&flg);CHKERRQ(ierr);
292*c4762a1bSJed Brown   if (!flg) {
293*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD," Error in MatDuplicate(), sA != sB \n");CHKERRQ(ierr);
294*c4762a1bSJed Brown   }
295*c4762a1bSJed Brown   ierr = MatMultEqual(sA,sB,5,&flg);CHKERRQ(ierr);
296*c4762a1bSJed Brown   if (!flg) {
297*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMult()\n",rank);CHKERRQ(ierr);
298*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
299*c4762a1bSJed Brown   }
300*c4762a1bSJed Brown   ierr = MatMultAddEqual(sA,sB,5,&flg);CHKERRQ(ierr);
301*c4762a1bSJed Brown   if (!flg) {
302*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMultAdd(()\n",rank);CHKERRQ(ierr);
303*c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
304*c4762a1bSJed Brown   }
305*c4762a1bSJed Brown   ierr = MatDestroy(&sB);CHKERRQ(ierr);
306*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
307*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
308*c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
309*c4762a1bSJed Brown   ierr = VecDestroy(&s1);CHKERRQ(ierr);
310*c4762a1bSJed Brown   ierr = VecDestroy(&s2);CHKERRQ(ierr);
311*c4762a1bSJed Brown   ierr = MatDestroy(&sA);CHKERRQ(ierr);
312*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
313*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
314*c4762a1bSJed Brown 
315*c4762a1bSJed Brown   ierr = PetscFinalize();
316*c4762a1bSJed Brown   return ierr;
317*c4762a1bSJed Brown }
318*c4762a1bSJed Brown 
319*c4762a1bSJed Brown /*TEST
320*c4762a1bSJed Brown 
321*c4762a1bSJed Brown    test:
322*c4762a1bSJed Brown       nsize: {{1 3}}
323*c4762a1bSJed Brown       args: -bs {{1 2 3  5  7 8}} -mat_ignore_lower_triangular
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown TEST*/
326