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; 114279555eSSatish Balay #if defined(PETSC_USE_INFO) 12580c7c76SPierre Jolivet const char *sops[] = {"MatMult", "MatMultAdd", "MatMultTranspose", "MatMultTransposeAdd", "MatMultHermitianTranspose", "MatMultHermitianTransposeAdd"}; 13447fed29SStefano Zampini const char *sop; 144279555eSSatish 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); 274279555eSSatish 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 */ 294279555eSSatish 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; 1084279555eSSatish Balay #if defined(PETSC_USE_INFO) 109580c7c76SPierre Jolivet const char *sops[] = {"MatMatMult", "MatTransposeMatMult", "MatMatTransposeMult", "MatTransposeMatTransposeMult"}; 110447fed29SStefano Zampini const char *sop; 1114279555eSSatish 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 1384279555eSSatish Balay #if defined(PETSC_USE_INFO) 139447fed29SStefano Zampini sop = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)]; 1404279555eSSatish 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: 207*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 20886a22c91SHong Zhang 20986a22c91SHong Zhang Level: intermediate 21086a22c91SHong Zhang 211*20f4b53cSBarry Smith .seealso: `Mat`, `MatMultAddEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()`, `MatIsLinear()` 21286a22c91SHong Zhang @*/ 213d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 214d71ae5a4SJacob Faibussowitsch { 21586a22c91SHong Zhang PetscFunctionBegin; 2169566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_FALSE)); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21886a22c91SHong Zhang } 21986a22c91SHong Zhang 22086a22c91SHong Zhang /*@ 221*20f4b53cSBarry Smith MatMultAddEqual - Compares matrix-vector product plus vector add of two matrices. 22286a22c91SHong Zhang 223c3339decSBarry Smith Collective 22486a22c91SHong Zhang 22586a22c91SHong Zhang Input Parameters: 22686a22c91SHong Zhang + A - the first matrix 2274222ddf1SHong Zhang . B - the second matrix 22886a22c91SHong Zhang - n - number of random vectors to be tested 22986a22c91SHong Zhang 23086a22c91SHong Zhang Output Parameter: 231*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 23286a22c91SHong Zhang 23386a22c91SHong Zhang Level: intermediate 23486a22c91SHong Zhang 235*20f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()` 23686a22c91SHong Zhang @*/ 237d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 238d71ae5a4SJacob Faibussowitsch { 23986a22c91SHong Zhang PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_TRUE)); 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24263879571SHong Zhang } 24363879571SHong Zhang 24463879571SHong Zhang /*@ 24563879571SHong Zhang MatMultTransposeEqual - Compares matrix-vector products of two matrices. 24663879571SHong Zhang 247c3339decSBarry Smith Collective 24863879571SHong Zhang 24963879571SHong Zhang Input Parameters: 25063879571SHong Zhang + A - the first matrix 2514222ddf1SHong Zhang . B - the second matrix 25263879571SHong Zhang - n - number of random vectors to be tested 25363879571SHong Zhang 25463879571SHong Zhang Output Parameter: 255*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 25663879571SHong Zhang 25763879571SHong Zhang Level: intermediate 25863879571SHong Zhang 259*20f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeAddEqual()` 26063879571SHong Zhang @*/ 261d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 262d71ae5a4SJacob Faibussowitsch { 26363879571SHong Zhang PetscFunctionBegin; 2649566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_FALSE)); 2653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26663879571SHong Zhang } 26763879571SHong Zhang 26863879571SHong Zhang /*@ 26963879571SHong Zhang MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 27063879571SHong Zhang 271c3339decSBarry Smith Collective 27263879571SHong Zhang 27363879571SHong Zhang Input Parameters: 27463879571SHong Zhang + A - the first matrix 2754222ddf1SHong Zhang . B - the second matrix 27663879571SHong Zhang - n - number of random vectors to be tested 27763879571SHong Zhang 27863879571SHong Zhang Output Parameter: 279*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 28063879571SHong Zhang 28163879571SHong Zhang Level: intermediate 28263879571SHong Zhang 283*20f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 28463879571SHong Zhang @*/ 285d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 286d71ae5a4SJacob Faibussowitsch { 28763879571SHong Zhang PetscFunctionBegin; 2889566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_TRUE)); 2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 290e573050aSPierre Jolivet } 291e573050aSPierre Jolivet 292e573050aSPierre Jolivet /*@ 293e573050aSPierre Jolivet MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices. 294e573050aSPierre Jolivet 295c3339decSBarry Smith Collective 296e573050aSPierre Jolivet 297e573050aSPierre Jolivet Input Parameters: 298e573050aSPierre Jolivet + A - the first matrix 299e573050aSPierre Jolivet . B - the second matrix 300e573050aSPierre Jolivet - n - number of random vectors to be tested 301e573050aSPierre Jolivet 302e573050aSPierre Jolivet Output Parameter: 303*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 304e573050aSPierre Jolivet 305e573050aSPierre Jolivet Level: intermediate 306e573050aSPierre Jolivet 307*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 308e573050aSPierre Jolivet @*/ 309d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 310d71ae5a4SJacob Faibussowitsch { 311e573050aSPierre Jolivet PetscFunctionBegin; 3129566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 2, PETSC_FALSE)); 3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 314e573050aSPierre Jolivet } 315e573050aSPierre Jolivet 316e573050aSPierre Jolivet /*@ 317e573050aSPierre Jolivet MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices. 318e573050aSPierre Jolivet 319c3339decSBarry Smith Collective 320e573050aSPierre Jolivet 321e573050aSPierre Jolivet Input Parameters: 322e573050aSPierre Jolivet + A - the first matrix 323e573050aSPierre Jolivet . B - the second matrix 324e573050aSPierre Jolivet - n - number of random vectors to be tested 325e573050aSPierre Jolivet 326e573050aSPierre Jolivet Output Parameter: 327*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 328e573050aSPierre Jolivet 329e573050aSPierre Jolivet Level: intermediate 330e573050aSPierre Jolivet 331*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 332e573050aSPierre Jolivet @*/ 333d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 334d71ae5a4SJacob Faibussowitsch { 335e573050aSPierre Jolivet PetscFunctionBegin; 3369566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 2, PETSC_TRUE)); 3373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33886a22c91SHong Zhang } 339a52676ecSHong Zhang 340a52676ecSHong Zhang /*@ 341a52676ecSHong Zhang MatMatMultEqual - Test A*B*x = C*x for n random vector x 342a52676ecSHong Zhang 343c3339decSBarry Smith Collective 344a52676ecSHong Zhang 345a52676ecSHong Zhang Input Parameters: 346a52676ecSHong Zhang + A - the first matrix 347c2916339SPierre Jolivet . B - the second matrix 348c2916339SPierre Jolivet . C - the third matrix 349a52676ecSHong Zhang - n - number of random vectors to be tested 350a52676ecSHong Zhang 351a52676ecSHong Zhang Output Parameter: 352*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 353a52676ecSHong Zhang 354a52676ecSHong Zhang Level: intermediate 355a52676ecSHong Zhang 356*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 357a52676ecSHong Zhang @*/ 358d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 359d71ae5a4SJacob Faibussowitsch { 360a52676ecSHong Zhang PetscFunctionBegin; 3619566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_FALSE)); 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 363a52676ecSHong Zhang } 364a52676ecSHong Zhang 365a52676ecSHong Zhang /*@ 366a52676ecSHong Zhang MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x 367a52676ecSHong Zhang 368c3339decSBarry Smith Collective 369a52676ecSHong Zhang 370a52676ecSHong Zhang Input Parameters: 371a52676ecSHong Zhang + A - the first matrix 3724222ddf1SHong Zhang . B - the second matrix 3734222ddf1SHong Zhang . C - the third matrix 374a52676ecSHong Zhang - n - number of random vectors to be tested 375a52676ecSHong Zhang 376a52676ecSHong Zhang Output Parameter: 377*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 378a52676ecSHong Zhang 379a52676ecSHong Zhang Level: intermediate 380a52676ecSHong Zhang 381*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 382a52676ecSHong Zhang @*/ 383d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTransposeMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 384d71ae5a4SJacob Faibussowitsch { 385a52676ecSHong Zhang PetscFunctionBegin; 3869566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_TRUE, PETSC_FALSE)); 3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 388a52676ecSHong Zhang } 38986919fd6SHong Zhang 39026546c1bSToby Isaac /*@ 39126546c1bSToby Isaac MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x 39226546c1bSToby Isaac 393c3339decSBarry Smith Collective 39426546c1bSToby Isaac 39526546c1bSToby Isaac Input Parameters: 39626546c1bSToby Isaac + A - the first matrix 3974222ddf1SHong Zhang . B - the second matrix 3984222ddf1SHong Zhang . C - the third matrix 39926546c1bSToby Isaac - n - number of random vectors to be tested 40026546c1bSToby Isaac 40126546c1bSToby Isaac Output Parameter: 402*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 40326546c1bSToby Isaac 40426546c1bSToby Isaac Level: intermediate 40526546c1bSToby Isaac 406*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 40726546c1bSToby Isaac @*/ 408d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatTransposeMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 409d71ae5a4SJacob Faibussowitsch { 410447fed29SStefano Zampini PetscFunctionBegin; 4119566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_TRUE)); 4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 413447fed29SStefano Zampini } 414447fed29SStefano Zampini 415d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProjMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool rart, PetscBool *flg) 416d71ae5a4SJacob Faibussowitsch { 417447fed29SStefano Zampini Vec x, v1, v2, v3, v4, Cx, Bx; 418447fed29SStefano Zampini PetscReal norm_abs, norm_rel, tol = PETSC_SQRT_MACHINE_EPSILON; 419447fed29SStefano Zampini PetscInt i, am, an, bm, bn, cm, cn; 420447fed29SStefano Zampini PetscRandom rdm; 421cc48ffa7SToby Isaac PetscScalar none = -1.0; 422cc48ffa7SToby Isaac 423cc48ffa7SToby Isaac PetscFunctionBegin; 4249566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 4259566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 4269371c9d4SSatish Balay if (rart) { 4279371c9d4SSatish Balay PetscInt t = bm; 4289371c9d4SSatish Balay bm = bn; 4299371c9d4SSatish Balay bn = t; 4309371c9d4SSatish Balay } 4319566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &cm, &cn)); 432aed4548fSBarry 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); 433cc48ffa7SToby Isaac 434447fed29SStefano Zampini /* Create left vector of A: v2 */ 4359566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Bx, &v2)); 436447fed29SStefano Zampini 437447fed29SStefano Zampini /* Create right vectors of B: x, v3, v4 */ 438447fed29SStefano Zampini if (rart) { 4399566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &v1, &x)); 440447fed29SStefano Zampini } else { 4419566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &x, &v1)); 442447fed29SStefano Zampini } 4439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &v3)); 444447fed29SStefano Zampini 4459566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &Cx, &v4)); 4469566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 4479566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 448cc48ffa7SToby Isaac 449cc48ffa7SToby Isaac *flg = PETSC_TRUE; 450447fed29SStefano Zampini for (i = 0; i < n; i++) { 4519566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 4529566063dSJacob Faibussowitsch PetscCall(VecCopy(x, Cx)); 4539566063dSJacob Faibussowitsch PetscCall(MatMult(C, Cx, v4)); /* v4 = C*x */ 454447fed29SStefano Zampini if (rart) { 4559566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, x, v1)); 456cc48ffa7SToby Isaac } else { 4579566063dSJacob Faibussowitsch PetscCall(MatMult(B, x, v1)); 458cc48ffa7SToby Isaac } 4599566063dSJacob Faibussowitsch PetscCall(VecCopy(v1, Bx)); 4609566063dSJacob Faibussowitsch PetscCall(MatMult(A, Bx, v2)); /* v2 = A*B*x */ 4619566063dSJacob Faibussowitsch PetscCall(VecCopy(v2, v1)); 462447fed29SStefano Zampini if (rart) { 4639566063dSJacob Faibussowitsch PetscCall(MatMult(B, v1, v3)); /* v3 = R*A*R^t*x */ 464447fed29SStefano Zampini } else { 4659566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, v1, v3)); /* v3 = Bt*A*B*x */ 466447fed29SStefano Zampini } 4679566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_abs)); 4689566063dSJacob Faibussowitsch PetscCall(VecAXPY(v4, none, v3)); 4699566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_rel)); 470447fed29SStefano Zampini 471447fed29SStefano Zampini if (norm_abs > tol) norm_rel /= norm_abs; 472447fed29SStefano Zampini if (norm_rel > tol) { 473cc48ffa7SToby Isaac *flg = PETSC_FALSE; 4749566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th Mat%sMult() %g\n", i, rart ? "RARt" : "PtAP", (double)norm_rel)); 475cc48ffa7SToby Isaac break; 476cc48ffa7SToby Isaac } 477cc48ffa7SToby Isaac } 478447fed29SStefano Zampini 4799566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 4819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 4829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Cx)); 4839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v1)); 4849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v2)); 4859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v3)); 4869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v4)); 4873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 488cc48ffa7SToby Isaac } 489cc48ffa7SToby Isaac 49086919fd6SHong Zhang /*@ 4914222ddf1SHong Zhang MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B 4924222ddf1SHong Zhang 493c3339decSBarry Smith Collective 4944222ddf1SHong Zhang 4954222ddf1SHong Zhang Input Parameters: 4964222ddf1SHong Zhang + A - the first matrix 4974222ddf1SHong Zhang . B - the second matrix 4984222ddf1SHong Zhang . C - the third matrix 4994222ddf1SHong Zhang - n - number of random vectors to be tested 5004222ddf1SHong Zhang 5014222ddf1SHong Zhang Output Parameter: 502*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 5034222ddf1SHong Zhang 5044222ddf1SHong Zhang Level: intermediate 5054222ddf1SHong Zhang 506*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 5074222ddf1SHong Zhang @*/ 508d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 509d71ae5a4SJacob Faibussowitsch { 5104222ddf1SHong Zhang PetscFunctionBegin; 5119566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_FALSE, flg)); 5123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5134222ddf1SHong Zhang } 5144222ddf1SHong Zhang 515447fed29SStefano Zampini /*@ 516447fed29SStefano Zampini MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t 517447fed29SStefano Zampini 518c3339decSBarry Smith Collective 519447fed29SStefano Zampini 520447fed29SStefano Zampini Input Parameters: 521447fed29SStefano Zampini + A - the first matrix 522447fed29SStefano Zampini . B - the second matrix 523447fed29SStefano Zampini . C - the third matrix 524447fed29SStefano Zampini - n - number of random vectors to be tested 525447fed29SStefano Zampini 526447fed29SStefano Zampini Output Parameter: 527*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 528447fed29SStefano Zampini 529447fed29SStefano Zampini Level: intermediate 530447fed29SStefano Zampini 531*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 532447fed29SStefano Zampini @*/ 533d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRARtMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 534d71ae5a4SJacob Faibussowitsch { 535447fed29SStefano Zampini PetscFunctionBegin; 5369566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_TRUE, flg)); 5373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5384222ddf1SHong Zhang } 5394222ddf1SHong Zhang 5404222ddf1SHong Zhang /*@ 541*20f4b53cSBarry Smith MatIsLinear - Check if a shell matrix `A` is a linear operator. 54286919fd6SHong Zhang 543c3339decSBarry Smith Collective 54486919fd6SHong Zhang 54586919fd6SHong Zhang Input Parameters: 54686919fd6SHong Zhang + A - the shell matrix 54786919fd6SHong Zhang - n - number of random vectors to be tested 54886919fd6SHong Zhang 54986919fd6SHong Zhang Output Parameter: 550*20f4b53cSBarry Smith . flg - `PETSC_TRUE` if the shell matrix is linear; `PETSC_FALSE` otherwise. 55186919fd6SHong Zhang 55286919fd6SHong Zhang Level: intermediate 553*20f4b53cSBarry Smith 554*20f4b53cSBarry Smith .seealso: `Mat`, `MatMatMultAddEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 55586919fd6SHong Zhang @*/ 556d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsLinear(Mat A, PetscInt n, PetscBool *flg) 557d71ae5a4SJacob Faibussowitsch { 55886919fd6SHong Zhang Vec x, y, s1, s2; 55986919fd6SHong Zhang PetscRandom rctx; 56086919fd6SHong Zhang PetscScalar a; 56186919fd6SHong Zhang PetscInt k; 56286919fd6SHong Zhang PetscReal norm, normA; 56386919fd6SHong Zhang MPI_Comm comm; 56486919fd6SHong Zhang PetscMPIInt rank; 56586919fd6SHong Zhang 56686919fd6SHong Zhang PetscFunctionBegin; 56786919fd6SHong Zhang PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 5699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 57086919fd6SHong Zhang 5719566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rctx)); 5729566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 5739566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &s1)); 5749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 5759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s1, &s2)); 57686919fd6SHong Zhang 57786919fd6SHong Zhang *flg = PETSC_TRUE; 57886919fd6SHong Zhang for (k = 0; k < n; k++) { 5799566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 5809566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 58148a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscRandomGetValue(rctx, &a)); 5829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm)); 58386919fd6SHong Zhang 58486919fd6SHong Zhang /* s2 = a*A*x + A*y */ 5859566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s2)); /* s2 = A*y */ 5869566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, s1)); /* s1 = A*x */ 5879566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, a, s1)); /* s2 = a s1 + s2 */ 58886919fd6SHong Zhang 58986919fd6SHong Zhang /* s1 = A * (a x + y) */ 5909566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, a, x)); /* y = a x + y */ 5919566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s1)); 5929566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_INFINITY, &normA)); 59386919fd6SHong Zhang 5949566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, -1.0, s1)); /* s2 = - s1 + s2 */ 5959566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &norm)); 5961b8dac88SHong Zhang if (norm / normA > 100. * PETSC_MACHINE_EPSILON) { 59786919fd6SHong Zhang *flg = PETSC_FALSE; 5989566063dSJacob 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))); 59986919fd6SHong Zhang break; 60086919fd6SHong Zhang } 60186919fd6SHong Zhang } 6029566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 6039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 6049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 6059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 6069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 6073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60886919fd6SHong Zhang } 609