xref: /petsc/src/mat/utils/multequal.c (revision 4279555e62e618147958e52682840638ab8c8be0)
186a22c91SHong Zhang 
2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I   "petscmat.h"  I*/
386a22c91SHong Zhang 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultEqual_Private(Mat A, Mat B, PetscInt n, PetscBool *flg, PetscInt t, PetscBool add)
5d71ae5a4SJacob Faibussowitsch {
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;
11*4279555eSSatish Balay #if defined(PETSC_USE_INFO)
12580c7c76SPierre Jolivet   const char *sops[] = {"MatMult", "MatMultAdd", "MatMultTranspose", "MatMultTransposeAdd", "MatMultHermitianTranspose", "MatMultHermitianTransposeAdd"};
13447fed29SStefano Zampini   const char *sop;
14*4279555eSSatish Balay #endif
15447fed29SStefano Zampini 
16447fed29SStefano Zampini   PetscFunctionBegin;
17447fed29SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
18447fed29SStefano Zampini   PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
19447fed29SStefano Zampini   PetscCheckSameComm(A, 1, B, 2);
20447fed29SStefano Zampini   PetscValidLogicalCollectiveInt(A, n, 3);
21dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(flg, 4);
22e573050aSPierre Jolivet   PetscValidLogicalCollectiveInt(A, t, 5);
23447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(A, add, 6);
249566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(B, &bm, &bn));
26aed4548fSBarry 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);
27*4279555eSSatish Balay #if defined(PETSC_USE_INFO)
28e573050aSPierre Jolivet   sop = sops[(add ? 1 : 0) + 2 * t]; /* t = 0 => no transpose, t = 1 => transpose, t = 2 => Hermitian transpose */
29*4279555eSSatish Balay #endif
309566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)A), &rctx));
319566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
32447fed29SStefano Zampini   if (t) {
339566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &s1, &Ax));
349566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &s2, &Bx));
35447fed29SStefano Zampini   } else {
369566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &Ax, &s1));
379566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &Bx, &s2));
38447fed29SStefano Zampini   }
39447fed29SStefano Zampini   if (add) {
409566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(s1, &Ay));
419566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(s2, &By));
42447fed29SStefano Zampini   }
43447fed29SStefano Zampini 
44447fed29SStefano Zampini   *flg = PETSC_TRUE;
45447fed29SStefano Zampini   for (k = 0; k < n; k++) {
469566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Ax, rctx));
479566063dSJacob Faibussowitsch     PetscCall(VecCopy(Ax, Bx));
48447fed29SStefano Zampini     if (add) {
499566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(Ay, rctx));
509566063dSJacob Faibussowitsch       PetscCall(VecCopy(Ay, By));
51447fed29SStefano Zampini     }
52e573050aSPierre Jolivet     if (t == 1) {
53447fed29SStefano Zampini       if (add) {
549566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(A, Ax, Ay, s1));
559566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(B, Bx, By, s2));
56447fed29SStefano Zampini       } else {
579566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(A, Ax, s1));
589566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(B, Bx, s2));
59447fed29SStefano Zampini       }
60e573050aSPierre Jolivet     } else if (t == 2) {
61e573050aSPierre Jolivet       if (add) {
629566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTransposeAdd(A, Ax, Ay, s1));
639566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTransposeAdd(B, Bx, By, s2));
64e573050aSPierre Jolivet       } else {
659566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTranspose(A, Ax, s1));
669566063dSJacob Faibussowitsch         PetscCall(MatMultHermitianTranspose(B, Bx, s2));
67e573050aSPierre Jolivet       }
68447fed29SStefano Zampini     } else {
69447fed29SStefano Zampini       if (add) {
709566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(A, Ax, Ay, s1));
719566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(B, Bx, By, s2));
72447fed29SStefano Zampini       } else {
739566063dSJacob Faibussowitsch         PetscCall(MatMult(A, Ax, s1));
749566063dSJacob Faibussowitsch         PetscCall(MatMult(B, Bx, s2));
75447fed29SStefano Zampini       }
76447fed29SStefano Zampini     }
779566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_INFINITY, &r2));
78447fed29SStefano Zampini     if (r2 < tol) {
799566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_INFINITY, &r1));
80447fed29SStefano Zampini     } else {
819566063dSJacob Faibussowitsch       PetscCall(VecAXPY(s2, none, s1));
829566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_INFINITY, &r1));
83447fed29SStefano Zampini       r1 /= r2;
84447fed29SStefano Zampini     }
85447fed29SStefano Zampini     if (r1 > tol) {
86447fed29SStefano Zampini       *flg = PETSC_FALSE;
879566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s() %g\n", k, sop, (double)r1));
88447fed29SStefano Zampini       break;
89447fed29SStefano Zampini     }
90447fed29SStefano Zampini   }
919566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Ax));
939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Bx));
949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Ay));
959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&By));
969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99447fed29SStefano Zampini }
100447fed29SStefano Zampini 
101d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg, PetscBool At, PetscBool Bt)
102d71ae5a4SJacob Faibussowitsch {
103447fed29SStefano Zampini   Vec         Ax, Bx, Cx, s1, s2, s3;
104447fed29SStefano Zampini   PetscRandom rctx;
105447fed29SStefano Zampini   PetscReal   r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON;
106447fed29SStefano Zampini   PetscInt    am, an, bm, bn, cm, cn, k;
107447fed29SStefano Zampini   PetscScalar none = -1.0;
108*4279555eSSatish Balay #if defined(PETSC_USE_INFO)
109580c7c76SPierre Jolivet   const char *sops[] = {"MatMatMult", "MatTransposeMatMult", "MatMatTransposeMult", "MatTransposeMatTransposeMult"};
110447fed29SStefano Zampini   const char *sop;
111*4279555eSSatish Balay #endif
112447fed29SStefano Zampini 
113447fed29SStefano Zampini   PetscFunctionBegin;
114447fed29SStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
115447fed29SStefano Zampini   PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
116447fed29SStefano Zampini   PetscCheckSameComm(A, 1, B, 2);
117447fed29SStefano Zampini   PetscValidHeaderSpecific(C, MAT_CLASSID, 3);
118447fed29SStefano Zampini   PetscCheckSameComm(A, 1, C, 3);
119447fed29SStefano Zampini   PetscValidLogicalCollectiveInt(A, n, 4);
120dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(flg, 5);
121447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(A, At, 6);
122447fed29SStefano Zampini   PetscValidLogicalCollectiveBool(B, Bt, 7);
1239566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
1249566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(B, &bm, &bn));
1259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &cm, &cn));
1269371c9d4SSatish Balay   if (At) {
1279371c9d4SSatish Balay     PetscInt tt = an;
1289371c9d4SSatish Balay     an          = am;
1299371c9d4SSatish Balay     am          = tt;
1309371c9d4SSatish Balay   };
1319371c9d4SSatish Balay   if (Bt) {
1329371c9d4SSatish Balay     PetscInt tt = bn;
1339371c9d4SSatish Balay     bn          = bm;
1349371c9d4SSatish Balay     bm          = tt;
1359371c9d4SSatish Balay   };
136aed4548fSBarry 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);
137447fed29SStefano Zampini 
138*4279555eSSatish Balay #if defined(PETSC_USE_INFO)
139447fed29SStefano Zampini   sop = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)];
140*4279555eSSatish Balay #endif
1419566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)C), &rctx));
1429566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
143447fed29SStefano Zampini   if (Bt) {
1449566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &s1, &Bx));
145447fed29SStefano Zampini   } else {
1469566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &Bx, &s1));
147447fed29SStefano Zampini   }
148447fed29SStefano Zampini   if (At) {
1499566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &s2, &Ax));
150447fed29SStefano Zampini   } else {
1519566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &Ax, &s2));
152447fed29SStefano Zampini   }
1539566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C, &Cx, &s3));
154447fed29SStefano Zampini 
155447fed29SStefano Zampini   *flg = PETSC_TRUE;
156447fed29SStefano Zampini   for (k = 0; k < n; k++) {
1579566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(Bx, rctx));
158447fed29SStefano Zampini     if (Bt) {
1599566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(B, Bx, s1));
160447fed29SStefano Zampini     } else {
1619566063dSJacob Faibussowitsch       PetscCall(MatMult(B, Bx, s1));
162447fed29SStefano Zampini     }
1639566063dSJacob Faibussowitsch     PetscCall(VecCopy(s1, Ax));
164447fed29SStefano Zampini     if (At) {
1659566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A, Ax, s2));
166447fed29SStefano Zampini     } else {
1679566063dSJacob Faibussowitsch       PetscCall(MatMult(A, Ax, s2));
168447fed29SStefano Zampini     }
1699566063dSJacob Faibussowitsch     PetscCall(VecCopy(Bx, Cx));
1709566063dSJacob Faibussowitsch     PetscCall(MatMult(C, Cx, s3));
171447fed29SStefano Zampini 
1729566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_INFINITY, &r2));
173447fed29SStefano Zampini     if (r2 < tol) {
1749566063dSJacob Faibussowitsch       PetscCall(VecNorm(s3, NORM_INFINITY, &r1));
175447fed29SStefano Zampini     } else {
1769566063dSJacob Faibussowitsch       PetscCall(VecAXPY(s2, none, s3));
1779566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_INFINITY, &r1));
178447fed29SStefano Zampini       r1 /= r2;
179447fed29SStefano Zampini     }
180447fed29SStefano Zampini     if (r1 > tol) {
181447fed29SStefano Zampini       *flg = PETSC_FALSE;
1829566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s %g\n", k, sop, (double)r1));
183447fed29SStefano Zampini       break;
184447fed29SStefano Zampini     }
185447fed29SStefano Zampini   }
1869566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
1879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Ax));
1889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Bx));
1899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Cx));
1909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
1919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
1929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s3));
1933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
194447fed29SStefano Zampini }
195447fed29SStefano Zampini 
19686a22c91SHong Zhang /*@
19786a22c91SHong Zhang    MatMultEqual - Compares matrix-vector products of two matrices.
19886a22c91SHong Zhang 
199c3339decSBarry Smith    Collective
20086a22c91SHong Zhang 
20186a22c91SHong Zhang    Input Parameters:
20286a22c91SHong Zhang +  A - the first matrix
2034222ddf1SHong Zhang .  B - the second matrix
20486a22c91SHong Zhang -  n - number of random vectors to be tested
20586a22c91SHong Zhang 
20686a22c91SHong Zhang    Output Parameter:
20786a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
20886a22c91SHong Zhang 
20986a22c91SHong Zhang    Level: intermediate
21086a22c91SHong Zhang 
21186a22c91SHong Zhang @*/
212d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultEqual(Mat A, Mat B, PetscInt n, PetscBool *flg)
213d71ae5a4SJacob Faibussowitsch {
21486a22c91SHong Zhang   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_FALSE));
2163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21786a22c91SHong Zhang }
21886a22c91SHong Zhang 
21986a22c91SHong Zhang /*@
22086a22c91SHong Zhang    MatMultAddEqual - Compares matrix-vector products of two matrices.
22186a22c91SHong Zhang 
222c3339decSBarry Smith    Collective
22386a22c91SHong Zhang 
22486a22c91SHong Zhang    Input Parameters:
22586a22c91SHong Zhang +  A - the first matrix
2264222ddf1SHong Zhang .  B - the second matrix
22786a22c91SHong Zhang -  n - number of random vectors to be tested
22886a22c91SHong Zhang 
22986a22c91SHong Zhang    Output Parameter:
23086a22c91SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
23186a22c91SHong Zhang 
23286a22c91SHong Zhang    Level: intermediate
23386a22c91SHong Zhang 
23486a22c91SHong Zhang @*/
235d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg)
236d71ae5a4SJacob Faibussowitsch {
23786a22c91SHong Zhang   PetscFunctionBegin;
2389566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_TRUE));
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24063879571SHong Zhang }
24163879571SHong Zhang 
24263879571SHong Zhang /*@
24363879571SHong Zhang    MatMultTransposeEqual - Compares matrix-vector products of two matrices.
24463879571SHong Zhang 
245c3339decSBarry Smith    Collective
24663879571SHong Zhang 
24763879571SHong Zhang    Input Parameters:
24863879571SHong Zhang +  A - the first matrix
2494222ddf1SHong Zhang .  B - the second matrix
25063879571SHong Zhang -  n - number of random vectors to be tested
25163879571SHong Zhang 
25263879571SHong Zhang    Output Parameter:
25363879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
25463879571SHong Zhang 
25563879571SHong Zhang    Level: intermediate
25663879571SHong Zhang 
25763879571SHong Zhang @*/
258d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg)
259d71ae5a4SJacob Faibussowitsch {
26063879571SHong Zhang   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_FALSE));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26363879571SHong Zhang }
26463879571SHong Zhang 
26563879571SHong Zhang /*@
26663879571SHong Zhang    MatMultTransposeAddEqual - Compares matrix-vector products of two matrices.
26763879571SHong Zhang 
268c3339decSBarry Smith    Collective
26963879571SHong Zhang 
27063879571SHong Zhang    Input Parameters:
27163879571SHong Zhang +  A - the first matrix
2724222ddf1SHong Zhang .  B - the second matrix
27363879571SHong Zhang -  n - number of random vectors to be tested
27463879571SHong Zhang 
27563879571SHong Zhang    Output Parameter:
27663879571SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
27763879571SHong Zhang 
27863879571SHong Zhang    Level: intermediate
27963879571SHong Zhang 
28063879571SHong Zhang @*/
281d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg)
282d71ae5a4SJacob Faibussowitsch {
28363879571SHong Zhang   PetscFunctionBegin;
2849566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_TRUE));
2853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
286e573050aSPierre Jolivet }
287e573050aSPierre Jolivet 
288e573050aSPierre Jolivet /*@
289e573050aSPierre Jolivet    MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices.
290e573050aSPierre Jolivet 
291c3339decSBarry Smith    Collective
292e573050aSPierre Jolivet 
293e573050aSPierre Jolivet    Input Parameters:
294e573050aSPierre Jolivet +  A - the first matrix
295e573050aSPierre Jolivet .  B - the second matrix
296e573050aSPierre Jolivet -  n - number of random vectors to be tested
297e573050aSPierre Jolivet 
298e573050aSPierre Jolivet    Output Parameter:
299e573050aSPierre Jolivet .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
300e573050aSPierre Jolivet 
301e573050aSPierre Jolivet    Level: intermediate
302e573050aSPierre Jolivet 
303e573050aSPierre Jolivet @*/
304d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg)
305d71ae5a4SJacob Faibussowitsch {
306e573050aSPierre Jolivet   PetscFunctionBegin;
3079566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 2, PETSC_FALSE));
3083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
309e573050aSPierre Jolivet }
310e573050aSPierre Jolivet 
311e573050aSPierre Jolivet /*@
312e573050aSPierre Jolivet    MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices.
313e573050aSPierre Jolivet 
314c3339decSBarry Smith    Collective
315e573050aSPierre Jolivet 
316e573050aSPierre Jolivet    Input Parameters:
317e573050aSPierre Jolivet +  A - the first matrix
318e573050aSPierre Jolivet .  B - the second matrix
319e573050aSPierre Jolivet -  n - number of random vectors to be tested
320e573050aSPierre Jolivet 
321e573050aSPierre Jolivet    Output Parameter:
322e573050aSPierre Jolivet .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
323e573050aSPierre Jolivet 
324e573050aSPierre Jolivet    Level: intermediate
325e573050aSPierre Jolivet 
326e573050aSPierre Jolivet @*/
327d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg)
328d71ae5a4SJacob Faibussowitsch {
329e573050aSPierre Jolivet   PetscFunctionBegin;
3309566063dSJacob Faibussowitsch   PetscCall(MatMultEqual_Private(A, B, n, flg, 2, PETSC_TRUE));
3313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33286a22c91SHong Zhang }
333a52676ecSHong Zhang 
334a52676ecSHong Zhang /*@
335a52676ecSHong Zhang    MatMatMultEqual - Test A*B*x = C*x for n random vector x
336a52676ecSHong Zhang 
337c3339decSBarry Smith    Collective
338a52676ecSHong Zhang 
339a52676ecSHong Zhang    Input Parameters:
340a52676ecSHong Zhang +  A - the first matrix
341c2916339SPierre Jolivet .  B - the second matrix
342c2916339SPierre Jolivet .  C - the third matrix
343a52676ecSHong Zhang -  n - number of random vectors to be tested
344a52676ecSHong Zhang 
345a52676ecSHong Zhang    Output Parameter:
346f0fc11ceSJed Brown .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
347a52676ecSHong Zhang 
348a52676ecSHong Zhang    Level: intermediate
349a52676ecSHong Zhang 
350a52676ecSHong Zhang @*/
351d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg)
352d71ae5a4SJacob Faibussowitsch {
353a52676ecSHong Zhang   PetscFunctionBegin;
3549566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_FALSE));
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
356a52676ecSHong Zhang }
357a52676ecSHong Zhang 
358a52676ecSHong Zhang /*@
359a52676ecSHong Zhang    MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x
360a52676ecSHong Zhang 
361c3339decSBarry Smith    Collective
362a52676ecSHong Zhang 
363a52676ecSHong Zhang    Input Parameters:
364a52676ecSHong Zhang +  A - the first matrix
3654222ddf1SHong Zhang .  B - the second matrix
3664222ddf1SHong Zhang .  C - the third matrix
367a52676ecSHong Zhang -  n - number of random vectors to be tested
368a52676ecSHong Zhang 
369a52676ecSHong Zhang    Output Parameter:
370a52676ecSHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
371a52676ecSHong Zhang 
372a52676ecSHong Zhang    Level: intermediate
373a52676ecSHong Zhang 
374a52676ecSHong Zhang @*/
375d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTransposeMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg)
376d71ae5a4SJacob Faibussowitsch {
377a52676ecSHong Zhang   PetscFunctionBegin;
3789566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_TRUE, PETSC_FALSE));
3793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
380a52676ecSHong Zhang }
38186919fd6SHong Zhang 
38226546c1bSToby Isaac /*@
38326546c1bSToby Isaac    MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x
38426546c1bSToby Isaac 
385c3339decSBarry Smith    Collective
38626546c1bSToby Isaac 
38726546c1bSToby Isaac    Input Parameters:
38826546c1bSToby Isaac +  A - the first matrix
3894222ddf1SHong Zhang .  B - the second matrix
3904222ddf1SHong Zhang .  C - the third matrix
39126546c1bSToby Isaac -  n - number of random vectors to be tested
39226546c1bSToby Isaac 
39326546c1bSToby Isaac    Output Parameter:
39426546c1bSToby Isaac .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
39526546c1bSToby Isaac 
39626546c1bSToby Isaac    Level: intermediate
39726546c1bSToby Isaac 
39826546c1bSToby Isaac @*/
399d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatTransposeMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg)
400d71ae5a4SJacob Faibussowitsch {
401447fed29SStefano Zampini   PetscFunctionBegin;
4029566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_TRUE));
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
404447fed29SStefano Zampini }
405447fed29SStefano Zampini 
406d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProjMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool rart, PetscBool *flg)
407d71ae5a4SJacob Faibussowitsch {
408447fed29SStefano Zampini   Vec         x, v1, v2, v3, v4, Cx, Bx;
409447fed29SStefano Zampini   PetscReal   norm_abs, norm_rel, tol = PETSC_SQRT_MACHINE_EPSILON;
410447fed29SStefano Zampini   PetscInt    i, am, an, bm, bn, cm, cn;
411447fed29SStefano Zampini   PetscRandom rdm;
412cc48ffa7SToby Isaac   PetscScalar none = -1.0;
413cc48ffa7SToby Isaac 
414cc48ffa7SToby Isaac   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
4169566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(B, &bm, &bn));
4179371c9d4SSatish Balay   if (rart) {
4189371c9d4SSatish Balay     PetscInt t = bm;
4199371c9d4SSatish Balay     bm         = bn;
4209371c9d4SSatish Balay     bn         = t;
4219371c9d4SSatish Balay   }
4229566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &cm, &cn));
423aed4548fSBarry 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);
424cc48ffa7SToby Isaac 
425447fed29SStefano Zampini   /* Create left vector of A: v2 */
4269566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &Bx, &v2));
427447fed29SStefano Zampini 
428447fed29SStefano Zampini   /* Create right vectors of B: x, v3, v4 */
429447fed29SStefano Zampini   if (rart) {
4309566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &v1, &x));
431447fed29SStefano Zampini   } else {
4329566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(B, &x, &v1));
433447fed29SStefano Zampini   }
4349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &v3));
435447fed29SStefano Zampini 
4369566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C, &Cx, &v4));
4379566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
4389566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
439cc48ffa7SToby Isaac 
440cc48ffa7SToby Isaac   *flg = PETSC_TRUE;
441447fed29SStefano Zampini   for (i = 0; i < n; i++) {
4429566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
4439566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, Cx));
4449566063dSJacob Faibussowitsch     PetscCall(MatMult(C, Cx, v4)); /* v4 = C*x   */
445447fed29SStefano Zampini     if (rart) {
4469566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(B, x, v1));
447cc48ffa7SToby Isaac     } else {
4489566063dSJacob Faibussowitsch       PetscCall(MatMult(B, x, v1));
449cc48ffa7SToby Isaac     }
4509566063dSJacob Faibussowitsch     PetscCall(VecCopy(v1, Bx));
4519566063dSJacob Faibussowitsch     PetscCall(MatMult(A, Bx, v2)); /* v2 = A*B*x */
4529566063dSJacob Faibussowitsch     PetscCall(VecCopy(v2, v1));
453447fed29SStefano Zampini     if (rart) {
4549566063dSJacob Faibussowitsch       PetscCall(MatMult(B, v1, v3)); /* v3 = R*A*R^t*x */
455447fed29SStefano Zampini     } else {
4569566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(B, v1, v3)); /* v3 = Bt*A*B*x */
457447fed29SStefano Zampini     }
4589566063dSJacob Faibussowitsch     PetscCall(VecNorm(v4, NORM_2, &norm_abs));
4599566063dSJacob Faibussowitsch     PetscCall(VecAXPY(v4, none, v3));
4609566063dSJacob Faibussowitsch     PetscCall(VecNorm(v4, NORM_2, &norm_rel));
461447fed29SStefano Zampini 
462447fed29SStefano Zampini     if (norm_abs > tol) norm_rel /= norm_abs;
463447fed29SStefano Zampini     if (norm_rel > tol) {
464cc48ffa7SToby Isaac       *flg = PETSC_FALSE;
4659566063dSJacob Faibussowitsch       PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th Mat%sMult() %g\n", i, rart ? "RARt" : "PtAP", (double)norm_rel));
466cc48ffa7SToby Isaac       break;
467cc48ffa7SToby Isaac     }
468cc48ffa7SToby Isaac   }
469447fed29SStefano Zampini 
4709566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
4719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
4729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Bx));
4739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Cx));
4749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v1));
4759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v2));
4769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v3));
4779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v4));
4783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
479cc48ffa7SToby Isaac }
480cc48ffa7SToby Isaac 
48186919fd6SHong Zhang /*@
4824222ddf1SHong Zhang    MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B
4834222ddf1SHong Zhang 
484c3339decSBarry Smith    Collective
4854222ddf1SHong Zhang 
4864222ddf1SHong Zhang    Input Parameters:
4874222ddf1SHong Zhang +  A - the first matrix
4884222ddf1SHong Zhang .  B - the second matrix
4894222ddf1SHong Zhang .  C - the third matrix
4904222ddf1SHong Zhang -  n - number of random vectors to be tested
4914222ddf1SHong Zhang 
4924222ddf1SHong Zhang    Output Parameter:
4934222ddf1SHong Zhang .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
4944222ddf1SHong Zhang 
4954222ddf1SHong Zhang    Level: intermediate
4964222ddf1SHong Zhang 
4974222ddf1SHong Zhang @*/
498d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg)
499d71ae5a4SJacob Faibussowitsch {
5004222ddf1SHong Zhang   PetscFunctionBegin;
5019566063dSJacob Faibussowitsch   PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_FALSE, flg));
5023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5034222ddf1SHong Zhang }
5044222ddf1SHong Zhang 
505447fed29SStefano Zampini /*@
506447fed29SStefano Zampini    MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t
507447fed29SStefano Zampini 
508c3339decSBarry Smith    Collective
509447fed29SStefano Zampini 
510447fed29SStefano Zampini    Input Parameters:
511447fed29SStefano Zampini +  A - the first matrix
512447fed29SStefano Zampini .  B - the second matrix
513447fed29SStefano Zampini .  C - the third matrix
514447fed29SStefano Zampini -  n - number of random vectors to be tested
515447fed29SStefano Zampini 
516447fed29SStefano Zampini    Output Parameter:
517447fed29SStefano Zampini .  flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise.
518447fed29SStefano Zampini 
519447fed29SStefano Zampini    Level: intermediate
520447fed29SStefano Zampini 
521447fed29SStefano Zampini @*/
522d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRARtMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg)
523d71ae5a4SJacob Faibussowitsch {
524447fed29SStefano Zampini   PetscFunctionBegin;
5259566063dSJacob Faibussowitsch   PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_TRUE, flg));
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5274222ddf1SHong Zhang }
5284222ddf1SHong Zhang 
5294222ddf1SHong Zhang /*@
53086919fd6SHong Zhang    MatIsLinear - Check if a shell matrix A is a linear operator.
53186919fd6SHong Zhang 
532c3339decSBarry Smith    Collective
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 @*/
543d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsLinear(Mat A, PetscInt n, PetscBool *flg)
544d71ae5a4SJacob Faibussowitsch {
54586919fd6SHong Zhang   Vec         x, y, s1, s2;
54686919fd6SHong Zhang   PetscRandom rctx;
54786919fd6SHong Zhang   PetscScalar a;
54886919fd6SHong Zhang   PetscInt    k;
54986919fd6SHong Zhang   PetscReal   norm, normA;
55086919fd6SHong Zhang   MPI_Comm    comm;
55186919fd6SHong Zhang   PetscMPIInt rank;
55286919fd6SHong Zhang 
55386919fd6SHong Zhang   PetscFunctionBegin;
55486919fd6SHong Zhang   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
5559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
5569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
55786919fd6SHong Zhang 
5589566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(comm, &rctx));
5599566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
5609566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &s1));
5619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
5629566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(s1, &s2));
56386919fd6SHong Zhang 
56486919fd6SHong Zhang   *flg = PETSC_TRUE;
56586919fd6SHong Zhang   for (k = 0; k < n; k++) {
5669566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rctx));
5679566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y, rctx));
56848a46eb9SPierre Jolivet     if (rank == 0) PetscCall(PetscRandomGetValue(rctx, &a));
5699566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm));
57086919fd6SHong Zhang 
57186919fd6SHong Zhang     /* s2 = a*A*x + A*y */
5729566063dSJacob Faibussowitsch     PetscCall(MatMult(A, y, s2));  /* s2 = A*y */
5739566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, s1));  /* s1 = A*x */
5749566063dSJacob Faibussowitsch     PetscCall(VecAXPY(s2, a, s1)); /* s2 = a s1 + s2 */
57586919fd6SHong Zhang 
57686919fd6SHong Zhang     /* s1 = A * (a x + y) */
5779566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, a, x)); /* y = a x + y */
5789566063dSJacob Faibussowitsch     PetscCall(MatMult(A, y, s1));
5799566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_INFINITY, &normA));
58086919fd6SHong Zhang 
5819566063dSJacob Faibussowitsch     PetscCall(VecAXPY(s2, -1.0, s1)); /* s2 = - s1 + s2 */
5829566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_INFINITY, &norm));
5831b8dac88SHong Zhang     if (norm / normA > 100. * PETSC_MACHINE_EPSILON) {
58486919fd6SHong Zhang       *flg = PETSC_FALSE;
5859566063dSJacob 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)));
58686919fd6SHong Zhang       break;
58786919fd6SHong Zhang     }
58886919fd6SHong Zhang   }
5899566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
5909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
5919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
5929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
5939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
5943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59586919fd6SHong Zhang }
596