xref: /petsc/src/mat/utils/multequal.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558) !
186a22c91SHong Zhang 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
386a22c91SHong Zhang 
49371c9d4SSatish Balay static PetscErrorCode MatMultEqual_Private(Mat A, Mat B, PetscInt n, PetscBool *flg, PetscInt t, PetscBool add) {
5b84f494bSStefano Zampini   Vec         Ax = NULL, Bx = NULL, s1 = NULL, s2 = NULL, Ay = NULL, By = NULL;
6447fed29SStefano Zampini   PetscRandom rctx;
7447fed29SStefano Zampini   PetscReal   r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON;
8447fed29SStefano Zampini   PetscInt    am, an, bm, bn, k;
9447fed29SStefano Zampini   PetscScalar none   = -1.0;
10580c7c76SPierre Jolivet   const char *sops[] = {"MatMult", "MatMultAdd", "MatMultTranspose", "MatMultTransposeAdd", "MatMultHermitianTranspose", "MatMultHermitianTransposeAdd"};
11447fed29SStefano Zampini   const char *sop;
12447fed29SStefano Zampini 
13447fed29SStefano Zampini   PetscFunctionBegin;
14447fed29SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
15447fed29SStefano Zampini   PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
16447fed29SStefano Zampini   PetscCheckSameComm(A, 1, B, 2);
17447fed29SStefano Zampini   PetscValidLogicalCollectiveInt(A, n, 3);
18dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(flg, 4);
19e573050aSPierre Jolivet   PetscValidLogicalCollectiveInt(A, t, 5);
20447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(A, add, 6);
219566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
229566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(B, &bm, &bn));
23aed4548fSBarry Smith   PetscCheck(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);
24e573050aSPierre Jolivet   sop = sops[(add ? 1 : 0) + 2 * t]; /* t = 0 => no transpose, t = 1 => transpose, t = 2 => Hermitian transpose */
259566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)A), &rctx));
269566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
27447fed29SStefano Zampini   if (t) {
289566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &s1, &Ax));
299566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &s2, &Bx));
30447fed29SStefano Zampini   } else {
319566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &Ax, &s1));
329566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &Bx, &s2));
33447fed29SStefano Zampini   }
34447fed29SStefano Zampini   if (add) {
359566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(s1, &Ay));
369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(s2, &By));
37447fed29SStefano Zampini   }
38447fed29SStefano Zampini 
39447fed29SStefano Zampini   *flg = PETSC_TRUE;
40447fed29SStefano Zampini   for (k = 0; k < n; k++) {
419566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Ax, rctx));
429566063dSJacob Faibussowitsch     PetscCall(VecCopy(Ax, Bx));
43447fed29SStefano Zampini     if (add) {
449566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(Ay, rctx));
459566063dSJacob Faibussowitsch       PetscCall(VecCopy(Ay, By));
46447fed29SStefano Zampini     }
47e573050aSPierre Jolivet     if (t == 1) {
48447fed29SStefano Zampini       if (add) {
499566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(A, Ax, Ay, s1));
509566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(B, Bx, By, s2));
51447fed29SStefano Zampini       } else {
529566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(A, Ax, s1));
539566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(B, Bx, s2));
54447fed29SStefano Zampini       }
55e573050aSPierre Jolivet     } else if (t == 2) {
56e573050aSPierre Jolivet       if (add) {
579566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTransposeAdd(A, Ax, Ay, s1));
589566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTransposeAdd(B, Bx, By, s2));
59e573050aSPierre Jolivet       } else {
609566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTranspose(A, Ax, s1));
619566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTranspose(B, Bx, s2));
62e573050aSPierre Jolivet       }
63447fed29SStefano Zampini     } else {
64447fed29SStefano Zampini       if (add) {
659566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(A, Ax, Ay, s1));
669566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(B, Bx, By, s2));
67447fed29SStefano Zampini       } else {
689566063dSJacob Faibussowitsch         PetscCall(MatMult(A, Ax, s1));
699566063dSJacob Faibussowitsch         PetscCall(MatMult(B, Bx, s2));
70447fed29SStefano Zampini       }
71447fed29SStefano Zampini     }
729566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_INFINITY, &r2));
73447fed29SStefano Zampini     if (r2 < tol) {
749566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_INFINITY, &r1));
75447fed29SStefano Zampini     } else {
769566063dSJacob Faibussowitsch       PetscCall(VecAXPY(s2, none, s1));
779566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_INFINITY, &r1));
78447fed29SStefano Zampini       r1 /= r2;
79447fed29SStefano Zampini     }
80447fed29SStefano Zampini     if (r1 > tol) {
81447fed29SStefano Zampini       *flg = PETSC_FALSE;
829566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s() %g\n", k, sop, (double)r1));
83447fed29SStefano Zampini       break;
84447fed29SStefano Zampini     }
85447fed29SStefano Zampini   }
869566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Ax));
889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Bx));
899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Ay));
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&By));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
93447fed29SStefano Zampini   PetscFunctionReturn(0);
94447fed29SStefano Zampini }
95447fed29SStefano Zampini 
969371c9d4SSatish Balay static PetscErrorCode MatMatMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg, PetscBool At, PetscBool Bt) {
97447fed29SStefano Zampini   Vec         Ax, Bx, Cx, s1, s2, s3;
98447fed29SStefano Zampini   PetscRandom rctx;
99447fed29SStefano Zampini   PetscReal   r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON;
100447fed29SStefano Zampini   PetscInt    am, an, bm, bn, cm, cn, k;
101447fed29SStefano Zampini   PetscScalar none   = -1.0;
102580c7c76SPierre Jolivet   const char *sops[] = {"MatMatMult", "MatTransposeMatMult", "MatMatTransposeMult", "MatTransposeMatTransposeMult"};
103447fed29SStefano Zampini   const char *sop;
104447fed29SStefano Zampini 
105447fed29SStefano Zampini   PetscFunctionBegin;
106447fed29SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
107447fed29SStefano Zampini   PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
108447fed29SStefano Zampini   PetscCheckSameComm(A, 1, B, 2);
109447fed29SStefano Zampini   PetscValidHeaderSpecific(C, MAT_CLASSID, 3);
110447fed29SStefano Zampini   PetscCheckSameComm(A, 1, C, 3);
111447fed29SStefano Zampini   PetscValidLogicalCollectiveInt(A, n, 4);
112dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(flg, 5);
113447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(A, At, 6);
114447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(B, Bt, 7);
1159566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
1169566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(B, &bm, &bn));
1179566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &cm, &cn));
1189371c9d4SSatish Balay   if (At) {
1199371c9d4SSatish Balay     PetscInt tt = an;
1209371c9d4SSatish Balay     an          = am;
1219371c9d4SSatish Balay     am          = tt;
1229371c9d4SSatish Balay   };
1239371c9d4SSatish Balay   if (Bt) {
1249371c9d4SSatish Balay     PetscInt tt = bn;
1259371c9d4SSatish Balay     bn          = bm;
1269371c9d4SSatish Balay     bm          = tt;
1279371c9d4SSatish Balay   };
128aed4548fSBarry Smith   PetscCheck(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);
129447fed29SStefano Zampini 
130447fed29SStefano Zampini   sop = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)];
1319566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)C), &rctx));
1329566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
133447fed29SStefano Zampini   if (Bt) {
1349566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &s1, &Bx));
135447fed29SStefano Zampini   } else {
1369566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &Bx, &s1));
137447fed29SStefano Zampini   }
138447fed29SStefano Zampini   if (At) {
1399566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &s2, &Ax));
140447fed29SStefano Zampini   } else {
1419566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &Ax, &s2));
142447fed29SStefano Zampini   }
1439566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C, &Cx, &s3));
144447fed29SStefano Zampini 
145447fed29SStefano Zampini   *flg = PETSC_TRUE;
146447fed29SStefano Zampini   for (k = 0; k < n; k++) {
1479566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Bx, rctx));
148447fed29SStefano Zampini     if (Bt) {
1499566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(B, Bx, s1));
150447fed29SStefano Zampini     } else {
1519566063dSJacob Faibussowitsch       PetscCall(MatMult(B, Bx, s1));
152447fed29SStefano Zampini     }
1539566063dSJacob Faibussowitsch     PetscCall(VecCopy(s1, Ax));
154447fed29SStefano Zampini     if (At) {
1559566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A, Ax, s2));
156447fed29SStefano Zampini     } else {
1579566063dSJacob Faibussowitsch       PetscCall(MatMult(A, Ax, s2));
158447fed29SStefano Zampini     }
1599566063dSJacob Faibussowitsch     PetscCall(VecCopy(Bx, Cx));
1609566063dSJacob Faibussowitsch     PetscCall(MatMult(C, Cx, s3));
161447fed29SStefano Zampini 
1629566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_INFINITY, &r2));
163447fed29SStefano Zampini     if (r2 < tol) {
1649566063dSJacob Faibussowitsch       PetscCall(VecNorm(s3, NORM_INFINITY, &r1));
165447fed29SStefano Zampini     } else {
1669566063dSJacob Faibussowitsch       PetscCall(VecAXPY(s2, none, s3));
1679566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_INFINITY, &r1));
168447fed29SStefano Zampini       r1 /= r2;
169447fed29SStefano Zampini     }
170447fed29SStefano Zampini     if (r1 > tol) {
171447fed29SStefano Zampini       *flg = PETSC_FALSE;
1729566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s %g\n", k, sop, (double)r1));
173447fed29SStefano Zampini       break;
174447fed29SStefano Zampini     }
175447fed29SStefano Zampini   }
1769566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Ax));
1789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Bx));
1799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Cx));
1809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
1819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
1829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s3));
183447fed29SStefano Zampini   PetscFunctionReturn(0);
184447fed29SStefano Zampini }
185447fed29SStefano Zampini 
18686a22c91SHong Zhang /*@
18786a22c91SHong Zhang    MatMultEqual - Compares matrix-vector products of two matrices.
18886a22c91SHong Zhang 
18986a22c91SHong Zhang    Collective on Mat
19086a22c91SHong Zhang 
19186a22c91SHong Zhang    Input Parameters:
19286a22c91SHong Zhang +  A - the first matrix
1934222ddf1SHong Zhang .  B - the second matrix
19486a22c91SHong Zhang -  n - number of random vectors to be tested
19586a22c91SHong Zhang 
19686a22c91SHong Zhang    Output Parameter:
19786a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
19886a22c91SHong Zhang 
19986a22c91SHong Zhang    Level: intermediate
20086a22c91SHong Zhang 
20186a22c91SHong Zhang @*/
2029371c9d4SSatish Balay PetscErrorCode MatMultEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) {
20386a22c91SHong Zhang   PetscFunctionBegin;
2049566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_FALSE));
20586a22c91SHong Zhang   PetscFunctionReturn(0);
20686a22c91SHong Zhang }
20786a22c91SHong Zhang 
20886a22c91SHong Zhang /*@
20986a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
21086a22c91SHong Zhang 
21186a22c91SHong Zhang    Collective on Mat
21286a22c91SHong Zhang 
21386a22c91SHong Zhang    Input Parameters:
21486a22c91SHong Zhang +  A - the first matrix
2154222ddf1SHong Zhang .  B - the second matrix
21686a22c91SHong Zhang -  n - number of random vectors to be tested
21786a22c91SHong Zhang 
21886a22c91SHong Zhang    Output Parameter:
21986a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
22086a22c91SHong Zhang 
22186a22c91SHong Zhang    Level: intermediate
22286a22c91SHong Zhang 
22386a22c91SHong Zhang @*/
2249371c9d4SSatish Balay PetscErrorCode MatMultAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) {
22586a22c91SHong Zhang   PetscFunctionBegin;
2269566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_TRUE));
22763879571SHong Zhang   PetscFunctionReturn(0);
22863879571SHong Zhang }
22963879571SHong Zhang 
23063879571SHong Zhang /*@
23163879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
23263879571SHong Zhang 
23363879571SHong Zhang    Collective on Mat
23463879571SHong Zhang 
23563879571SHong Zhang    Input Parameters:
23663879571SHong Zhang +  A - the first matrix
2374222ddf1SHong Zhang .  B - the second matrix
23863879571SHong Zhang -  n - number of random vectors to be tested
23963879571SHong Zhang 
24063879571SHong Zhang    Output Parameter:
24163879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
24263879571SHong Zhang 
24363879571SHong Zhang    Level: intermediate
24463879571SHong Zhang 
24563879571SHong Zhang @*/
2469371c9d4SSatish Balay PetscErrorCode MatMultTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) {
24763879571SHong Zhang   PetscFunctionBegin;
2489566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_FALSE));
24963879571SHong Zhang   PetscFunctionReturn(0);
25063879571SHong Zhang }
25163879571SHong Zhang 
25263879571SHong Zhang /*@
25363879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
25463879571SHong Zhang 
25563879571SHong Zhang    Collective on Mat
25663879571SHong Zhang 
25763879571SHong Zhang    Input Parameters:
25863879571SHong Zhang +  A - the first matrix
2594222ddf1SHong Zhang .  B - the second matrix
26063879571SHong Zhang -  n - number of random vectors to be tested
26163879571SHong Zhang 
26263879571SHong Zhang    Output Parameter:
26363879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
26463879571SHong Zhang 
26563879571SHong Zhang    Level: intermediate
26663879571SHong Zhang 
26763879571SHong Zhang @*/
2689371c9d4SSatish Balay PetscErrorCode MatMultTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) {
26963879571SHong Zhang   PetscFunctionBegin;
2709566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_TRUE));
271e573050aSPierre Jolivet   PetscFunctionReturn(0);
272e573050aSPierre Jolivet }
273e573050aSPierre Jolivet 
274e573050aSPierre Jolivet /*@
275e573050aSPierre Jolivet    MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices.
276e573050aSPierre Jolivet 
277e573050aSPierre Jolivet    Collective on Mat
278e573050aSPierre Jolivet 
279e573050aSPierre Jolivet    Input Parameters:
280e573050aSPierre Jolivet +  A - the first matrix
281e573050aSPierre Jolivet .  B - the second matrix
282e573050aSPierre Jolivet -  n - number of random vectors to be tested
283e573050aSPierre Jolivet 
284e573050aSPierre Jolivet    Output Parameter:
285e573050aSPierre Jolivet .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
286e573050aSPierre Jolivet 
287e573050aSPierre Jolivet    Level: intermediate
288e573050aSPierre Jolivet 
289e573050aSPierre Jolivet @*/
2909371c9d4SSatish Balay PetscErrorCode MatMultHermitianTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) {
291e573050aSPierre Jolivet   PetscFunctionBegin;
2929566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 2, PETSC_FALSE));
293e573050aSPierre Jolivet   PetscFunctionReturn(0);
294e573050aSPierre Jolivet }
295e573050aSPierre Jolivet 
296e573050aSPierre Jolivet /*@
297e573050aSPierre Jolivet    MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices.
298e573050aSPierre Jolivet 
299e573050aSPierre Jolivet    Collective on Mat
300e573050aSPierre Jolivet 
301e573050aSPierre Jolivet    Input Parameters:
302e573050aSPierre Jolivet +  A - the first matrix
303e573050aSPierre Jolivet .  B - the second matrix
304e573050aSPierre Jolivet -  n - number of random vectors to be tested
305e573050aSPierre Jolivet 
306e573050aSPierre Jolivet    Output Parameter:
307e573050aSPierre Jolivet .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
308e573050aSPierre Jolivet 
309e573050aSPierre Jolivet    Level: intermediate
310e573050aSPierre Jolivet 
311e573050aSPierre Jolivet @*/
3129371c9d4SSatish Balay PetscErrorCode MatMultHermitianTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) {
313e573050aSPierre Jolivet   PetscFunctionBegin;
3149566063dSJacob Faibussowitsch   PetscCall(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 @*/
3359371c9d4SSatish Balay PetscErrorCode MatMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) {
336a52676ecSHong Zhang   PetscFunctionBegin;
3379566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_FALSE));
338a52676ecSHong Zhang   PetscFunctionReturn(0);
339a52676ecSHong Zhang }
340a52676ecSHong Zhang 
341a52676ecSHong Zhang /*@
342a52676ecSHong Zhang    MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x
343a52676ecSHong Zhang 
344a52676ecSHong Zhang    Collective on Mat
345a52676ecSHong Zhang 
346a52676ecSHong Zhang    Input Parameters:
347a52676ecSHong Zhang +  A - the first matrix
3484222ddf1SHong Zhang .  B - the second matrix
3494222ddf1SHong Zhang .  C - the third matrix
350a52676ecSHong Zhang -  n - number of random vectors to be tested
351a52676ecSHong Zhang 
352a52676ecSHong Zhang    Output Parameter:
353a52676ecSHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
354a52676ecSHong Zhang 
355a52676ecSHong Zhang    Level: intermediate
356a52676ecSHong Zhang 
357a52676ecSHong Zhang @*/
3589371c9d4SSatish Balay PetscErrorCode MatTransposeMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) {
359a52676ecSHong Zhang   PetscFunctionBegin;
3609566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_TRUE, PETSC_FALSE));
361a52676ecSHong Zhang   PetscFunctionReturn(0);
362a52676ecSHong Zhang }
36386919fd6SHong Zhang 
36426546c1bSToby Isaac /*@
36526546c1bSToby Isaac    MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x
36626546c1bSToby Isaac 
36726546c1bSToby Isaac    Collective on Mat
36826546c1bSToby Isaac 
36926546c1bSToby Isaac    Input Parameters:
37026546c1bSToby Isaac +  A - the first matrix
3714222ddf1SHong Zhang .  B - the second matrix
3724222ddf1SHong Zhang .  C - the third matrix
37326546c1bSToby Isaac -  n - number of random vectors to be tested
37426546c1bSToby Isaac 
37526546c1bSToby Isaac    Output Parameter:
37626546c1bSToby Isaac .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
37726546c1bSToby Isaac 
37826546c1bSToby Isaac    Level: intermediate
37926546c1bSToby Isaac 
38026546c1bSToby Isaac @*/
3819371c9d4SSatish Balay PetscErrorCode MatMatTransposeMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) {
382447fed29SStefano Zampini   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_TRUE));
384447fed29SStefano Zampini   PetscFunctionReturn(0);
385447fed29SStefano Zampini }
386447fed29SStefano Zampini 
3879371c9d4SSatish Balay static PetscErrorCode MatProjMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool rart, PetscBool *flg) {
388447fed29SStefano Zampini   Vec         x, v1, v2, v3, v4, Cx, Bx;
389447fed29SStefano Zampini   PetscReal   norm_abs, norm_rel, tol = PETSC_SQRT_MACHINE_EPSILON;
390447fed29SStefano Zampini   PetscInt    i, am, an, bm, bn, cm, cn;
391447fed29SStefano Zampini   PetscRandom rdm;
392cc48ffa7SToby Isaac   PetscScalar none = -1.0;
393cc48ffa7SToby Isaac 
394cc48ffa7SToby Isaac   PetscFunctionBegin;
3959566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
3969566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(B, &bm, &bn));
3979371c9d4SSatish Balay   if (rart) {
3989371c9d4SSatish Balay     PetscInt t = bm;
3999371c9d4SSatish Balay     bm         = bn;
4009371c9d4SSatish Balay     bn         = t;
4019371c9d4SSatish Balay   }
4029566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &cm, &cn));
403aed4548fSBarry Smith   PetscCheck(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 */
4069566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &Bx, &v2));
407447fed29SStefano Zampini 
408447fed29SStefano Zampini   /* Create right vectors of B: x, v3, v4 */
409447fed29SStefano Zampini   if (rart) {
4109566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &v1, &x));
411447fed29SStefano Zampini   } else {
4129566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &x, &v1));
413447fed29SStefano Zampini   }
4149566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &v3));
415447fed29SStefano Zampini 
4169566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C, &Cx, &v4));
4179566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
4189566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
419cc48ffa7SToby Isaac 
420cc48ffa7SToby Isaac   *flg = PETSC_TRUE;
421447fed29SStefano Zampini   for (i = 0; i < n; i++) {
4229566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
4239566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, Cx));
4249566063dSJacob Faibussowitsch     PetscCall(MatMult(C, Cx, v4)); /* v4 = C*x   */
425447fed29SStefano Zampini     if (rart) {
4269566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(B, x, v1));
427cc48ffa7SToby Isaac     } else {
4289566063dSJacob Faibussowitsch       PetscCall(MatMult(B, x, v1));
429cc48ffa7SToby Isaac     }
4309566063dSJacob Faibussowitsch     PetscCall(VecCopy(v1, Bx));
4319566063dSJacob Faibussowitsch     PetscCall(MatMult(A, Bx, v2)); /* v2 = A*B*x */
4329566063dSJacob Faibussowitsch     PetscCall(VecCopy(v2, v1));
433447fed29SStefano Zampini     if (rart) {
4349566063dSJacob Faibussowitsch       PetscCall(MatMult(B, v1, v3)); /* v3 = R*A*R^t*x */
435447fed29SStefano Zampini     } else {
4369566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(B, v1, v3)); /* v3 = Bt*A*B*x */
437447fed29SStefano Zampini     }
4389566063dSJacob Faibussowitsch     PetscCall(VecNorm(v4, NORM_2, &norm_abs));
4399566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v4, none, v3));
4409566063dSJacob Faibussowitsch     PetscCall(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;
4459566063dSJacob Faibussowitsch       PetscCall(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 
4509566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
4519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
4529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Bx));
4539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Cx));
4549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v1));
4559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v2));
4569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v3));
4579566063dSJacob Faibussowitsch   PetscCall(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 @*/
4789371c9d4SSatish Balay PetscErrorCode MatPtAPMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) {
4794222ddf1SHong Zhang   PetscFunctionBegin;
4809566063dSJacob Faibussowitsch   PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_FALSE, flg));
481447fed29SStefano Zampini   PetscFunctionReturn(0);
4824222ddf1SHong Zhang }
4834222ddf1SHong Zhang 
484447fed29SStefano Zampini /*@
485447fed29SStefano Zampini    MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t
486447fed29SStefano Zampini 
487447fed29SStefano Zampini    Collective on Mat
488447fed29SStefano Zampini 
489447fed29SStefano Zampini    Input Parameters:
490447fed29SStefano Zampini +  A - the first matrix
491447fed29SStefano Zampini .  B - the second matrix
492447fed29SStefano Zampini .  C - the third matrix
493447fed29SStefano Zampini -  n - number of random vectors to be tested
494447fed29SStefano Zampini 
495447fed29SStefano Zampini    Output Parameter:
496447fed29SStefano Zampini .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
497447fed29SStefano Zampini 
498447fed29SStefano Zampini    Level: intermediate
499447fed29SStefano Zampini 
500447fed29SStefano Zampini @*/
5019371c9d4SSatish Balay PetscErrorCode MatRARtMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) {
502447fed29SStefano Zampini   PetscFunctionBegin;
5039566063dSJacob Faibussowitsch   PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_TRUE, flg));
5044222ddf1SHong Zhang   PetscFunctionReturn(0);
5054222ddf1SHong Zhang }
5064222ddf1SHong Zhang 
5074222ddf1SHong Zhang /*@
50886919fd6SHong Zhang    MatIsLinear - Check if a shell matrix A is a linear operator.
50986919fd6SHong Zhang 
51086919fd6SHong Zhang    Collective on Mat
51186919fd6SHong Zhang 
51286919fd6SHong Zhang    Input Parameters:
51386919fd6SHong Zhang +  A - the shell matrix
51486919fd6SHong Zhang -  n - number of random vectors to be tested
51586919fd6SHong Zhang 
51686919fd6SHong Zhang    Output Parameter:
51786919fd6SHong Zhang .  flg - PETSC_TRUE if the shell matrix is linear; PETSC_FALSE otherwise.
51886919fd6SHong Zhang 
51986919fd6SHong Zhang    Level: intermediate
52086919fd6SHong Zhang @*/
5219371c9d4SSatish Balay PetscErrorCode MatIsLinear(Mat A, PetscInt n, PetscBool *flg) {
52286919fd6SHong Zhang   Vec         x, y, s1, s2;
52386919fd6SHong Zhang   PetscRandom rctx;
52486919fd6SHong Zhang   PetscScalar a;
52586919fd6SHong Zhang   PetscInt    k;
52686919fd6SHong Zhang   PetscReal   norm, normA;
52786919fd6SHong Zhang   MPI_Comm    comm;
52886919fd6SHong Zhang   PetscMPIInt rank;
52986919fd6SHong Zhang 
53086919fd6SHong Zhang   PetscFunctionBegin;
53186919fd6SHong Zhang   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
5329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
53486919fd6SHong Zhang 
5359566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(comm, &rctx));
5369566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
5379566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &s1));
5389566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
5399566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s1, &s2));
54086919fd6SHong Zhang 
54186919fd6SHong Zhang   *flg = PETSC_TRUE;
54286919fd6SHong Zhang   for (k = 0; k < n; k++) {
5439566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rctx));
5449566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y, rctx));
545*48a46eb9SPierre Jolivet     if (rank == 0) PetscCall(PetscRandomGetValue(rctx, &a));
5469566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm));
54786919fd6SHong Zhang 
54886919fd6SHong Zhang     /* s2 = a*A*x + A*y */
5499566063dSJacob Faibussowitsch     PetscCall(MatMult(A, y, s2));  /* s2 = A*y */
5509566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, s1));  /* s1 = A*x */
5519566063dSJacob Faibussowitsch     PetscCall(VecAXPY(s2, a, s1)); /* s2 = a s1 + s2 */
55286919fd6SHong Zhang 
55386919fd6SHong Zhang     /* s1 = A * (a x + y) */
5549566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, a, x)); /* y = a x + y */
5559566063dSJacob Faibussowitsch     PetscCall(MatMult(A, y, s1));
5569566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_INFINITY, &normA));
55786919fd6SHong Zhang 
5589566063dSJacob Faibussowitsch     PetscCall(VecAXPY(s2, -1.0, s1)); /* s2 = - s1 + s2 */
5599566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_INFINITY, &norm));
5601b8dac88SHong Zhang     if (norm / normA > 100. * PETSC_MACHINE_EPSILON) {
56186919fd6SHong Zhang       *flg = PETSC_FALSE;
5629566063dSJacob Faibussowitsch       PetscCall(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)));
56386919fd6SHong Zhang       break;
56486919fd6SHong Zhang     }
56586919fd6SHong Zhang   }
5669566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
5679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
5689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
5699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
5709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
57186919fd6SHong Zhang   PetscFunctionReturn(0);
57286919fd6SHong Zhang }
573