1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 286a22c91SHong Zhang 3274439ffSJunchao Zhang /* 4*d7c1f440SPierre Jolivet n; try the MatMult variant n times 5274439ffSJunchao Zhang flg: return the boolean result, equal or not 6274439ffSJunchao Zhang t: 0 => no transpose; 1 => transpose; 2 => Hermitian transpose 7274439ffSJunchao Zhang add: 0 => no add (e.g., y = Ax); 1 => add third vector (e.g., z = Ax + y); 2 => add update (e.g., y = Ax + y) 8274439ffSJunchao Zhang */ 936eb99e2SToby Isaac static PetscErrorCode MatMultEqual_Private(Mat A, Mat B, PetscInt n, PetscBool *flg, PetscInt t, PetscInt add) 10d71ae5a4SJacob Faibussowitsch { 11b84f494bSStefano Zampini Vec Ax = NULL, Bx = NULL, s1 = NULL, s2 = NULL, Ay = NULL, By = NULL; 12447fed29SStefano Zampini PetscRandom rctx; 13447fed29SStefano Zampini PetscReal r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON; 14447fed29SStefano Zampini PetscInt am, an, bm, bn, k; 15447fed29SStefano Zampini PetscScalar none = -1.0; 164279555eSSatish Balay #if defined(PETSC_USE_INFO) 1736eb99e2SToby Isaac const char *sops[] = {"MatMult", "MatMultAdd", "MatMultAdd (update)", "MatMultTranspose", "MatMultTransposeAdd", "MatMultTransposeAdd (update)", "MatMultHermitianTranspose", "MatMultHermitianTransposeAdd", "MatMultHermitianTransposeAdd (update)"}; 18447fed29SStefano Zampini const char *sop; 194279555eSSatish Balay #endif 20447fed29SStefano Zampini 21447fed29SStefano Zampini PetscFunctionBegin; 22447fed29SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 23447fed29SStefano Zampini PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 24447fed29SStefano Zampini PetscCheckSameComm(A, 1, B, 2); 25447fed29SStefano Zampini PetscValidLogicalCollectiveInt(A, n, 3); 264f572ea9SToby Isaac PetscAssertPointer(flg, 4); 27e573050aSPierre Jolivet PetscValidLogicalCollectiveInt(A, t, 5); 280d6f747bSJacob Faibussowitsch PetscValidLogicalCollectiveInt(A, add, 6); 299566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 309566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 31aed4548fSBarry 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); 324279555eSSatish Balay #if defined(PETSC_USE_INFO) 33274439ffSJunchao Zhang sop = sops[add + 3 * t]; 344279555eSSatish Balay #endif 359566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)A), &rctx)); 369566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 37447fed29SStefano Zampini if (t) { 389566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &s1, &Ax)); 399566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &s2, &Bx)); 40447fed29SStefano Zampini } else { 419566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Ax, &s1)); 429566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &Bx, &s2)); 43447fed29SStefano Zampini } 44447fed29SStefano Zampini if (add) { 459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s1, &Ay)); 469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s2, &By)); 47447fed29SStefano Zampini } 48447fed29SStefano Zampini 49447fed29SStefano Zampini *flg = PETSC_TRUE; 50447fed29SStefano Zampini for (k = 0; k < n; k++) { 5136eb99e2SToby Isaac Vec Aadd = NULL, Badd = NULL; 5236eb99e2SToby Isaac 539566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ax, rctx)); 549566063dSJacob Faibussowitsch PetscCall(VecCopy(Ax, Bx)); 55447fed29SStefano Zampini if (add) { 569566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ay, rctx)); 579566063dSJacob Faibussowitsch PetscCall(VecCopy(Ay, By)); 5836eb99e2SToby Isaac Aadd = Ay; 5936eb99e2SToby Isaac Badd = By; 6036eb99e2SToby Isaac if (add == 2) { 6136eb99e2SToby Isaac PetscCall(VecCopy(Ay, s1)); 6236eb99e2SToby Isaac PetscCall(VecCopy(By, s2)); 6336eb99e2SToby Isaac Aadd = s1; 6436eb99e2SToby Isaac Badd = s2; 6536eb99e2SToby Isaac } 66447fed29SStefano Zampini } 67e573050aSPierre Jolivet if (t == 1) { 68447fed29SStefano Zampini if (add) { 6936eb99e2SToby Isaac PetscCall(MatMultTransposeAdd(A, Ax, Aadd, s1)); 7036eb99e2SToby Isaac PetscCall(MatMultTransposeAdd(B, Bx, Badd, s2)); 71447fed29SStefano Zampini } else { 729566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Ax, s1)); 739566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Bx, s2)); 74447fed29SStefano Zampini } 75e573050aSPierre Jolivet } else if (t == 2) { 76e573050aSPierre Jolivet if (add) { 7736eb99e2SToby Isaac PetscCall(MatMultHermitianTransposeAdd(A, Ax, Aadd, s1)); 7836eb99e2SToby Isaac PetscCall(MatMultHermitianTransposeAdd(B, Bx, Badd, s2)); 79e573050aSPierre Jolivet } else { 809566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(A, Ax, s1)); 819566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(B, Bx, s2)); 82e573050aSPierre Jolivet } 83447fed29SStefano Zampini } else { 84447fed29SStefano Zampini if (add) { 8536eb99e2SToby Isaac PetscCall(MatMultAdd(A, Ax, Aadd, s1)); 8636eb99e2SToby Isaac PetscCall(MatMultAdd(B, Bx, Badd, s2)); 87447fed29SStefano Zampini } else { 889566063dSJacob Faibussowitsch PetscCall(MatMult(A, Ax, s1)); 899566063dSJacob Faibussowitsch PetscCall(MatMult(B, Bx, s2)); 90447fed29SStefano Zampini } 91447fed29SStefano Zampini } 929566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 93447fed29SStefano Zampini if (r2 < tol) { 949566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_INFINITY, &r1)); 95447fed29SStefano Zampini } else { 969566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, none, s1)); 979566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 98447fed29SStefano Zampini r1 /= r2; 99447fed29SStefano Zampini } 100447fed29SStefano Zampini if (r1 > tol) { 101447fed29SStefano Zampini *flg = PETSC_FALSE; 1029566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s() %g\n", k, sop, (double)r1)); 103447fed29SStefano Zampini break; 104447fed29SStefano Zampini } 105447fed29SStefano Zampini } 1069566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 1079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 1099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ay)); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&By)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114447fed29SStefano Zampini } 115447fed29SStefano Zampini 116d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg, PetscBool At, PetscBool Bt) 117d71ae5a4SJacob Faibussowitsch { 118447fed29SStefano Zampini Vec Ax, Bx, Cx, s1, s2, s3; 119447fed29SStefano Zampini PetscRandom rctx; 120447fed29SStefano Zampini PetscReal r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON; 121447fed29SStefano Zampini PetscInt am, an, bm, bn, cm, cn, k; 122447fed29SStefano Zampini PetscScalar none = -1.0; 1234279555eSSatish Balay #if defined(PETSC_USE_INFO) 124580c7c76SPierre Jolivet const char *sops[] = {"MatMatMult", "MatTransposeMatMult", "MatMatTransposeMult", "MatTransposeMatTransposeMult"}; 125447fed29SStefano Zampini const char *sop; 1264279555eSSatish Balay #endif 127447fed29SStefano Zampini 128447fed29SStefano Zampini PetscFunctionBegin; 129447fed29SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 130447fed29SStefano Zampini PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 131447fed29SStefano Zampini PetscCheckSameComm(A, 1, B, 2); 132447fed29SStefano Zampini PetscValidHeaderSpecific(C, MAT_CLASSID, 3); 133447fed29SStefano Zampini PetscCheckSameComm(A, 1, C, 3); 134447fed29SStefano Zampini PetscValidLogicalCollectiveInt(A, n, 4); 1354f572ea9SToby Isaac PetscAssertPointer(flg, 5); 136447fed29SStefano Zampini PetscValidLogicalCollectiveBool(A, At, 6); 137447fed29SStefano Zampini PetscValidLogicalCollectiveBool(B, Bt, 7); 1389566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 1399566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 1409566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &cm, &cn)); 1419371c9d4SSatish Balay if (At) { 1429371c9d4SSatish Balay PetscInt tt = an; 1439371c9d4SSatish Balay an = am; 1449371c9d4SSatish Balay am = tt; 145a8f51744SPierre Jolivet } 1469371c9d4SSatish Balay if (Bt) { 1479371c9d4SSatish Balay PetscInt tt = bn; 1489371c9d4SSatish Balay bn = bm; 1499371c9d4SSatish Balay bm = tt; 150a8f51744SPierre Jolivet } 151aed4548fSBarry 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); 152447fed29SStefano Zampini 1534279555eSSatish Balay #if defined(PETSC_USE_INFO) 154447fed29SStefano Zampini sop = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)]; 1554279555eSSatish Balay #endif 1569566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)C), &rctx)); 1579566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 158447fed29SStefano Zampini if (Bt) { 1599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &s1, &Bx)); 160447fed29SStefano Zampini } else { 1619566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &Bx, &s1)); 162447fed29SStefano Zampini } 163447fed29SStefano Zampini if (At) { 1649566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &s2, &Ax)); 165447fed29SStefano Zampini } else { 1669566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Ax, &s2)); 167447fed29SStefano Zampini } 1689566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &Cx, &s3)); 169447fed29SStefano Zampini 170447fed29SStefano Zampini *flg = PETSC_TRUE; 171447fed29SStefano Zampini for (k = 0; k < n; k++) { 1729566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Bx, rctx)); 173447fed29SStefano Zampini if (Bt) { 1749566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Bx, s1)); 175447fed29SStefano Zampini } else { 1769566063dSJacob Faibussowitsch PetscCall(MatMult(B, Bx, s1)); 177447fed29SStefano Zampini } 1789566063dSJacob Faibussowitsch PetscCall(VecCopy(s1, Ax)); 179447fed29SStefano Zampini if (At) { 1809566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Ax, s2)); 181447fed29SStefano Zampini } else { 1829566063dSJacob Faibussowitsch PetscCall(MatMult(A, Ax, s2)); 183447fed29SStefano Zampini } 1849566063dSJacob Faibussowitsch PetscCall(VecCopy(Bx, Cx)); 1859566063dSJacob Faibussowitsch PetscCall(MatMult(C, Cx, s3)); 186447fed29SStefano Zampini 1879566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 188447fed29SStefano Zampini if (r2 < tol) { 1899566063dSJacob Faibussowitsch PetscCall(VecNorm(s3, NORM_INFINITY, &r1)); 190447fed29SStefano Zampini } else { 1919566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, none, s3)); 1929566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 193447fed29SStefano Zampini r1 /= r2; 194447fed29SStefano Zampini } 195447fed29SStefano Zampini if (r1 > tol) { 196447fed29SStefano Zampini *flg = PETSC_FALSE; 1979566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s %g\n", k, sop, (double)r1)); 198447fed29SStefano Zampini break; 199447fed29SStefano Zampini } 200447fed29SStefano Zampini } 2019566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 2029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 2039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 2049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Cx)); 2059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s3)); 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 209447fed29SStefano Zampini } 210447fed29SStefano Zampini 21186a22c91SHong Zhang /*@ 21286a22c91SHong Zhang MatMultEqual - Compares matrix-vector products of two matrices. 21386a22c91SHong Zhang 214c3339decSBarry Smith Collective 21586a22c91SHong Zhang 21686a22c91SHong Zhang Input Parameters: 21786a22c91SHong Zhang + A - the first matrix 2184222ddf1SHong Zhang . B - the second matrix 21986a22c91SHong Zhang - n - number of random vectors to be tested 22086a22c91SHong Zhang 22186a22c91SHong Zhang Output Parameter: 22220f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 22386a22c91SHong Zhang 22486a22c91SHong Zhang Level: intermediate 22586a22c91SHong Zhang 22620f4b53cSBarry Smith .seealso: `Mat`, `MatMultAddEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()`, `MatIsLinear()` 22786a22c91SHong Zhang @*/ 228d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 229d71ae5a4SJacob Faibussowitsch { 23086a22c91SHong Zhang PetscFunctionBegin; 23136eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 0)); 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23386a22c91SHong Zhang } 23486a22c91SHong Zhang 23586a22c91SHong Zhang /*@ 23620f4b53cSBarry Smith MatMultAddEqual - Compares matrix-vector product plus vector add of two matrices. 23786a22c91SHong Zhang 238c3339decSBarry Smith Collective 23986a22c91SHong Zhang 24086a22c91SHong Zhang Input Parameters: 24186a22c91SHong Zhang + A - the first matrix 2424222ddf1SHong Zhang . B - the second matrix 24386a22c91SHong Zhang - n - number of random vectors to be tested 24486a22c91SHong Zhang 24586a22c91SHong Zhang Output Parameter: 24620f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 24786a22c91SHong Zhang 24886a22c91SHong Zhang Level: intermediate 24986a22c91SHong Zhang 25020f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()` 25186a22c91SHong Zhang @*/ 252d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 253d71ae5a4SJacob Faibussowitsch { 25486a22c91SHong Zhang PetscFunctionBegin; 25536eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 1)); 25636eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 2)); 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25863879571SHong Zhang } 25963879571SHong Zhang 26063879571SHong Zhang /*@ 26163879571SHong Zhang MatMultTransposeEqual - Compares matrix-vector products of two matrices. 26263879571SHong Zhang 263c3339decSBarry Smith Collective 26463879571SHong Zhang 26563879571SHong Zhang Input Parameters: 26663879571SHong Zhang + A - the first matrix 2674222ddf1SHong Zhang . B - the second matrix 26863879571SHong Zhang - n - number of random vectors to be tested 26963879571SHong Zhang 27063879571SHong Zhang Output Parameter: 27120f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 27263879571SHong Zhang 27363879571SHong Zhang Level: intermediate 27463879571SHong Zhang 27520f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeAddEqual()` 27663879571SHong Zhang @*/ 277d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 278d71ae5a4SJacob Faibussowitsch { 27963879571SHong Zhang PetscFunctionBegin; 28036eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 0)); 2813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28263879571SHong Zhang } 28363879571SHong Zhang 28463879571SHong Zhang /*@ 28563879571SHong Zhang MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 28663879571SHong Zhang 287c3339decSBarry Smith Collective 28863879571SHong Zhang 28963879571SHong Zhang Input Parameters: 29063879571SHong Zhang + A - the first matrix 2914222ddf1SHong Zhang . B - the second matrix 29263879571SHong Zhang - n - number of random vectors to be tested 29363879571SHong Zhang 29463879571SHong Zhang Output Parameter: 29520f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 29663879571SHong Zhang 29763879571SHong Zhang Level: intermediate 29863879571SHong Zhang 29920f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 30063879571SHong Zhang @*/ 301d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 302d71ae5a4SJacob Faibussowitsch { 30363879571SHong Zhang PetscFunctionBegin; 30436eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 1)); 30536eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 2)); 3063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 307e573050aSPierre Jolivet } 308e573050aSPierre Jolivet 309e573050aSPierre Jolivet /*@ 310e573050aSPierre Jolivet MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices. 311e573050aSPierre Jolivet 312c3339decSBarry Smith Collective 313e573050aSPierre Jolivet 314e573050aSPierre Jolivet Input Parameters: 315e573050aSPierre Jolivet + A - the first matrix 316e573050aSPierre Jolivet . B - the second matrix 317e573050aSPierre Jolivet - n - number of random vectors to be tested 318e573050aSPierre Jolivet 319e573050aSPierre Jolivet Output Parameter: 32020f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 321e573050aSPierre Jolivet 322e573050aSPierre Jolivet Level: intermediate 323e573050aSPierre Jolivet 32452cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 325e573050aSPierre Jolivet @*/ 326d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 327d71ae5a4SJacob Faibussowitsch { 328e573050aSPierre Jolivet PetscFunctionBegin; 32936eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 0)); 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 331e573050aSPierre Jolivet } 332e573050aSPierre Jolivet 333e573050aSPierre Jolivet /*@ 334e573050aSPierre Jolivet MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices. 335e573050aSPierre Jolivet 336c3339decSBarry Smith Collective 337e573050aSPierre Jolivet 338e573050aSPierre Jolivet Input Parameters: 339e573050aSPierre Jolivet + A - the first matrix 340e573050aSPierre Jolivet . B - the second matrix 341e573050aSPierre Jolivet - n - number of random vectors to be tested 342e573050aSPierre Jolivet 343e573050aSPierre Jolivet Output Parameter: 34420f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 345e573050aSPierre Jolivet 346e573050aSPierre Jolivet Level: intermediate 347e573050aSPierre Jolivet 34852cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 349e573050aSPierre Jolivet @*/ 350d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 351d71ae5a4SJacob Faibussowitsch { 352e573050aSPierre Jolivet PetscFunctionBegin; 35336eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 1)); 35436eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 2)); 3553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35686a22c91SHong Zhang } 357a52676ecSHong Zhang 358a52676ecSHong Zhang /*@ 359a52676ecSHong Zhang MatMatMultEqual - Test A*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 365c2916339SPierre Jolivet . B - the second matrix 366c2916339SPierre Jolivet . C - the third matrix 367a52676ecSHong Zhang - n - number of random vectors to be tested 368a52676ecSHong Zhang 369a52676ecSHong Zhang Output Parameter: 37020f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 371a52676ecSHong Zhang 372a52676ecSHong Zhang Level: intermediate 373a52676ecSHong Zhang 37452cd20c4SPierre Jolivet .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 375a52676ecSHong Zhang @*/ 376d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 377d71ae5a4SJacob Faibussowitsch { 378a52676ecSHong Zhang PetscFunctionBegin; 3799566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_FALSE)); 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381a52676ecSHong Zhang } 382a52676ecSHong Zhang 383a52676ecSHong Zhang /*@ 384a52676ecSHong Zhang MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x 385a52676ecSHong Zhang 386c3339decSBarry Smith Collective 387a52676ecSHong Zhang 388a52676ecSHong Zhang Input Parameters: 389a52676ecSHong Zhang + A - the first matrix 3904222ddf1SHong Zhang . B - the second matrix 3914222ddf1SHong Zhang . C - the third matrix 392a52676ecSHong Zhang - n - number of random vectors to be tested 393a52676ecSHong Zhang 394a52676ecSHong Zhang Output Parameter: 39520f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 396a52676ecSHong Zhang 397a52676ecSHong Zhang Level: intermediate 398a52676ecSHong Zhang 39952cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 400a52676ecSHong Zhang @*/ 401d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTransposeMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 402d71ae5a4SJacob Faibussowitsch { 403a52676ecSHong Zhang PetscFunctionBegin; 4049566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_TRUE, PETSC_FALSE)); 4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 406a52676ecSHong Zhang } 40786919fd6SHong Zhang 40826546c1bSToby Isaac /*@ 40926546c1bSToby Isaac MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x 41026546c1bSToby Isaac 411c3339decSBarry Smith Collective 41226546c1bSToby Isaac 41326546c1bSToby Isaac Input Parameters: 41426546c1bSToby Isaac + A - the first matrix 4154222ddf1SHong Zhang . B - the second matrix 4164222ddf1SHong Zhang . C - the third matrix 41726546c1bSToby Isaac - n - number of random vectors to be tested 41826546c1bSToby Isaac 41926546c1bSToby Isaac Output Parameter: 42020f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 42126546c1bSToby Isaac 42226546c1bSToby Isaac Level: intermediate 42326546c1bSToby Isaac 42452cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 42526546c1bSToby Isaac @*/ 426d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatTransposeMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 427d71ae5a4SJacob Faibussowitsch { 428447fed29SStefano Zampini PetscFunctionBegin; 4299566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_TRUE)); 4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 431447fed29SStefano Zampini } 432447fed29SStefano Zampini 433d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProjMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool rart, PetscBool *flg) 434d71ae5a4SJacob Faibussowitsch { 435447fed29SStefano Zampini Vec x, v1, v2, v3, v4, Cx, Bx; 436447fed29SStefano Zampini PetscReal norm_abs, norm_rel, tol = PETSC_SQRT_MACHINE_EPSILON; 437447fed29SStefano Zampini PetscInt i, am, an, bm, bn, cm, cn; 438447fed29SStefano Zampini PetscRandom rdm; 439cc48ffa7SToby Isaac PetscScalar none = -1.0; 440cc48ffa7SToby Isaac 441cc48ffa7SToby Isaac PetscFunctionBegin; 4429566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 4439566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 4449371c9d4SSatish Balay if (rart) { 4459371c9d4SSatish Balay PetscInt t = bm; 4469371c9d4SSatish Balay bm = bn; 4479371c9d4SSatish Balay bn = t; 4489371c9d4SSatish Balay } 4499566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &cm, &cn)); 450aed4548fSBarry 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); 451cc48ffa7SToby Isaac 452447fed29SStefano Zampini /* Create left vector of A: v2 */ 4539566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Bx, &v2)); 454447fed29SStefano Zampini 455447fed29SStefano Zampini /* Create right vectors of B: x, v3, v4 */ 456447fed29SStefano Zampini if (rart) { 4579566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &v1, &x)); 458447fed29SStefano Zampini } else { 4599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &x, &v1)); 460447fed29SStefano Zampini } 4619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &v3)); 462447fed29SStefano Zampini 4639566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &Cx, &v4)); 4649566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 4659566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 466cc48ffa7SToby Isaac 467cc48ffa7SToby Isaac *flg = PETSC_TRUE; 468447fed29SStefano Zampini for (i = 0; i < n; i++) { 4699566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 4709566063dSJacob Faibussowitsch PetscCall(VecCopy(x, Cx)); 4719566063dSJacob Faibussowitsch PetscCall(MatMult(C, Cx, v4)); /* v4 = C*x */ 472447fed29SStefano Zampini if (rart) { 4739566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, x, v1)); 474cc48ffa7SToby Isaac } else { 4759566063dSJacob Faibussowitsch PetscCall(MatMult(B, x, v1)); 476cc48ffa7SToby Isaac } 4779566063dSJacob Faibussowitsch PetscCall(VecCopy(v1, Bx)); 4789566063dSJacob Faibussowitsch PetscCall(MatMult(A, Bx, v2)); /* v2 = A*B*x */ 4799566063dSJacob Faibussowitsch PetscCall(VecCopy(v2, v1)); 480447fed29SStefano Zampini if (rart) { 4819566063dSJacob Faibussowitsch PetscCall(MatMult(B, v1, v3)); /* v3 = R*A*R^t*x */ 482447fed29SStefano Zampini } else { 4839566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, v1, v3)); /* v3 = Bt*A*B*x */ 484447fed29SStefano Zampini } 4859566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_abs)); 4869566063dSJacob Faibussowitsch PetscCall(VecAXPY(v4, none, v3)); 4879566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_rel)); 488447fed29SStefano Zampini 489447fed29SStefano Zampini if (norm_abs > tol) norm_rel /= norm_abs; 490447fed29SStefano Zampini if (norm_rel > tol) { 491cc48ffa7SToby Isaac *flg = PETSC_FALSE; 4929566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th Mat%sMult() %g\n", i, rart ? "RARt" : "PtAP", (double)norm_rel)); 493cc48ffa7SToby Isaac break; 494cc48ffa7SToby Isaac } 495cc48ffa7SToby Isaac } 496447fed29SStefano Zampini 4979566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 4989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 4999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 5009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Cx)); 5019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v1)); 5029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v2)); 5039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v3)); 5049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v4)); 5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 506cc48ffa7SToby Isaac } 507cc48ffa7SToby Isaac 50886919fd6SHong Zhang /*@ 5094222ddf1SHong Zhang MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B 5104222ddf1SHong Zhang 511c3339decSBarry Smith Collective 5124222ddf1SHong Zhang 5134222ddf1SHong Zhang Input Parameters: 5144222ddf1SHong Zhang + A - the first matrix 5154222ddf1SHong Zhang . B - the second matrix 5164222ddf1SHong Zhang . C - the third matrix 5174222ddf1SHong Zhang - n - number of random vectors to be tested 5184222ddf1SHong Zhang 5194222ddf1SHong Zhang Output Parameter: 52020f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 5214222ddf1SHong Zhang 5224222ddf1SHong Zhang Level: intermediate 5234222ddf1SHong Zhang 52452cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 5254222ddf1SHong Zhang @*/ 526d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 527d71ae5a4SJacob Faibussowitsch { 5284222ddf1SHong Zhang PetscFunctionBegin; 5299566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_FALSE, flg)); 5303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5314222ddf1SHong Zhang } 5324222ddf1SHong Zhang 533447fed29SStefano Zampini /*@ 534447fed29SStefano Zampini MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t 535447fed29SStefano Zampini 536c3339decSBarry Smith Collective 537447fed29SStefano Zampini 538447fed29SStefano Zampini Input Parameters: 539447fed29SStefano Zampini + A - the first matrix 540447fed29SStefano Zampini . B - the second matrix 541447fed29SStefano Zampini . C - the third matrix 542447fed29SStefano Zampini - n - number of random vectors to be tested 543447fed29SStefano Zampini 544447fed29SStefano Zampini Output Parameter: 54520f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 546447fed29SStefano Zampini 547447fed29SStefano Zampini Level: intermediate 548447fed29SStefano Zampini 54952cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 550447fed29SStefano Zampini @*/ 551d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRARtMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 552d71ae5a4SJacob Faibussowitsch { 553447fed29SStefano Zampini PetscFunctionBegin; 5549566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_TRUE, flg)); 5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5564222ddf1SHong Zhang } 5574222ddf1SHong Zhang 5584222ddf1SHong Zhang /*@ 55920f4b53cSBarry Smith MatIsLinear - Check if a shell matrix `A` is a linear operator. 56086919fd6SHong Zhang 561c3339decSBarry Smith Collective 56286919fd6SHong Zhang 56386919fd6SHong Zhang Input Parameters: 56486919fd6SHong Zhang + A - the shell matrix 56586919fd6SHong Zhang - n - number of random vectors to be tested 56686919fd6SHong Zhang 56786919fd6SHong Zhang Output Parameter: 56820f4b53cSBarry Smith . flg - `PETSC_TRUE` if the shell matrix is linear; `PETSC_FALSE` otherwise. 56986919fd6SHong Zhang 57086919fd6SHong Zhang Level: intermediate 57120f4b53cSBarry Smith 57252cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 57386919fd6SHong Zhang @*/ 574d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsLinear(Mat A, PetscInt n, PetscBool *flg) 575d71ae5a4SJacob Faibussowitsch { 57686919fd6SHong Zhang Vec x, y, s1, s2; 57786919fd6SHong Zhang PetscRandom rctx; 57886919fd6SHong Zhang PetscScalar a; 57986919fd6SHong Zhang PetscInt k; 58086919fd6SHong Zhang PetscReal norm, normA; 58186919fd6SHong Zhang MPI_Comm comm; 58286919fd6SHong Zhang PetscMPIInt rank; 58386919fd6SHong Zhang 58486919fd6SHong Zhang PetscFunctionBegin; 58586919fd6SHong Zhang PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5869566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 5879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 58886919fd6SHong Zhang 5899566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rctx)); 5909566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 5919566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &s1)); 5929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 5939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s1, &s2)); 59486919fd6SHong Zhang 59586919fd6SHong Zhang *flg = PETSC_TRUE; 59686919fd6SHong Zhang for (k = 0; k < n; k++) { 5979566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 5989566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 59948a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscRandomGetValue(rctx, &a)); 6009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm)); 60186919fd6SHong Zhang 60286919fd6SHong Zhang /* s2 = a*A*x + A*y */ 6039566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s2)); /* s2 = A*y */ 6049566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, s1)); /* s1 = A*x */ 6059566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, a, s1)); /* s2 = a s1 + s2 */ 60686919fd6SHong Zhang 60786919fd6SHong Zhang /* s1 = A * (a x + y) */ 6089566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, a, x)); /* y = a x + y */ 6099566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s1)); 6109566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_INFINITY, &normA)); 61186919fd6SHong Zhang 6129566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, -1.0, s1)); /* s2 = - s1 + s2 */ 6139566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &norm)); 6141b8dac88SHong Zhang if (norm / normA > 100. * PETSC_MACHINE_EPSILON) { 61586919fd6SHong Zhang *flg = PETSC_FALSE; 6169566063dSJacob 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))); 61786919fd6SHong Zhang break; 61886919fd6SHong Zhang } 61986919fd6SHong Zhang } 6209566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 6219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 6229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 6239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 6249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 6253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62686919fd6SHong Zhang } 627