xref: /petsc/src/mat/utils/multequal.c (revision dadcf80911fb48939c55327431ae8d7e47dbe367)
186a22c91SHong Zhang 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>  /*I   "petscmat.h"  I*/
386a22c91SHong Zhang 
4e573050aSPierre Jolivet static PetscErrorCode MatMultEqual_Private(Mat A,Mat B,PetscInt n,PetscBool *flg,PetscInt t,PetscBool add)
5447fed29SStefano Zampini {
6b84f494bSStefano Zampini   Vec            Ax = NULL,Bx = NULL,s1 = NULL,s2 = NULL,Ay = NULL, By = NULL;
7447fed29SStefano Zampini   PetscRandom    rctx;
8447fed29SStefano Zampini   PetscReal      r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON;
9447fed29SStefano Zampini   PetscInt       am,an,bm,bn,k;
10447fed29SStefano Zampini   PetscScalar    none = -1.0;
11580c7c76SPierre Jolivet   const char*    sops[] = {"MatMult","MatMultAdd","MatMultTranspose","MatMultTransposeAdd","MatMultHermitianTranspose","MatMultHermitianTransposeAdd"};
12447fed29SStefano Zampini   const char*    sop;
13447fed29SStefano Zampini 
14447fed29SStefano Zampini   PetscFunctionBegin;
15447fed29SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
16447fed29SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,2);
17447fed29SStefano Zampini   PetscCheckSameComm(A,1,B,2);
18447fed29SStefano Zampini   PetscValidLogicalCollectiveInt(A,n,3);
19*dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(flg,4);
20e573050aSPierre Jolivet   PetscValidLogicalCollectiveInt(A,t,5);
21447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(A,add,6);
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&am,&an));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(B,&bm,&bn));
242c71b3e2SJacob Faibussowitsch   PetscCheckFalse(am != bm || an != bn,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT,am,bm,an,bn);
25e573050aSPierre Jolivet   sop  = sops[(add ? 1 : 0) + 2 * t]; /* t = 0 => no transpose, t = 1 => transpose, t = 2 => Hermitian transpose */
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject)A),&rctx));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rctx));
28447fed29SStefano Zampini   if (t) {
295f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(A,&s1,&Ax));
305f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,&s2,&Bx));
31447fed29SStefano Zampini   } else {
325f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(A,&Ax,&s1));
335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,&Bx,&s2));
34447fed29SStefano Zampini   }
35447fed29SStefano Zampini   if (add) {
365f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(s1,&Ay));
375f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(s2,&By));
38447fed29SStefano Zampini   }
39447fed29SStefano Zampini 
40447fed29SStefano Zampini   *flg = PETSC_TRUE;
41447fed29SStefano Zampini   for (k=0; k<n; k++) {
425f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(Ax,rctx));
435f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(Ax,Bx));
44447fed29SStefano Zampini     if (add) {
455f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(Ay,rctx));
465f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(Ay,By));
47447fed29SStefano Zampini     }
48e573050aSPierre Jolivet     if (t == 1) {
49447fed29SStefano Zampini       if (add) {
505f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultTransposeAdd(A,Ax,Ay,s1));
515f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultTransposeAdd(B,Bx,By,s2));
52447fed29SStefano Zampini       } else {
535f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultTranspose(A,Ax,s1));
545f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultTranspose(B,Bx,s2));
55447fed29SStefano Zampini       }
56e573050aSPierre Jolivet     } else if (t == 2) {
57e573050aSPierre Jolivet       if (add) {
585f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultHermitianTransposeAdd(A,Ax,Ay,s1));
595f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultHermitianTransposeAdd(B,Bx,By,s2));
60e573050aSPierre Jolivet       } else {
615f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultHermitianTranspose(A,Ax,s1));
625f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultHermitianTranspose(B,Bx,s2));
63e573050aSPierre Jolivet       }
64447fed29SStefano Zampini     } else {
65447fed29SStefano Zampini       if (add) {
665f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultAdd(A,Ax,Ay,s1));
675f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMultAdd(B,Bx,By,s2));
68447fed29SStefano Zampini       } else {
695f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMult(A,Ax,s1));
705f80ce2aSJacob Faibussowitsch         CHKERRQ(MatMult(B,Bx,s2));
71447fed29SStefano Zampini       }
72447fed29SStefano Zampini     }
735f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_INFINITY,&r2));
74447fed29SStefano Zampini     if (r2 < tol) {
755f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s1,NORM_INFINITY,&r1));
76447fed29SStefano Zampini     } else {
775f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(s2,none,s1));
785f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_INFINITY,&r1));
79447fed29SStefano Zampini       r1  /= r2;
80447fed29SStefano Zampini     }
81447fed29SStefano Zampini     if (r1 > tol) {
82447fed29SStefano Zampini       *flg = PETSC_FALSE;
835f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"Error: %" PetscInt_FMT "-th %s() %g\n",k,sop,(double)r1));
84447fed29SStefano Zampini       break;
85447fed29SStefano Zampini     }
86447fed29SStefano Zampini   }
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rctx));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Ax));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Bx));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Ay));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&By));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s1));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s2));
94447fed29SStefano Zampini   PetscFunctionReturn(0);
95447fed29SStefano Zampini }
96447fed29SStefano Zampini 
97447fed29SStefano Zampini static PetscErrorCode MatMatMultEqual_Private(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg,PetscBool At,PetscBool Bt)
98447fed29SStefano Zampini {
99447fed29SStefano Zampini   Vec            Ax,Bx,Cx,s1,s2,s3;
100447fed29SStefano Zampini   PetscRandom    rctx;
101447fed29SStefano Zampini   PetscReal      r1,r2,tol=PETSC_SQRT_MACHINE_EPSILON;
102447fed29SStefano Zampini   PetscInt       am,an,bm,bn,cm,cn,k;
103447fed29SStefano Zampini   PetscScalar    none = -1.0;
104580c7c76SPierre Jolivet   const char*    sops[] = {"MatMatMult","MatTransposeMatMult","MatMatTransposeMult","MatTransposeMatTransposeMult"};
105447fed29SStefano Zampini   const char*    sop;
106447fed29SStefano Zampini 
107447fed29SStefano Zampini   PetscFunctionBegin;
108447fed29SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
109447fed29SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,2);
110447fed29SStefano Zampini   PetscCheckSameComm(A,1,B,2);
111447fed29SStefano Zampini   PetscValidHeaderSpecific(C,MAT_CLASSID,3);
112447fed29SStefano Zampini   PetscCheckSameComm(A,1,C,3);
113447fed29SStefano Zampini   PetscValidLogicalCollectiveInt(A,n,4);
114*dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(flg,5);
115447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(A,At,6);
116447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(B,Bt,7);
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&am,&an));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(B,&bm,&bn));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(C,&cm,&cn));
120447fed29SStefano Zampini   if (At) { PetscInt tt = an; an = am; am = tt; };
121447fed29SStefano Zampini   if (Bt) { PetscInt tt = bn; bn = bm; bm = tt; };
1222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(an != bm || am != cm || bn != cn,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT,am,an,bm,bn,cm,cn);
123447fed29SStefano Zampini 
124447fed29SStefano Zampini   sop  = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)];
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject)C),&rctx));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rctx));
127447fed29SStefano Zampini   if (Bt) {
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,&s1,&Bx));
129447fed29SStefano Zampini   } else {
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,&Bx,&s1));
131447fed29SStefano Zampini   }
132447fed29SStefano Zampini   if (At) {
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(A,&s2,&Ax));
134447fed29SStefano Zampini   } else {
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(A,&Ax,&s2));
136447fed29SStefano Zampini   }
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(C,&Cx,&s3));
138447fed29SStefano Zampini 
139447fed29SStefano Zampini   *flg = PETSC_TRUE;
140447fed29SStefano Zampini   for (k=0; k<n; k++) {
1415f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(Bx,rctx));
142447fed29SStefano Zampini     if (Bt) {
1435f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultTranspose(B,Bx,s1));
144447fed29SStefano Zampini     } else {
1455f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(B,Bx,s1));
146447fed29SStefano Zampini     }
1475f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(s1,Ax));
148447fed29SStefano Zampini     if (At) {
1495f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultTranspose(A,Ax,s2));
150447fed29SStefano Zampini     } else {
1515f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(A,Ax,s2));
152447fed29SStefano Zampini     }
1535f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(Bx,Cx));
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(C,Cx,s3));
155447fed29SStefano Zampini 
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_INFINITY,&r2));
157447fed29SStefano Zampini     if (r2 < tol) {
1585f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s3,NORM_INFINITY,&r1));
159447fed29SStefano Zampini     } else {
1605f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(s2,none,s3));
1615f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(s2,NORM_INFINITY,&r1));
162447fed29SStefano Zampini       r1  /= r2;
163447fed29SStefano Zampini     }
164447fed29SStefano Zampini     if (r1 > tol) {
165447fed29SStefano Zampini       *flg = PETSC_FALSE;
1665f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"Error: %" PetscInt_FMT "-th %s %g\n",k,sop,(double)r1));
167447fed29SStefano Zampini       break;
168447fed29SStefano Zampini     }
169447fed29SStefano Zampini   }
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rctx));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Ax));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Bx));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Cx));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s1));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s2));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s3));
177447fed29SStefano Zampini   PetscFunctionReturn(0);
178447fed29SStefano Zampini }
179447fed29SStefano Zampini 
18086a22c91SHong Zhang /*@
18186a22c91SHong Zhang    MatMultEqual - Compares matrix-vector products of two matrices.
18286a22c91SHong Zhang 
18386a22c91SHong Zhang    Collective on Mat
18486a22c91SHong Zhang 
18586a22c91SHong Zhang    Input Parameters:
18686a22c91SHong Zhang +  A - the first matrix
1874222ddf1SHong Zhang .  B - the second matrix
18886a22c91SHong Zhang -  n - number of random vectors to be tested
18986a22c91SHong Zhang 
19086a22c91SHong Zhang    Output Parameter:
19186a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
19286a22c91SHong Zhang 
19386a22c91SHong Zhang    Level: intermediate
19486a22c91SHong Zhang 
19586a22c91SHong Zhang @*/
1967087cfbeSBarry Smith PetscErrorCode MatMultEqual(Mat A,Mat B,PetscInt n,PetscBool *flg)
19786a22c91SHong Zhang {
19886a22c91SHong Zhang   PetscFunctionBegin;
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual_Private(A,B,n,flg,0,PETSC_FALSE));
20086a22c91SHong Zhang   PetscFunctionReturn(0);
20186a22c91SHong Zhang }
20286a22c91SHong Zhang 
20386a22c91SHong Zhang /*@
20486a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
20586a22c91SHong Zhang 
20686a22c91SHong Zhang    Collective on Mat
20786a22c91SHong Zhang 
20886a22c91SHong Zhang    Input Parameters:
20986a22c91SHong Zhang +  A - the first matrix
2104222ddf1SHong Zhang .  B - the second matrix
21186a22c91SHong Zhang -  n - number of random vectors to be tested
21286a22c91SHong Zhang 
21386a22c91SHong Zhang    Output Parameter:
21486a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
21586a22c91SHong Zhang 
21686a22c91SHong Zhang    Level: intermediate
21786a22c91SHong Zhang 
21886a22c91SHong Zhang @*/
2197087cfbeSBarry Smith PetscErrorCode  MatMultAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg)
22086a22c91SHong Zhang {
22186a22c91SHong Zhang   PetscFunctionBegin;
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual_Private(A,B,n,flg,0,PETSC_TRUE));
22363879571SHong Zhang   PetscFunctionReturn(0);
22463879571SHong Zhang }
22563879571SHong Zhang 
22663879571SHong Zhang /*@
22763879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
22863879571SHong Zhang 
22963879571SHong Zhang    Collective on Mat
23063879571SHong Zhang 
23163879571SHong Zhang    Input Parameters:
23263879571SHong Zhang +  A - the first matrix
2334222ddf1SHong Zhang .  B - the second matrix
23463879571SHong Zhang -  n - number of random vectors to be tested
23563879571SHong Zhang 
23663879571SHong Zhang    Output Parameter:
23763879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
23863879571SHong Zhang 
23963879571SHong Zhang    Level: intermediate
24063879571SHong Zhang 
24163879571SHong Zhang @*/
2427087cfbeSBarry Smith PetscErrorCode  MatMultTransposeEqual(Mat A,Mat B,PetscInt n,PetscBool *flg)
24363879571SHong Zhang {
24463879571SHong Zhang   PetscFunctionBegin;
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual_Private(A,B,n,flg,1,PETSC_FALSE));
24663879571SHong Zhang   PetscFunctionReturn(0);
24763879571SHong Zhang }
24863879571SHong Zhang 
24963879571SHong Zhang /*@
25063879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
25163879571SHong Zhang 
25263879571SHong Zhang    Collective on Mat
25363879571SHong Zhang 
25463879571SHong Zhang    Input Parameters:
25563879571SHong Zhang +  A - the first matrix
2564222ddf1SHong Zhang .  B - the second matrix
25763879571SHong Zhang -  n - number of random vectors to be tested
25863879571SHong Zhang 
25963879571SHong Zhang    Output Parameter:
26063879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
26163879571SHong Zhang 
26263879571SHong Zhang    Level: intermediate
26363879571SHong Zhang 
26463879571SHong Zhang @*/
2657087cfbeSBarry Smith PetscErrorCode  MatMultTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg)
26663879571SHong Zhang {
26763879571SHong Zhang   PetscFunctionBegin;
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual_Private(A,B,n,flg,1,PETSC_TRUE));
269e573050aSPierre Jolivet   PetscFunctionReturn(0);
270e573050aSPierre Jolivet }
271e573050aSPierre Jolivet 
272e573050aSPierre Jolivet /*@
273e573050aSPierre Jolivet    MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices.
274e573050aSPierre Jolivet 
275e573050aSPierre Jolivet    Collective on Mat
276e573050aSPierre Jolivet 
277e573050aSPierre Jolivet    Input Parameters:
278e573050aSPierre Jolivet +  A - the first matrix
279e573050aSPierre Jolivet .  B - the second matrix
280e573050aSPierre Jolivet -  n - number of random vectors to be tested
281e573050aSPierre Jolivet 
282e573050aSPierre Jolivet    Output Parameter:
283e573050aSPierre Jolivet .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
284e573050aSPierre Jolivet 
285e573050aSPierre Jolivet    Level: intermediate
286e573050aSPierre Jolivet 
287e573050aSPierre Jolivet @*/
288e573050aSPierre Jolivet PetscErrorCode  MatMultHermitianTransposeEqual(Mat A,Mat B,PetscInt n,PetscBool *flg)
289e573050aSPierre Jolivet {
290e573050aSPierre Jolivet   PetscFunctionBegin;
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual_Private(A,B,n,flg,2,PETSC_FALSE));
292e573050aSPierre Jolivet   PetscFunctionReturn(0);
293e573050aSPierre Jolivet }
294e573050aSPierre Jolivet 
295e573050aSPierre Jolivet /*@
296e573050aSPierre Jolivet    MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices.
297e573050aSPierre Jolivet 
298e573050aSPierre Jolivet    Collective on Mat
299e573050aSPierre Jolivet 
300e573050aSPierre Jolivet    Input Parameters:
301e573050aSPierre Jolivet +  A - the first matrix
302e573050aSPierre Jolivet .  B - the second matrix
303e573050aSPierre Jolivet -  n - number of random vectors to be tested
304e573050aSPierre Jolivet 
305e573050aSPierre Jolivet    Output Parameter:
306e573050aSPierre Jolivet .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
307e573050aSPierre Jolivet 
308e573050aSPierre Jolivet    Level: intermediate
309e573050aSPierre Jolivet 
310e573050aSPierre Jolivet @*/
311e573050aSPierre Jolivet PetscErrorCode  MatMultHermitianTransposeAddEqual(Mat A,Mat B,PetscInt n,PetscBool *flg)
312e573050aSPierre Jolivet {
313e573050aSPierre Jolivet   PetscFunctionBegin;
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual_Private(A,B,n,flg,2,PETSC_TRUE));
31586a22c91SHong Zhang   PetscFunctionReturn(0);
31686a22c91SHong Zhang }
317a52676ecSHong Zhang 
318a52676ecSHong Zhang /*@
319a52676ecSHong Zhang    MatMatMultEqual - Test A*B*x = C*x for n random vector x
320a52676ecSHong Zhang 
321a52676ecSHong Zhang    Collective on Mat
322a52676ecSHong Zhang 
323a52676ecSHong Zhang    Input Parameters:
324a52676ecSHong Zhang +  A - the first matrix
325c2916339SPierre Jolivet .  B - the second matrix
326c2916339SPierre Jolivet .  C - the third matrix
327a52676ecSHong Zhang -  n - number of random vectors to be tested
328a52676ecSHong Zhang 
329a52676ecSHong Zhang    Output Parameter:
330f0fc11ceSJed Brown .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
331a52676ecSHong Zhang 
332a52676ecSHong Zhang    Level: intermediate
333a52676ecSHong Zhang 
334a52676ecSHong Zhang @*/
335a52676ecSHong Zhang PetscErrorCode MatMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
336a52676ecSHong Zhang {
337a52676ecSHong Zhang   PetscFunctionBegin;
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual_Private(A,B,C,n,flg,PETSC_FALSE,PETSC_FALSE));
339a52676ecSHong Zhang   PetscFunctionReturn(0);
340a52676ecSHong Zhang }
341a52676ecSHong Zhang 
342a52676ecSHong Zhang /*@
343a52676ecSHong Zhang    MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x
344a52676ecSHong Zhang 
345a52676ecSHong Zhang    Collective on Mat
346a52676ecSHong Zhang 
347a52676ecSHong Zhang    Input Parameters:
348a52676ecSHong Zhang +  A - the first matrix
3494222ddf1SHong Zhang .  B - the second matrix
3504222ddf1SHong Zhang .  C - the third matrix
351a52676ecSHong Zhang -  n - number of random vectors to be tested
352a52676ecSHong Zhang 
353a52676ecSHong Zhang    Output Parameter:
354a52676ecSHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
355a52676ecSHong Zhang 
356a52676ecSHong Zhang    Level: intermediate
357a52676ecSHong Zhang 
358a52676ecSHong Zhang @*/
359a52676ecSHong Zhang PetscErrorCode MatTransposeMatMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
360a52676ecSHong Zhang {
361a52676ecSHong Zhang   PetscFunctionBegin;
3625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual_Private(A,B,C,n,flg,PETSC_TRUE,PETSC_FALSE));
363a52676ecSHong Zhang   PetscFunctionReturn(0);
364a52676ecSHong Zhang }
36586919fd6SHong Zhang 
36626546c1bSToby Isaac /*@
36726546c1bSToby Isaac    MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x
36826546c1bSToby Isaac 
36926546c1bSToby Isaac    Collective on Mat
37026546c1bSToby Isaac 
37126546c1bSToby Isaac    Input Parameters:
37226546c1bSToby Isaac +  A - the first matrix
3734222ddf1SHong Zhang .  B - the second matrix
3744222ddf1SHong Zhang .  C - the third matrix
37526546c1bSToby Isaac -  n - number of random vectors to be tested
37626546c1bSToby Isaac 
37726546c1bSToby Isaac    Output Parameter:
37826546c1bSToby Isaac .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
37926546c1bSToby Isaac 
38026546c1bSToby Isaac    Level: intermediate
38126546c1bSToby Isaac 
38226546c1bSToby Isaac @*/
383cc48ffa7SToby Isaac PetscErrorCode MatMatTransposeMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
384cc48ffa7SToby Isaac {
385447fed29SStefano Zampini   PetscFunctionBegin;
3865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual_Private(A,B,C,n,flg,PETSC_FALSE,PETSC_TRUE));
387447fed29SStefano Zampini   PetscFunctionReturn(0);
388447fed29SStefano Zampini }
389447fed29SStefano Zampini 
390447fed29SStefano Zampini static PetscErrorCode MatProjMultEqual_Private(Mat A,Mat B,Mat C,PetscInt n,PetscBool rart,PetscBool *flg)
391447fed29SStefano Zampini {
392447fed29SStefano Zampini   Vec            x,v1,v2,v3,v4,Cx,Bx;
393447fed29SStefano Zampini   PetscReal      norm_abs,norm_rel,tol=PETSC_SQRT_MACHINE_EPSILON;
394447fed29SStefano Zampini   PetscInt       i,am,an,bm,bn,cm,cn;
395447fed29SStefano Zampini   PetscRandom    rdm;
396cc48ffa7SToby Isaac   PetscScalar    none = -1.0;
397cc48ffa7SToby Isaac 
398cc48ffa7SToby Isaac   PetscFunctionBegin;
3995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&am,&an));
4005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(B,&bm,&bn));
401447fed29SStefano Zampini   if (rart) { PetscInt t = bm; bm = bn; bn = t; }
4025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(C,&cm,&cn));
4032c71b3e2SJacob Faibussowitsch   PetscCheckFalse(an != bm || bn != cm || bn != cn,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A, B, C local dim %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT,am,an,bm,bn,cm,cn);
404cc48ffa7SToby Isaac 
405447fed29SStefano Zampini   /* Create left vector of A: v2 */
4065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&Bx,&v2));
407447fed29SStefano Zampini 
408447fed29SStefano Zampini   /* Create right vectors of B: x, v3, v4 */
409447fed29SStefano Zampini   if (rart) {
4105f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,&v1,&x));
411447fed29SStefano Zampini   } else {
4125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(B,&x,&v1));
413447fed29SStefano Zampini   }
4145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&v3));
415447fed29SStefano Zampini 
4165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(C,&Cx,&v4));
4175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
4185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
419cc48ffa7SToby Isaac 
420cc48ffa7SToby Isaac   *flg = PETSC_TRUE;
421447fed29SStefano Zampini   for (i=0; i<n; i++) {
4225f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(x,rdm));
4235f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(x,Cx));
4245f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(C,Cx,v4));           /* v4 = C*x   */
425447fed29SStefano Zampini     if (rart) {
4265f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultTranspose(B,x,v1));
427cc48ffa7SToby Isaac     } else {
4285f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(B,x,v1));
429cc48ffa7SToby Isaac     }
4305f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(v1,Bx));
4315f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,Bx,v2));          /* v2 = A*B*x */
4325f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(v2,v1));
433447fed29SStefano Zampini     if (rart) {
4345f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(B,v1,v3)); /* v3 = R*A*R^t*x */
435447fed29SStefano Zampini     } else {
4365f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultTranspose(B,v1,v3)); /* v3 = Bt*A*B*x */
437447fed29SStefano Zampini     }
4385f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(v4,NORM_2,&norm_abs));
4395f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(v4,none,v3));
4405f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(v4,NORM_2,&norm_rel));
441447fed29SStefano Zampini 
442447fed29SStefano Zampini     if (norm_abs > tol) norm_rel /= norm_abs;
443447fed29SStefano Zampini     if (norm_rel > tol) {
444cc48ffa7SToby Isaac       *flg = PETSC_FALSE;
4455f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"Error: %" PetscInt_FMT "-th Mat%sMult() %g\n",i,rart ? "RARt" : "PtAP",(double)norm_rel));
446cc48ffa7SToby Isaac       break;
447cc48ffa7SToby Isaac     }
448cc48ffa7SToby Isaac   }
449447fed29SStefano Zampini 
4505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
4515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
4525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Bx));
4535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Cx));
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v1));
4555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v2));
4565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v3));
4575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v4));
458cc48ffa7SToby Isaac   PetscFunctionReturn(0);
459cc48ffa7SToby Isaac }
460cc48ffa7SToby Isaac 
46186919fd6SHong Zhang /*@
4624222ddf1SHong Zhang    MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B
4634222ddf1SHong Zhang 
4644222ddf1SHong Zhang    Collective on Mat
4654222ddf1SHong Zhang 
4664222ddf1SHong Zhang    Input Parameters:
4674222ddf1SHong Zhang +  A - the first matrix
4684222ddf1SHong Zhang .  B - the second matrix
4694222ddf1SHong Zhang .  C - the third matrix
4704222ddf1SHong Zhang -  n - number of random vectors to be tested
4714222ddf1SHong Zhang 
4724222ddf1SHong Zhang    Output Parameter:
4734222ddf1SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
4744222ddf1SHong Zhang 
4754222ddf1SHong Zhang    Level: intermediate
4764222ddf1SHong Zhang 
4774222ddf1SHong Zhang @*/
4784222ddf1SHong Zhang PetscErrorCode MatPtAPMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
4794222ddf1SHong Zhang {
4804222ddf1SHong Zhang   PetscFunctionBegin;
4815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProjMultEqual_Private(A,B,C,n,PETSC_FALSE,flg));
482447fed29SStefano Zampini   PetscFunctionReturn(0);
4834222ddf1SHong Zhang }
4844222ddf1SHong Zhang 
485447fed29SStefano Zampini /*@
486447fed29SStefano Zampini    MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t
487447fed29SStefano Zampini 
488447fed29SStefano Zampini    Collective on Mat
489447fed29SStefano Zampini 
490447fed29SStefano Zampini    Input Parameters:
491447fed29SStefano Zampini +  A - the first matrix
492447fed29SStefano Zampini .  B - the second matrix
493447fed29SStefano Zampini .  C - the third matrix
494447fed29SStefano Zampini -  n - number of random vectors to be tested
495447fed29SStefano Zampini 
496447fed29SStefano Zampini    Output Parameter:
497447fed29SStefano Zampini .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
498447fed29SStefano Zampini 
499447fed29SStefano Zampini    Level: intermediate
500447fed29SStefano Zampini 
501447fed29SStefano Zampini @*/
502447fed29SStefano Zampini PetscErrorCode MatRARtMultEqual(Mat A,Mat B,Mat C,PetscInt n,PetscBool *flg)
503447fed29SStefano Zampini {
504447fed29SStefano Zampini   PetscFunctionBegin;
5055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatProjMultEqual_Private(A,B,C,n,PETSC_TRUE,flg));
5064222ddf1SHong Zhang   PetscFunctionReturn(0);
5074222ddf1SHong Zhang }
5084222ddf1SHong Zhang 
5094222ddf1SHong Zhang /*@
51086919fd6SHong Zhang    MatIsLinear - Check if a shell matrix A is a linear operator.
51186919fd6SHong Zhang 
51286919fd6SHong Zhang    Collective on Mat
51386919fd6SHong Zhang 
51486919fd6SHong Zhang    Input Parameters:
51586919fd6SHong Zhang +  A - the shell matrix
51686919fd6SHong Zhang -  n - number of random vectors to be tested
51786919fd6SHong Zhang 
51886919fd6SHong Zhang    Output Parameter:
51986919fd6SHong Zhang .  flg - PETSC_TRUE if the shell matrix is linear; PETSC_FALSE otherwise.
52086919fd6SHong Zhang 
52186919fd6SHong Zhang    Level: intermediate
52286919fd6SHong Zhang @*/
52386919fd6SHong Zhang PetscErrorCode MatIsLinear(Mat A,PetscInt n,PetscBool *flg)
52486919fd6SHong Zhang {
52586919fd6SHong Zhang   Vec            x,y,s1,s2;
52686919fd6SHong Zhang   PetscRandom    rctx;
52786919fd6SHong Zhang   PetscScalar    a;
52886919fd6SHong Zhang   PetscInt       k;
52986919fd6SHong Zhang   PetscReal      norm,normA;
53086919fd6SHong Zhang   MPI_Comm       comm;
53186919fd6SHong Zhang   PetscMPIInt    rank;
53286919fd6SHong Zhang 
53386919fd6SHong Zhang   PetscFunctionBegin;
53486919fd6SHong Zhang   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
5355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)A,&comm));
5365f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
53786919fd6SHong Zhang 
5385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(comm,&rctx));
5395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rctx));
5405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,&s1));
5415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&y));
5425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(s1,&s2));
54386919fd6SHong Zhang 
54486919fd6SHong Zhang   *flg = PETSC_TRUE;
54586919fd6SHong Zhang   for (k=0; k<n; k++) {
5465f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(x,rctx));
5475f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(y,rctx));
548dd400576SPatrick Sanan     if (rank == 0) {
5495f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rctx,&a));
55086919fd6SHong Zhang     }
5515f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm));
55286919fd6SHong Zhang 
55386919fd6SHong Zhang     /* s2 = a*A*x + A*y */
5545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,y,s2)); /* s2 = A*y */
5555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,x,s1)); /* s1 = A*x */
5565f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(s2,a,s1)); /* s2 = a s1 + s2 */
55786919fd6SHong Zhang 
55886919fd6SHong Zhang     /* s1 = A * (a x + y) */
5595f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(y,a,x)); /* y = a x + y */
5605f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,y,s1));
5615f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s1,NORM_INFINITY,&normA));
56286919fd6SHong Zhang 
5635f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(s2,-1.0,s1)); /* s2 = - s1 + s2 */
5645f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(s2,NORM_INFINITY,&norm));
5651b8dac88SHong Zhang     if (norm/normA > 100.*PETSC_MACHINE_EPSILON) {
56686919fd6SHong Zhang       *flg = PETSC_FALSE;
5675f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"Error: %" PetscInt_FMT "-th |A*(ax+y) - (a*A*x+A*y)|/|A(ax+y)| %g > tol %g\n",k,(double)(norm/normA),(double)(100.*PETSC_MACHINE_EPSILON)));
56886919fd6SHong Zhang       break;
56986919fd6SHong Zhang     }
57086919fd6SHong Zhang   }
5715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rctx));
5725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
5735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
5745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s1));
5755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&s2));
57686919fd6SHong Zhang   PetscFunctionReturn(0);
57786919fd6SHong Zhang }
578