xref: /petsc/src/mat/utils/multequal.c (revision f0fc11cebb1bb284829732915f9e84cabc170c2f) !
186a22c91SHong Zhang 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>  /*I   "petscmat.h"  I*/
386a22c91SHong Zhang 
486a22c91SHong Zhang /*@
586a22c91SHong Zhang    MatMultEqual - Compares matrix-vector products of two matrices.
686a22c91SHong Zhang 
786a22c91SHong Zhang    Collective on Mat
886a22c91SHong Zhang 
986a22c91SHong Zhang    Input Parameters:
1086a22c91SHong Zhang +  A - the first matrix
114222ddf1SHong Zhang .  B - the second matrix
1286a22c91SHong Zhang -  n - number of random vectors to be tested
1386a22c91SHong Zhang 
1486a22c91SHong Zhang    Output Parameter:
1586a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
1686a22c91SHong Zhang 
1786a22c91SHong Zhang    Level: intermediate
1886a22c91SHong Zhang 
1986a22c91SHong Zhang @*/
207087cfbeSBarry Smith PetscErrorCode MatMultEqual(Mat A,Mat B,PetscInt n,PetscBool  *flg)
2186a22c91SHong Zhang {
2286a22c91SHong Zhang   PetscErrorCode ierr;
2386a22c91SHong Zhang   Vec            x,s1,s2;
2486a22c91SHong Zhang   PetscRandom    rctx;
252a0b5429SBarry Smith   PetscReal      r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON;
2686a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
27e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
2886a22c91SHong Zhang 
2986a22c91SHong Zhang   PetscFunctionBegin;
300700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
310700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,2);
3286a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
3386a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
34e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
35cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
36ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
37c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
382a7a6963SBarry Smith   ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr);
39cb5d8e9eSHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
4086a22c91SHong Zhang 
414eb6d288SHong Zhang   *flg = PETSC_TRUE;
4286a22c91SHong Zhang   for (k=0; k<n; k++) {
43abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
4486a22c91SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr);
4586a22c91SHong Zhang     ierr = MatMult(B,x,s2);CHKERRQ(ierr);
46e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
47e6e6b1bdSHong Zhang     if (r2 < tol) {
48e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
49e6e6b1bdSHong Zhang     } else {
502dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
51e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
52e6e6b1bdSHong Zhang       r1  /= r2;
53e6e6b1bdSHong Zhang     }
54e6e6b1bdSHong Zhang     if (r1 > tol) {
5586a22c91SHong Zhang       *flg = PETSC_FALSE;
5657622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %D-th MatMult() %g\n",k,(double)r1);CHKERRQ(ierr);
574eb6d288SHong Zhang       break;
5886a22c91SHong Zhang     }
5986a22c91SHong Zhang   }
606bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
616bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
626bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
636bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
6486a22c91SHong Zhang   PetscFunctionReturn(0);
6586a22c91SHong Zhang }
6686a22c91SHong Zhang 
6786a22c91SHong Zhang /*@
6886a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
6986a22c91SHong Zhang 
7086a22c91SHong Zhang    Collective on Mat
7186a22c91SHong Zhang 
7286a22c91SHong Zhang    Input Parameters:
7386a22c91SHong Zhang +  A - the first matrix
744222ddf1SHong Zhang .  B - the second matrix
7586a22c91SHong Zhang -  n - number of random vectors to be tested
7686a22c91SHong Zhang 
7786a22c91SHong Zhang    Output Parameter:
7886a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
7986a22c91SHong Zhang 
8086a22c91SHong Zhang    Level: intermediate
8186a22c91SHong Zhang 
8286a22c91SHong Zhang @*/
837087cfbeSBarry Smith PetscErrorCode  MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscBool  *flg)
8486a22c91SHong Zhang {
8586a22c91SHong Zhang   PetscErrorCode ierr;
8686a22c91SHong Zhang   Vec            x,y,s1,s2;
8786a22c91SHong Zhang   PetscRandom    rctx;
882a0b5429SBarry Smith   PetscReal      r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON;
8986a22c91SHong Zhang   PetscInt       am,an,bm,bn,k;
90e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
9186a22c91SHong Zhang 
9286a22c91SHong Zhang   PetscFunctionBegin;
9386a22c91SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
9486a22c91SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
95e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
96cb5d8e9eSHong Zhang   PetscCheckSameComm(A,1,B,2);
97ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
98c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
992a7a6963SBarry Smith   ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr);
10063879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
10163879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
10286a22c91SHong Zhang 
1034eb6d288SHong Zhang   *flg = PETSC_TRUE;
10486a22c91SHong Zhang   for (k=0; k<n; k++) {
105abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
106abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
10786a22c91SHong Zhang     ierr = MatMultAdd(A,x,y,s1);CHKERRQ(ierr);
10886a22c91SHong Zhang     ierr = MatMultAdd(B,x,y,s2);CHKERRQ(ierr);
109e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
110e6e6b1bdSHong Zhang     if (r2 < tol) {
111e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
112e6e6b1bdSHong Zhang     } else {
1132dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
114e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
115e6e6b1bdSHong Zhang       r1  /= r2;
116e6e6b1bdSHong Zhang     }
117e6e6b1bdSHong Zhang     if (r1 > tol) {
11886a22c91SHong Zhang       *flg = PETSC_FALSE;
11957622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultAdd() %g\n",k,(double)r1);CHKERRQ(ierr);
12063879571SHong Zhang       break;
12163879571SHong Zhang     }
12263879571SHong Zhang   }
1236bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
1246bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1256bf464f9SBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
1266bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
1276bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
12863879571SHong Zhang   PetscFunctionReturn(0);
12963879571SHong Zhang }
13063879571SHong Zhang 
13163879571SHong Zhang /*@
13263879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
13363879571SHong Zhang 
13463879571SHong Zhang    Collective on Mat
13563879571SHong Zhang 
13663879571SHong Zhang    Input Parameters:
13763879571SHong Zhang +  A - the first matrix
1384222ddf1SHong Zhang .  B - the second matrix
13963879571SHong Zhang -  n - number of random vectors to be tested
14063879571SHong Zhang 
14163879571SHong Zhang    Output Parameter:
14263879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
14363879571SHong Zhang 
14463879571SHong Zhang    Level: intermediate
14563879571SHong Zhang 
14663879571SHong Zhang @*/
1477087cfbeSBarry Smith PetscErrorCode  MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscBool  *flg)
14863879571SHong Zhang {
14963879571SHong Zhang   PetscErrorCode ierr;
15063879571SHong Zhang   Vec            x,s1,s2;
15163879571SHong Zhang   PetscRandom    rctx;
1522a0b5429SBarry Smith   PetscReal      r1,r2,tol= PETSC_SQRT_MACHINE_EPSILON;
15363879571SHong Zhang   PetscInt       am,an,bm,bn,k;
154e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
15563879571SHong Zhang 
15663879571SHong Zhang   PetscFunctionBegin;
15763879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
15863879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
159e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
16063879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
161ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
162c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
1632a7a6963SBarry Smith   ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr);
16463879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
16563879571SHong Zhang 
16663879571SHong Zhang   *flg = PETSC_TRUE;
16763879571SHong Zhang   for (k=0; k<n; k++) {
168abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
16963879571SHong Zhang     ierr = MatMultTranspose(A,x,s1);CHKERRQ(ierr);
17063879571SHong Zhang     ierr = MatMultTranspose(B,x,s2);CHKERRQ(ierr);
171e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
172e6e6b1bdSHong Zhang     if (r2 < tol) {
173e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
174e6e6b1bdSHong Zhang     } else {
1752dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
176e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
177e6e6b1bdSHong Zhang       r1  /= r2;
178e6e6b1bdSHong Zhang     }
179e6e6b1bdSHong Zhang     if (r1 > tol) {
18063879571SHong Zhang       *flg = PETSC_FALSE;
18157622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTranspose() %g\n",k,(double)r1);CHKERRQ(ierr);
18263879571SHong Zhang       break;
18363879571SHong Zhang     }
18463879571SHong Zhang   }
1856bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
1866bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1876bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
1886bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
18963879571SHong Zhang   PetscFunctionReturn(0);
19063879571SHong Zhang }
19163879571SHong Zhang 
19263879571SHong Zhang /*@
19363879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
19463879571SHong Zhang 
19563879571SHong Zhang    Collective on Mat
19663879571SHong Zhang 
19763879571SHong Zhang    Input Parameters:
19863879571SHong Zhang +  A - the first matrix
1994222ddf1SHong Zhang .  B - the second matrix
20063879571SHong Zhang -  n - number of random vectors to be tested
20163879571SHong Zhang 
20263879571SHong Zhang    Output Parameter:
20363879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
20463879571SHong Zhang 
20563879571SHong Zhang    Level: intermediate
20663879571SHong Zhang 
20763879571SHong Zhang @*/
2087087cfbeSBarry Smith PetscErrorCode  MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscBool  *flg)
20963879571SHong Zhang {
21063879571SHong Zhang   PetscErrorCode ierr;
21163879571SHong Zhang   Vec            x,y,s1,s2;
21263879571SHong Zhang   PetscRandom    rctx;
2132a0b5429SBarry Smith   PetscReal      r1,r2,tol = PETSC_SQRT_MACHINE_EPSILON;
21463879571SHong Zhang   PetscInt       am,an,bm,bn,k;
215e6e6b1bdSHong Zhang   PetscScalar    none = -1.0;
21663879571SHong Zhang 
21763879571SHong Zhang   PetscFunctionBegin;
21863879571SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
21963879571SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
220e32f2f54SBarry Smith   if (am != bm || an != bn) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D %D %D",am,bm,an,bn);
22163879571SHong Zhang   PetscCheckSameComm(A,1,B,2);
222ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx);CHKERRQ(ierr);
223c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
2242a7a6963SBarry Smith   ierr = MatCreateVecs(A,&s1,&x);CHKERRQ(ierr);
22563879571SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
22663879571SHong Zhang   ierr = VecDuplicate(s1,&y);CHKERRQ(ierr);
22763879571SHong Zhang 
22863879571SHong Zhang   *flg = PETSC_TRUE;
22963879571SHong Zhang   for (k=0; k<n; k++) {
230abb0e124SMatthew Knepley     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
231abb0e124SMatthew Knepley     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
23263879571SHong Zhang     ierr = MatMultTransposeAdd(A,x,y,s1);CHKERRQ(ierr);
23363879571SHong Zhang     ierr = MatMultTransposeAdd(B,x,y,s2);CHKERRQ(ierr);
234e6e6b1bdSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
235e6e6b1bdSHong Zhang     if (r2 < tol) {
236e6e6b1bdSHong Zhang       ierr = VecNorm(s1,NORM_INFINITY,&r1);CHKERRQ(ierr);
237e6e6b1bdSHong Zhang     } else {
2382dcb1b2aSMatthew Knepley       ierr = VecAXPY(s2,none,s1);CHKERRQ(ierr);
239e6e6b1bdSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
240e6e6b1bdSHong Zhang       r1  /= r2;
241e6e6b1bdSHong Zhang     }
242e6e6b1bdSHong Zhang     if (r1 > tol) {
24363879571SHong Zhang       *flg = PETSC_FALSE;
24457622a8eSBarry Smith       ierr = PetscInfo2(A,"Error: %d-th MatMultTransposeAdd() %g\n",k,(double)r1);CHKERRQ(ierr);
2454eb6d288SHong Zhang       break;
24686a22c91SHong Zhang     }
24786a22c91SHong Zhang   }
2486bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
2496bf464f9SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
2506bf464f9SBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
2516bf464f9SBarry Smith   ierr = VecDestroy(&s1);CHKERRQ(ierr);
2526bf464f9SBarry Smith   ierr = VecDestroy(&s2);CHKERRQ(ierr);
25386a22c91SHong Zhang   PetscFunctionReturn(0);
25486a22c91SHong Zhang }
255a52676ecSHong Zhang 
256a52676ecSHong Zhang /*@
257a52676ecSHong Zhang    MatMatMultEqual - Test A*B*x = C*x for n random vector x
258a52676ecSHong Zhang 
259a52676ecSHong Zhang    Collective on Mat
260a52676ecSHong Zhang 
261a52676ecSHong Zhang    Input Parameters:
262a52676ecSHong Zhang +  A - the first matrix
263c2916339SPierre Jolivet .  B - the second matrix
264c2916339SPierre Jolivet .  C - the third matrix
265a52676ecSHong Zhang -  n - number of random vectors to be tested
266a52676ecSHong Zhang 
267a52676ecSHong Zhang    Output Parameter:
268*f0fc11ceSJed Brown .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
269a52676ecSHong Zhang 
270a52676ecSHong Zhang    Level: intermediate
271a52676ecSHong Zhang 
272a52676ecSHong Zhang @*/
273a52676ecSHong Zhang PetscErrorCode MatMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
274a52676ecSHong Zhang {
275a52676ecSHong Zhang   PetscErrorCode ierr;
276a52676ecSHong Zhang   Vec            x,s1,s2,s3;
277a52676ecSHong Zhang   PetscRandom    rctx;
278a52676ecSHong Zhang   PetscReal      r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON;
279a52676ecSHong Zhang   PetscInt       am,an,bm,bn,cm,cn,k;
280a52676ecSHong Zhang   PetscScalar    none = -1.0;
281a52676ecSHong Zhang 
282a52676ecSHong Zhang   PetscFunctionBegin;
283a52676ecSHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
284a52676ecSHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
285a52676ecSHong Zhang   ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr);
286a52676ecSHong Zhang   if (an != bm || am != cm || bn != cn) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %D %D %D %D",am,an,bm,bn,cm, cn);
287a52676ecSHong Zhang 
288a52676ecSHong Zhang   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx);CHKERRQ(ierr);
289a52676ecSHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
290a52676ecSHong Zhang   ierr = MatCreateVecs(B,&x,&s1);CHKERRQ(ierr);
291a52676ecSHong Zhang   ierr = MatCreateVecs(C,NULL,&s2);CHKERRQ(ierr);
292a52676ecSHong Zhang   ierr = VecDuplicate(s2,&s3);CHKERRQ(ierr);
293a52676ecSHong Zhang 
294a52676ecSHong Zhang   *flg = PETSC_TRUE;
295a52676ecSHong Zhang   for (k=0; k<n; k++) {
296a52676ecSHong Zhang     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
297a52676ecSHong Zhang     ierr = MatMult(B,x,s1);CHKERRQ(ierr);
298a52676ecSHong Zhang     ierr = MatMult(A,s1,s2);CHKERRQ(ierr);
299a52676ecSHong Zhang     ierr = MatMult(C,x,s3);CHKERRQ(ierr);
300a52676ecSHong Zhang 
301a52676ecSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
302a52676ecSHong Zhang     if (r2 < tol) {
303a52676ecSHong Zhang       ierr = VecNorm(s3,NORM_INFINITY,&r1);CHKERRQ(ierr);
304a52676ecSHong Zhang     } else {
305a52676ecSHong Zhang       ierr = VecAXPY(s2,none,s3);CHKERRQ(ierr);
306a52676ecSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
307a52676ecSHong Zhang       r1  /= r2;
308a52676ecSHong Zhang     }
309a52676ecSHong Zhang     if (r1 > tol) {
310a52676ecSHong Zhang       *flg = PETSC_FALSE;
311a52676ecSHong Zhang       ierr = PetscInfo2(A,"Error: %D-th MatMatMult() %g\n",k,(double)r1);CHKERRQ(ierr);
312a52676ecSHong Zhang       break;
313a52676ecSHong Zhang     }
314a52676ecSHong Zhang   }
315a52676ecSHong Zhang   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
316a52676ecSHong Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
317a52676ecSHong Zhang   ierr = VecDestroy(&s1);CHKERRQ(ierr);
318a52676ecSHong Zhang   ierr = VecDestroy(&s2);CHKERRQ(ierr);
319a52676ecSHong Zhang   ierr = VecDestroy(&s3);CHKERRQ(ierr);
320a52676ecSHong Zhang   PetscFunctionReturn(0);
321a52676ecSHong Zhang }
322a52676ecSHong Zhang 
323a52676ecSHong Zhang /*@
324a52676ecSHong Zhang    MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x
325a52676ecSHong Zhang 
326a52676ecSHong Zhang    Collective on Mat
327a52676ecSHong Zhang 
328a52676ecSHong Zhang    Input Parameters:
329a52676ecSHong Zhang +  A - the first matrix
3304222ddf1SHong Zhang .  B - the second matrix
3314222ddf1SHong Zhang .  C - the third matrix
332a52676ecSHong Zhang -  n - number of random vectors to be tested
333a52676ecSHong Zhang 
334a52676ecSHong Zhang    Output Parameter:
335a52676ecSHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
336a52676ecSHong Zhang 
337a52676ecSHong Zhang    Level: intermediate
338a52676ecSHong Zhang 
339a52676ecSHong Zhang @*/
340a52676ecSHong Zhang PetscErrorCode MatTransposeMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
341a52676ecSHong Zhang {
342a52676ecSHong Zhang   PetscErrorCode ierr;
343a52676ecSHong Zhang   Vec            x,s1,s2,s3;
344a52676ecSHong Zhang   PetscRandom    rctx;
345a52676ecSHong Zhang   PetscReal      r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON;
346a52676ecSHong Zhang   PetscInt       am,an,bm,bn,cm,cn,k;
347a52676ecSHong Zhang   PetscScalar    none = -1.0;
348a52676ecSHong Zhang 
349a52676ecSHong Zhang   PetscFunctionBegin;
350a52676ecSHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
351a52676ecSHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
352a52676ecSHong Zhang   ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr);
353a52676ecSHong Zhang   if (am != bm || an != cm || bn != cn) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %D %D %D %D",am,an,bm,bn,cm, cn);
354a52676ecSHong Zhang 
355a52676ecSHong Zhang   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx);CHKERRQ(ierr);
356a52676ecSHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
357a52676ecSHong Zhang   ierr = MatCreateVecs(B,&x,&s1);CHKERRQ(ierr);
358a52676ecSHong Zhang   ierr = MatCreateVecs(C,NULL,&s2);CHKERRQ(ierr);
359a52676ecSHong Zhang   ierr = VecDuplicate(s2,&s3);CHKERRQ(ierr);
360a52676ecSHong Zhang 
361a52676ecSHong Zhang   *flg = PETSC_TRUE;
362a52676ecSHong Zhang   for (k=0; k<n; k++) {
363a52676ecSHong Zhang     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
364a52676ecSHong Zhang     ierr = MatMult(B,x,s1);CHKERRQ(ierr);
365a52676ecSHong Zhang     ierr = MatMultTranspose(A,s1,s2);CHKERRQ(ierr);
366a52676ecSHong Zhang     ierr = MatMult(C,x,s3);CHKERRQ(ierr);
367a52676ecSHong Zhang 
368a52676ecSHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
369a52676ecSHong Zhang     if (r2 < tol) {
370a52676ecSHong Zhang       ierr = VecNorm(s3,NORM_INFINITY,&r1);CHKERRQ(ierr);
371a52676ecSHong Zhang     } else {
372a52676ecSHong Zhang       ierr = VecAXPY(s2,none,s3);CHKERRQ(ierr);
373a52676ecSHong Zhang       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
374a52676ecSHong Zhang       r1  /= r2;
375a52676ecSHong Zhang     }
376a52676ecSHong Zhang     if (r1 > tol) {
377a52676ecSHong Zhang       *flg = PETSC_FALSE;
378a52676ecSHong Zhang       ierr = PetscInfo2(A,"Error: %D-th MatTransposeMatMult() %g\n",k,(double)r1);CHKERRQ(ierr);
379a52676ecSHong Zhang       break;
380a52676ecSHong Zhang     }
381a52676ecSHong Zhang   }
382a52676ecSHong Zhang   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
383a52676ecSHong Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
384a52676ecSHong Zhang   ierr = VecDestroy(&s1);CHKERRQ(ierr);
385a52676ecSHong Zhang   ierr = VecDestroy(&s2);CHKERRQ(ierr);
386a52676ecSHong Zhang   ierr = VecDestroy(&s3);CHKERRQ(ierr);
387a52676ecSHong Zhang   PetscFunctionReturn(0);
388a52676ecSHong Zhang }
38986919fd6SHong Zhang 
39026546c1bSToby Isaac /*@
39126546c1bSToby Isaac    MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x
39226546c1bSToby Isaac 
39326546c1bSToby Isaac    Collective on Mat
39426546c1bSToby Isaac 
39526546c1bSToby Isaac    Input Parameters:
39626546c1bSToby Isaac +  A - the first matrix
3974222ddf1SHong Zhang .  B - the second matrix
3984222ddf1SHong Zhang .  C - the third matrix
39926546c1bSToby Isaac -  n - number of random vectors to be tested
40026546c1bSToby Isaac 
40126546c1bSToby Isaac    Output Parameter:
40226546c1bSToby Isaac .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
40326546c1bSToby Isaac 
40426546c1bSToby Isaac    Level: intermediate
40526546c1bSToby Isaac 
40626546c1bSToby Isaac @*/
407cc48ffa7SToby Isaac PetscErrorCode MatMatTransposeMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
408cc48ffa7SToby Isaac {
409cc48ffa7SToby Isaac   PetscErrorCode ierr;
410cc48ffa7SToby Isaac   Vec            x,s1,s2,s3;
411cc48ffa7SToby Isaac   PetscRandom    rctx;
412cc48ffa7SToby Isaac   PetscReal      r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON;
413cc48ffa7SToby Isaac   PetscInt       am,an,bm,bn,cm,cn,k;
414cc48ffa7SToby Isaac   PetscScalar    none = -1.0;
415cc48ffa7SToby Isaac 
416cc48ffa7SToby Isaac   PetscFunctionBegin;
417cc48ffa7SToby Isaac   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
418cc48ffa7SToby Isaac   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
419cc48ffa7SToby Isaac   ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr);
420cc48ffa7SToby Isaac   if (an != bn || am != cm || bm != cn) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %D %D %D %D",am,an,bm,bn,cm, cn);
421cc48ffa7SToby Isaac 
422cc48ffa7SToby Isaac   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx);CHKERRQ(ierr);
423cc48ffa7SToby Isaac   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
424cc48ffa7SToby Isaac   ierr = MatCreateVecs(B,&s1,&x);CHKERRQ(ierr);
425cc48ffa7SToby Isaac   ierr = MatCreateVecs(C,NULL,&s2);CHKERRQ(ierr);
426cc48ffa7SToby Isaac   ierr = VecDuplicate(s2,&s3);CHKERRQ(ierr);
427cc48ffa7SToby Isaac 
428cc48ffa7SToby Isaac   *flg = PETSC_TRUE;
429cc48ffa7SToby Isaac   for (k=0; k<n; k++) {
430cc48ffa7SToby Isaac     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
431cc48ffa7SToby Isaac     ierr = MatMultTranspose(B,x,s1);CHKERRQ(ierr);
432cc48ffa7SToby Isaac     ierr = MatMult(A,s1,s2);CHKERRQ(ierr);
433cc48ffa7SToby Isaac     ierr = MatMult(C,x,s3);CHKERRQ(ierr);
434cc48ffa7SToby Isaac 
435cc48ffa7SToby Isaac     ierr = VecNorm(s2,NORM_INFINITY,&r2);CHKERRQ(ierr);
436cc48ffa7SToby Isaac     if (r2 < tol) {
437cc48ffa7SToby Isaac       ierr = VecNorm(s3,NORM_INFINITY,&r1);CHKERRQ(ierr);
438cc48ffa7SToby Isaac     } else {
439cc48ffa7SToby Isaac       ierr = VecAXPY(s2,none,s3);CHKERRQ(ierr);
440cc48ffa7SToby Isaac       ierr = VecNorm(s2,NORM_INFINITY,&r1);CHKERRQ(ierr);
441cc48ffa7SToby Isaac       r1  /= r2;
442cc48ffa7SToby Isaac     }
443cc48ffa7SToby Isaac     if (r1 > tol) {
444cc48ffa7SToby Isaac       *flg = PETSC_FALSE;
445cc48ffa7SToby Isaac       ierr = PetscInfo2(A,"Error: %D-th MatMatTransposeMult() %g\n",k,(double)r1);CHKERRQ(ierr);
446cc48ffa7SToby Isaac       break;
447cc48ffa7SToby Isaac     }
448cc48ffa7SToby Isaac   }
449cc48ffa7SToby Isaac   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
450cc48ffa7SToby Isaac   ierr = VecDestroy(&x);CHKERRQ(ierr);
451cc48ffa7SToby Isaac   ierr = VecDestroy(&s1);CHKERRQ(ierr);
452cc48ffa7SToby Isaac   ierr = VecDestroy(&s2);CHKERRQ(ierr);
453cc48ffa7SToby Isaac   ierr = VecDestroy(&s3);CHKERRQ(ierr);
454cc48ffa7SToby Isaac   PetscFunctionReturn(0);
455cc48ffa7SToby Isaac }
456cc48ffa7SToby Isaac 
45786919fd6SHong Zhang /*@
4584222ddf1SHong Zhang    MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B
4594222ddf1SHong Zhang 
4604222ddf1SHong Zhang    Collective on Mat
4614222ddf1SHong Zhang 
4624222ddf1SHong Zhang    Input Parameters:
4634222ddf1SHong Zhang +  A - the first matrix
4644222ddf1SHong Zhang .  B - the second matrix
4654222ddf1SHong Zhang .  C - the third matrix
4664222ddf1SHong Zhang -  n - number of random vectors to be tested
4674222ddf1SHong Zhang 
4684222ddf1SHong Zhang    Output Parameter:
4694222ddf1SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
4704222ddf1SHong Zhang 
4714222ddf1SHong Zhang    Level: intermediate
4724222ddf1SHong Zhang 
4734222ddf1SHong Zhang @*/
4744222ddf1SHong Zhang PetscErrorCode MatPtAPMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
4754222ddf1SHong Zhang {
4764222ddf1SHong Zhang   PetscErrorCode ierr;
4774222ddf1SHong Zhang   Vec            x,v1,v2,v3,v4;
4784222ddf1SHong Zhang   PetscReal      norm_abs,norm_rel,tol=PETSC_SQRT_MACHINE_EPSILON;
4794222ddf1SHong Zhang   PetscInt       i,am,an,bm,bn,cm,cn;
4804222ddf1SHong Zhang   PetscRandom    rdm;
4814222ddf1SHong Zhang   PetscScalar    none = -1.0;
4824222ddf1SHong Zhang 
4834222ddf1SHong Zhang   PetscFunctionBegin;
4844222ddf1SHong Zhang   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
4854222ddf1SHong Zhang   ierr = MatGetLocalSize(B,&bm,&bn);CHKERRQ(ierr);
4864222ddf1SHong Zhang   ierr = MatGetLocalSize(C,&cm,&cn);CHKERRQ(ierr);
4874222ddf1SHong Zhang   if (an != bm || bn != cm || bn != cn) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %D %D %D %D %D %D",am,an,bm,bn,cm,cn);
4884222ddf1SHong Zhang 
4894222ddf1SHong Zhang   /* Create left vector of A: v2 */
4904222ddf1SHong Zhang   ierr = MatCreateVecs(A,NULL,&v2);CHKERRQ(ierr);
4914222ddf1SHong Zhang 
4924222ddf1SHong Zhang   /* Create right vectors of B: x, v3, v4 */
4934222ddf1SHong Zhang   ierr = MatCreateVecs(B,&x,&v1);CHKERRQ(ierr);
4944222ddf1SHong Zhang   ierr = VecDuplicate(x,&v3);CHKERRQ(ierr);
4954222ddf1SHong Zhang   ierr = VecDuplicate(x,&v4);CHKERRQ(ierr);
4964222ddf1SHong Zhang 
4974222ddf1SHong Zhang   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr);
4984222ddf1SHong Zhang   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
4994222ddf1SHong Zhang 
5004222ddf1SHong Zhang   *flg = PETSC_TRUE;
5014222ddf1SHong Zhang   for (i=0; i<n; i++) {
5024222ddf1SHong Zhang     ierr = VecSetRandom(x,rdm);CHKERRQ(ierr);
5034222ddf1SHong Zhang     ierr = MatMult(B,x,v1);CHKERRQ(ierr);
5044222ddf1SHong Zhang     ierr = MatMult(A,v1,v2);CHKERRQ(ierr);          /* v2 = A*B*x */
5054222ddf1SHong Zhang 
5064222ddf1SHong Zhang     ierr = MatMultTranspose(B,v2,v3);CHKERRQ(ierr); /* v3 = Bt*A*B*x */
5074222ddf1SHong Zhang     ierr = MatMult(C,x,v4);CHKERRQ(ierr);           /* v4 = C*x   */
5084222ddf1SHong Zhang     ierr = VecNorm(v4,NORM_2,&norm_abs);CHKERRQ(ierr);
5094222ddf1SHong Zhang     ierr = VecAXPY(v4,none,v3);CHKERRQ(ierr);
5104222ddf1SHong Zhang     ierr = VecNorm(v4,NORM_2,&norm_rel);CHKERRQ(ierr);
5114222ddf1SHong Zhang 
5124222ddf1SHong Zhang     if (norm_abs > tol) norm_rel /= norm_abs;
5134222ddf1SHong Zhang     if (norm_rel > tol) {
5144222ddf1SHong Zhang       *flg = PETSC_FALSE;
5154222ddf1SHong Zhang       ierr = PetscInfo2(A,"Error: %D-th MatPtAPMult() %g\n",i,(double)norm_rel);CHKERRQ(ierr);
5164222ddf1SHong Zhang       break;
5174222ddf1SHong Zhang     }
5184222ddf1SHong Zhang   }
5194222ddf1SHong Zhang 
5204222ddf1SHong Zhang   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
5214222ddf1SHong Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
5224222ddf1SHong Zhang   ierr = VecDestroy(&v1);CHKERRQ(ierr);
5234222ddf1SHong Zhang   ierr = VecDestroy(&v2);CHKERRQ(ierr);
5244222ddf1SHong Zhang   ierr = VecDestroy(&v3);CHKERRQ(ierr);
5254222ddf1SHong Zhang   ierr = VecDestroy(&v4);CHKERRQ(ierr);
5264222ddf1SHong Zhang   PetscFunctionReturn(0);
5274222ddf1SHong Zhang }
5284222ddf1SHong Zhang 
5294222ddf1SHong Zhang /*@
53086919fd6SHong Zhang    MatIsLinear - Check if a shell matrix A is a linear operator.
53186919fd6SHong Zhang 
53286919fd6SHong Zhang    Collective on Mat
53386919fd6SHong Zhang 
53486919fd6SHong Zhang    Input Parameters:
53586919fd6SHong Zhang +  A - the shell matrix
53686919fd6SHong Zhang -  n - number of random vectors to be tested
53786919fd6SHong Zhang 
53886919fd6SHong Zhang    Output Parameter:
53986919fd6SHong Zhang .  flg - PETSC_TRUE if the shell matrix is linear; PETSC_FALSE otherwise.
54086919fd6SHong Zhang 
54186919fd6SHong Zhang    Level: intermediate
54286919fd6SHong Zhang @*/
54386919fd6SHong Zhang PetscErrorCode MatIsLinear(Mat A,PetscInt n,PetscBool  *flg)
54486919fd6SHong Zhang {
54586919fd6SHong Zhang   PetscErrorCode ierr;
54686919fd6SHong Zhang   Vec            x,y,s1,s2;
54786919fd6SHong Zhang   PetscRandom    rctx;
54886919fd6SHong Zhang   PetscScalar    a;
54986919fd6SHong Zhang   PetscInt       k;
55086919fd6SHong Zhang   PetscReal      norm,normA;
55186919fd6SHong Zhang   MPI_Comm       comm;
55286919fd6SHong Zhang   PetscMPIInt    rank;
55386919fd6SHong Zhang 
55486919fd6SHong Zhang   PetscFunctionBegin;
55586919fd6SHong Zhang   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
55686919fd6SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
55786919fd6SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
55886919fd6SHong Zhang 
55986919fd6SHong Zhang   ierr = PetscRandomCreate(comm,&rctx);CHKERRQ(ierr);
56086919fd6SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
56186919fd6SHong Zhang   ierr = MatCreateVecs(A,&x,&s1);CHKERRQ(ierr);
56286919fd6SHong Zhang   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
56386919fd6SHong Zhang   ierr = VecDuplicate(s1,&s2);CHKERRQ(ierr);
56486919fd6SHong Zhang 
56586919fd6SHong Zhang   *flg = PETSC_TRUE;
56686919fd6SHong Zhang   for (k=0; k<n; k++) {
56786919fd6SHong Zhang     ierr = VecSetRandom(x,rctx);CHKERRQ(ierr);
5686d5db48cSHong Zhang     ierr = VecSetRandom(y,rctx);CHKERRQ(ierr);
56986919fd6SHong Zhang     if (!rank) {
57086919fd6SHong Zhang       ierr = PetscRandomGetValue(rctx,&a);CHKERRQ(ierr);
57186919fd6SHong Zhang     }
57286919fd6SHong Zhang     ierr = MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm);CHKERRQ(ierr);
57386919fd6SHong Zhang 
57486919fd6SHong Zhang     /* s2 = a*A*x + A*y */
57586919fd6SHong Zhang     ierr = MatMult(A,y,s2);CHKERRQ(ierr); /* s2 = A*y */
57686919fd6SHong Zhang     ierr = MatMult(A,x,s1);CHKERRQ(ierr); /* s1 = A*x */
57786919fd6SHong Zhang     ierr = VecAXPY(s2,a,s1);CHKERRQ(ierr); /* s2 = a s1 + s2 */
57886919fd6SHong Zhang 
57986919fd6SHong Zhang     /* s1 = A * (a x + y) */
58086919fd6SHong Zhang     ierr = VecAXPY(y,a,x);CHKERRQ(ierr); /* y = a x + y */
58186919fd6SHong Zhang     ierr = MatMult(A,y,s1);CHKERRQ(ierr);
58286919fd6SHong Zhang     ierr = VecNorm(s1,NORM_INFINITY,&normA);CHKERRQ(ierr);
58386919fd6SHong Zhang 
58486919fd6SHong Zhang     ierr = VecAXPY(s2,-1.0,s1);CHKERRQ(ierr); /* s2 = - s1 + s2 */
58586919fd6SHong Zhang     ierr = VecNorm(s2,NORM_INFINITY,&norm);CHKERRQ(ierr);
5861b8dac88SHong Zhang     if (norm/normA > 100.*PETSC_MACHINE_EPSILON) {
58786919fd6SHong Zhang       *flg = PETSC_FALSE;
5881b8dac88SHong Zhang       ierr = PetscInfo3(A,"Error: %D-th |A*(ax+y) - (a*A*x+A*y)|/|A(ax+y)| %g > tol %g\n",k,(double)norm/normA,100.*PETSC_MACHINE_EPSILON);CHKERRQ(ierr);
58986919fd6SHong Zhang       break;
59086919fd6SHong Zhang     }
59186919fd6SHong Zhang   }
59286919fd6SHong Zhang   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
59386919fd6SHong Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
59486919fd6SHong Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
59586919fd6SHong Zhang   ierr = VecDestroy(&s1);CHKERRQ(ierr);
59686919fd6SHong Zhang   ierr = VecDestroy(&s2);CHKERRQ(ierr);
59786919fd6SHong Zhang   PetscFunctionReturn(0);
59886919fd6SHong Zhang }
599