186a22c91SHong Zhang 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 386a22c91SHong Zhang 436eb99e2SToby Isaac static PetscErrorCode MatMultEqual_Private(Mat A, Mat B, PetscInt n, PetscBool *flg, PetscInt t, PetscInt 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) 1236eb99e2SToby Isaac const char *sops[] = {"MatMult", "MatMultAdd", "MatMultAdd (update)", "MatMultTranspose", "MatMultTransposeAdd", "MatMultTransposeAdd (update)", "MatMultHermitianTranspose", "MatMultHermitianTransposeAdd", "MatMultHermitianTransposeAdd (update)"}; 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); 230d6f747bSJacob Faibussowitsch PetscValidLogicalCollectiveInt(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) 2836eb99e2SToby Isaac sop = sops[add + 3 * t]; /* add = 0 => no add, add = 1 => add third vector, add = 2 => add update, 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++) { 4636eb99e2SToby Isaac Vec Aadd = NULL, Badd = NULL; 4736eb99e2SToby Isaac 489566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ax, rctx)); 499566063dSJacob Faibussowitsch PetscCall(VecCopy(Ax, Bx)); 50447fed29SStefano Zampini if (add) { 519566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ay, rctx)); 529566063dSJacob Faibussowitsch PetscCall(VecCopy(Ay, By)); 5336eb99e2SToby Isaac Aadd = Ay; 5436eb99e2SToby Isaac Badd = By; 5536eb99e2SToby Isaac if (add == 2) { 5636eb99e2SToby Isaac PetscCall(VecCopy(Ay, s1)); 5736eb99e2SToby Isaac PetscCall(VecCopy(By, s2)); 5836eb99e2SToby Isaac Aadd = s1; 5936eb99e2SToby Isaac Badd = s2; 6036eb99e2SToby Isaac } 61447fed29SStefano Zampini } 62e573050aSPierre Jolivet if (t == 1) { 63447fed29SStefano Zampini if (add) { 6436eb99e2SToby Isaac PetscCall(MatMultTransposeAdd(A, Ax, Aadd, s1)); 6536eb99e2SToby Isaac PetscCall(MatMultTransposeAdd(B, Bx, Badd, s2)); 66447fed29SStefano Zampini } else { 679566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Ax, s1)); 689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Bx, s2)); 69447fed29SStefano Zampini } 70e573050aSPierre Jolivet } else if (t == 2) { 71e573050aSPierre Jolivet if (add) { 7236eb99e2SToby Isaac PetscCall(MatMultHermitianTransposeAdd(A, Ax, Aadd, s1)); 7336eb99e2SToby Isaac PetscCall(MatMultHermitianTransposeAdd(B, Bx, Badd, s2)); 74e573050aSPierre Jolivet } else { 759566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(A, Ax, s1)); 769566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(B, Bx, s2)); 77e573050aSPierre Jolivet } 78447fed29SStefano Zampini } else { 79447fed29SStefano Zampini if (add) { 8036eb99e2SToby Isaac PetscCall(MatMultAdd(A, Ax, Aadd, s1)); 8136eb99e2SToby Isaac PetscCall(MatMultAdd(B, Bx, Badd, s2)); 82447fed29SStefano Zampini } else { 839566063dSJacob Faibussowitsch PetscCall(MatMult(A, Ax, s1)); 849566063dSJacob Faibussowitsch PetscCall(MatMult(B, Bx, s2)); 85447fed29SStefano Zampini } 86447fed29SStefano Zampini } 879566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 88447fed29SStefano Zampini if (r2 < tol) { 899566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_INFINITY, &r1)); 90447fed29SStefano Zampini } else { 919566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, none, s1)); 929566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 93447fed29SStefano Zampini r1 /= r2; 94447fed29SStefano Zampini } 95447fed29SStefano Zampini if (r1 > tol) { 96447fed29SStefano Zampini *flg = PETSC_FALSE; 979566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s() %g\n", k, sop, (double)r1)); 98447fed29SStefano Zampini break; 99447fed29SStefano Zampini } 100447fed29SStefano Zampini } 1019566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 1029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 1039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 1049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ay)); 1059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&By)); 1069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109447fed29SStefano Zampini } 110447fed29SStefano Zampini 111d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg, PetscBool At, PetscBool Bt) 112d71ae5a4SJacob Faibussowitsch { 113447fed29SStefano Zampini Vec Ax, Bx, Cx, s1, s2, s3; 114447fed29SStefano Zampini PetscRandom rctx; 115447fed29SStefano Zampini PetscReal r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON; 116447fed29SStefano Zampini PetscInt am, an, bm, bn, cm, cn, k; 117447fed29SStefano Zampini PetscScalar none = -1.0; 1184279555eSSatish Balay #if defined(PETSC_USE_INFO) 119580c7c76SPierre Jolivet const char *sops[] = {"MatMatMult", "MatTransposeMatMult", "MatMatTransposeMult", "MatTransposeMatTransposeMult"}; 120447fed29SStefano Zampini const char *sop; 1214279555eSSatish Balay #endif 122447fed29SStefano Zampini 123447fed29SStefano Zampini PetscFunctionBegin; 124447fed29SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 125447fed29SStefano Zampini PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 126447fed29SStefano Zampini PetscCheckSameComm(A, 1, B, 2); 127447fed29SStefano Zampini PetscValidHeaderSpecific(C, MAT_CLASSID, 3); 128447fed29SStefano Zampini PetscCheckSameComm(A, 1, C, 3); 129447fed29SStefano Zampini PetscValidLogicalCollectiveInt(A, n, 4); 130dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 5); 131447fed29SStefano Zampini PetscValidLogicalCollectiveBool(A, At, 6); 132447fed29SStefano Zampini PetscValidLogicalCollectiveBool(B, Bt, 7); 1339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 1349566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 1359566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &cm, &cn)); 1369371c9d4SSatish Balay if (At) { 1379371c9d4SSatish Balay PetscInt tt = an; 1389371c9d4SSatish Balay an = am; 1399371c9d4SSatish Balay am = tt; 1409371c9d4SSatish Balay }; 1419371c9d4SSatish Balay if (Bt) { 1429371c9d4SSatish Balay PetscInt tt = bn; 1439371c9d4SSatish Balay bn = bm; 1449371c9d4SSatish Balay bm = tt; 1459371c9d4SSatish Balay }; 146aed4548fSBarry 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); 147447fed29SStefano Zampini 1484279555eSSatish Balay #if defined(PETSC_USE_INFO) 149447fed29SStefano Zampini sop = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)]; 1504279555eSSatish Balay #endif 1519566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)C), &rctx)); 1529566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 153447fed29SStefano Zampini if (Bt) { 1549566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &s1, &Bx)); 155447fed29SStefano Zampini } else { 1569566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &Bx, &s1)); 157447fed29SStefano Zampini } 158447fed29SStefano Zampini if (At) { 1599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &s2, &Ax)); 160447fed29SStefano Zampini } else { 1619566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Ax, &s2)); 162447fed29SStefano Zampini } 1639566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &Cx, &s3)); 164447fed29SStefano Zampini 165447fed29SStefano Zampini *flg = PETSC_TRUE; 166447fed29SStefano Zampini for (k = 0; k < n; k++) { 1679566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Bx, rctx)); 168447fed29SStefano Zampini if (Bt) { 1699566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Bx, s1)); 170447fed29SStefano Zampini } else { 1719566063dSJacob Faibussowitsch PetscCall(MatMult(B, Bx, s1)); 172447fed29SStefano Zampini } 1739566063dSJacob Faibussowitsch PetscCall(VecCopy(s1, Ax)); 174447fed29SStefano Zampini if (At) { 1759566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Ax, s2)); 176447fed29SStefano Zampini } else { 1779566063dSJacob Faibussowitsch PetscCall(MatMult(A, Ax, s2)); 178447fed29SStefano Zampini } 1799566063dSJacob Faibussowitsch PetscCall(VecCopy(Bx, Cx)); 1809566063dSJacob Faibussowitsch PetscCall(MatMult(C, Cx, s3)); 181447fed29SStefano Zampini 1829566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 183447fed29SStefano Zampini if (r2 < tol) { 1849566063dSJacob Faibussowitsch PetscCall(VecNorm(s3, NORM_INFINITY, &r1)); 185447fed29SStefano Zampini } else { 1869566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, none, s3)); 1879566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 188447fed29SStefano Zampini r1 /= r2; 189447fed29SStefano Zampini } 190447fed29SStefano Zampini if (r1 > tol) { 191447fed29SStefano Zampini *flg = PETSC_FALSE; 1929566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s %g\n", k, sop, (double)r1)); 193447fed29SStefano Zampini break; 194447fed29SStefano Zampini } 195447fed29SStefano Zampini } 1969566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 1999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Cx)); 2009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 2019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 2029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s3)); 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 204447fed29SStefano Zampini } 205447fed29SStefano Zampini 20686a22c91SHong Zhang /*@ 20786a22c91SHong Zhang MatMultEqual - Compares matrix-vector products of two matrices. 20886a22c91SHong Zhang 209c3339decSBarry Smith Collective 21086a22c91SHong Zhang 21186a22c91SHong Zhang Input Parameters: 21286a22c91SHong Zhang + A - the first matrix 2134222ddf1SHong Zhang . B - the second matrix 21486a22c91SHong Zhang - n - number of random vectors to be tested 21586a22c91SHong Zhang 21686a22c91SHong Zhang Output Parameter: 21720f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 21886a22c91SHong Zhang 21986a22c91SHong Zhang Level: intermediate 22086a22c91SHong Zhang 22120f4b53cSBarry Smith .seealso: `Mat`, `MatMultAddEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()`, `MatIsLinear()` 22286a22c91SHong Zhang @*/ 223d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 224d71ae5a4SJacob Faibussowitsch { 22586a22c91SHong Zhang PetscFunctionBegin; 22636eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 0)); 2273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22886a22c91SHong Zhang } 22986a22c91SHong Zhang 23086a22c91SHong Zhang /*@ 23120f4b53cSBarry Smith MatMultAddEqual - Compares matrix-vector product plus vector add of two matrices. 23286a22c91SHong Zhang 233c3339decSBarry Smith Collective 23486a22c91SHong Zhang 23586a22c91SHong Zhang Input Parameters: 23686a22c91SHong Zhang + A - the first matrix 2374222ddf1SHong Zhang . B - the second matrix 23886a22c91SHong Zhang - n - number of random vectors to be tested 23986a22c91SHong Zhang 24086a22c91SHong Zhang Output Parameter: 24120f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 24286a22c91SHong Zhang 24386a22c91SHong Zhang Level: intermediate 24486a22c91SHong Zhang 24520f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()` 24686a22c91SHong Zhang @*/ 247d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 248d71ae5a4SJacob Faibussowitsch { 24986a22c91SHong Zhang PetscFunctionBegin; 25036eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 1)); 25136eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 2)); 2523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25363879571SHong Zhang } 25463879571SHong Zhang 25563879571SHong Zhang /*@ 25663879571SHong Zhang MatMultTransposeEqual - Compares matrix-vector products of two matrices. 25763879571SHong Zhang 258c3339decSBarry Smith Collective 25963879571SHong Zhang 26063879571SHong Zhang Input Parameters: 26163879571SHong Zhang + A - the first matrix 2624222ddf1SHong Zhang . B - the second matrix 26363879571SHong Zhang - n - number of random vectors to be tested 26463879571SHong Zhang 26563879571SHong Zhang Output Parameter: 26620f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 26763879571SHong Zhang 26863879571SHong Zhang Level: intermediate 26963879571SHong Zhang 27020f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeAddEqual()` 27163879571SHong Zhang @*/ 272d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 273d71ae5a4SJacob Faibussowitsch { 27463879571SHong Zhang PetscFunctionBegin; 27536eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 0)); 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27763879571SHong Zhang } 27863879571SHong Zhang 27963879571SHong Zhang /*@ 28063879571SHong Zhang MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 28163879571SHong Zhang 282c3339decSBarry Smith Collective 28363879571SHong Zhang 28463879571SHong Zhang Input Parameters: 28563879571SHong Zhang + A - the first matrix 2864222ddf1SHong Zhang . B - the second matrix 28763879571SHong Zhang - n - number of random vectors to be tested 28863879571SHong Zhang 28963879571SHong Zhang Output Parameter: 29020f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 29163879571SHong Zhang 29263879571SHong Zhang Level: intermediate 29363879571SHong Zhang 29420f4b53cSBarry Smith .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 29563879571SHong Zhang @*/ 296d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 297d71ae5a4SJacob Faibussowitsch { 29863879571SHong Zhang PetscFunctionBegin; 29936eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 1)); 30036eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 2)); 3013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 302e573050aSPierre Jolivet } 303e573050aSPierre Jolivet 304e573050aSPierre Jolivet /*@ 305e573050aSPierre Jolivet MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices. 306e573050aSPierre Jolivet 307c3339decSBarry Smith Collective 308e573050aSPierre Jolivet 309e573050aSPierre Jolivet Input Parameters: 310e573050aSPierre Jolivet + A - the first matrix 311e573050aSPierre Jolivet . B - the second matrix 312e573050aSPierre Jolivet - n - number of random vectors to be tested 313e573050aSPierre Jolivet 314e573050aSPierre Jolivet Output Parameter: 31520f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 316e573050aSPierre Jolivet 317e573050aSPierre Jolivet Level: intermediate 318e573050aSPierre Jolivet 319*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 320e573050aSPierre Jolivet @*/ 321d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 322d71ae5a4SJacob Faibussowitsch { 323e573050aSPierre Jolivet PetscFunctionBegin; 32436eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 0)); 3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 326e573050aSPierre Jolivet } 327e573050aSPierre Jolivet 328e573050aSPierre Jolivet /*@ 329e573050aSPierre Jolivet MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices. 330e573050aSPierre Jolivet 331c3339decSBarry Smith Collective 332e573050aSPierre Jolivet 333e573050aSPierre Jolivet Input Parameters: 334e573050aSPierre Jolivet + A - the first matrix 335e573050aSPierre Jolivet . B - the second matrix 336e573050aSPierre Jolivet - n - number of random vectors to be tested 337e573050aSPierre Jolivet 338e573050aSPierre Jolivet Output Parameter: 33920f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 340e573050aSPierre Jolivet 341e573050aSPierre Jolivet Level: intermediate 342e573050aSPierre Jolivet 343*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 344e573050aSPierre Jolivet @*/ 345d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 346d71ae5a4SJacob Faibussowitsch { 347e573050aSPierre Jolivet PetscFunctionBegin; 34836eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 1)); 34936eb99e2SToby Isaac PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 2)); 3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35186a22c91SHong Zhang } 352a52676ecSHong Zhang 353a52676ecSHong Zhang /*@ 354a52676ecSHong Zhang MatMatMultEqual - Test A*B*x = C*x for n random vector x 355a52676ecSHong Zhang 356c3339decSBarry Smith Collective 357a52676ecSHong Zhang 358a52676ecSHong Zhang Input Parameters: 359a52676ecSHong Zhang + A - the first matrix 360c2916339SPierre Jolivet . B - the second matrix 361c2916339SPierre Jolivet . C - the third matrix 362a52676ecSHong Zhang - n - number of random vectors to be tested 363a52676ecSHong Zhang 364a52676ecSHong Zhang Output Parameter: 36520f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 366a52676ecSHong Zhang 367a52676ecSHong Zhang Level: intermediate 368a52676ecSHong Zhang 369*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 370a52676ecSHong Zhang @*/ 371d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 372d71ae5a4SJacob Faibussowitsch { 373a52676ecSHong Zhang PetscFunctionBegin; 3749566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_FALSE)); 3753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 376a52676ecSHong Zhang } 377a52676ecSHong Zhang 378a52676ecSHong Zhang /*@ 379a52676ecSHong Zhang MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x 380a52676ecSHong Zhang 381c3339decSBarry Smith Collective 382a52676ecSHong Zhang 383a52676ecSHong Zhang Input Parameters: 384a52676ecSHong Zhang + A - the first matrix 3854222ddf1SHong Zhang . B - the second matrix 3864222ddf1SHong Zhang . C - the third matrix 387a52676ecSHong Zhang - n - number of random vectors to be tested 388a52676ecSHong Zhang 389a52676ecSHong Zhang Output Parameter: 39020f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 391a52676ecSHong Zhang 392a52676ecSHong Zhang Level: intermediate 393a52676ecSHong Zhang 394*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 395a52676ecSHong Zhang @*/ 396d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTransposeMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 397d71ae5a4SJacob Faibussowitsch { 398a52676ecSHong Zhang PetscFunctionBegin; 3999566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_TRUE, PETSC_FALSE)); 4003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 401a52676ecSHong Zhang } 40286919fd6SHong Zhang 40326546c1bSToby Isaac /*@ 40426546c1bSToby Isaac MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x 40526546c1bSToby Isaac 406c3339decSBarry Smith Collective 40726546c1bSToby Isaac 40826546c1bSToby Isaac Input Parameters: 40926546c1bSToby Isaac + A - the first matrix 4104222ddf1SHong Zhang . B - the second matrix 4114222ddf1SHong Zhang . C - the third matrix 41226546c1bSToby Isaac - n - number of random vectors to be tested 41326546c1bSToby Isaac 41426546c1bSToby Isaac Output Parameter: 41520f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 41626546c1bSToby Isaac 41726546c1bSToby Isaac Level: intermediate 41826546c1bSToby Isaac 419*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 42026546c1bSToby Isaac @*/ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatTransposeMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 422d71ae5a4SJacob Faibussowitsch { 423447fed29SStefano Zampini PetscFunctionBegin; 4249566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_TRUE)); 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 426447fed29SStefano Zampini } 427447fed29SStefano Zampini 428d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProjMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool rart, PetscBool *flg) 429d71ae5a4SJacob Faibussowitsch { 430447fed29SStefano Zampini Vec x, v1, v2, v3, v4, Cx, Bx; 431447fed29SStefano Zampini PetscReal norm_abs, norm_rel, tol = PETSC_SQRT_MACHINE_EPSILON; 432447fed29SStefano Zampini PetscInt i, am, an, bm, bn, cm, cn; 433447fed29SStefano Zampini PetscRandom rdm; 434cc48ffa7SToby Isaac PetscScalar none = -1.0; 435cc48ffa7SToby Isaac 436cc48ffa7SToby Isaac PetscFunctionBegin; 4379566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 4389566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 4399371c9d4SSatish Balay if (rart) { 4409371c9d4SSatish Balay PetscInt t = bm; 4419371c9d4SSatish Balay bm = bn; 4429371c9d4SSatish Balay bn = t; 4439371c9d4SSatish Balay } 4449566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &cm, &cn)); 445aed4548fSBarry 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); 446cc48ffa7SToby Isaac 447447fed29SStefano Zampini /* Create left vector of A: v2 */ 4489566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Bx, &v2)); 449447fed29SStefano Zampini 450447fed29SStefano Zampini /* Create right vectors of B: x, v3, v4 */ 451447fed29SStefano Zampini if (rart) { 4529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &v1, &x)); 453447fed29SStefano Zampini } else { 4549566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &x, &v1)); 455447fed29SStefano Zampini } 4569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &v3)); 457447fed29SStefano Zampini 4589566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &Cx, &v4)); 4599566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 4609566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 461cc48ffa7SToby Isaac 462cc48ffa7SToby Isaac *flg = PETSC_TRUE; 463447fed29SStefano Zampini for (i = 0; i < n; i++) { 4649566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 4659566063dSJacob Faibussowitsch PetscCall(VecCopy(x, Cx)); 4669566063dSJacob Faibussowitsch PetscCall(MatMult(C, Cx, v4)); /* v4 = C*x */ 467447fed29SStefano Zampini if (rart) { 4689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, x, v1)); 469cc48ffa7SToby Isaac } else { 4709566063dSJacob Faibussowitsch PetscCall(MatMult(B, x, v1)); 471cc48ffa7SToby Isaac } 4729566063dSJacob Faibussowitsch PetscCall(VecCopy(v1, Bx)); 4739566063dSJacob Faibussowitsch PetscCall(MatMult(A, Bx, v2)); /* v2 = A*B*x */ 4749566063dSJacob Faibussowitsch PetscCall(VecCopy(v2, v1)); 475447fed29SStefano Zampini if (rart) { 4769566063dSJacob Faibussowitsch PetscCall(MatMult(B, v1, v3)); /* v3 = R*A*R^t*x */ 477447fed29SStefano Zampini } else { 4789566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, v1, v3)); /* v3 = Bt*A*B*x */ 479447fed29SStefano Zampini } 4809566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_abs)); 4819566063dSJacob Faibussowitsch PetscCall(VecAXPY(v4, none, v3)); 4829566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_rel)); 483447fed29SStefano Zampini 484447fed29SStefano Zampini if (norm_abs > tol) norm_rel /= norm_abs; 485447fed29SStefano Zampini if (norm_rel > tol) { 486cc48ffa7SToby Isaac *flg = PETSC_FALSE; 4879566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th Mat%sMult() %g\n", i, rart ? "RARt" : "PtAP", (double)norm_rel)); 488cc48ffa7SToby Isaac break; 489cc48ffa7SToby Isaac } 490cc48ffa7SToby Isaac } 491447fed29SStefano Zampini 4929566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 4959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Cx)); 4969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v1)); 4979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v2)); 4989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v3)); 4999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v4)); 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 501cc48ffa7SToby Isaac } 502cc48ffa7SToby Isaac 50386919fd6SHong Zhang /*@ 5044222ddf1SHong Zhang MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B 5054222ddf1SHong Zhang 506c3339decSBarry Smith Collective 5074222ddf1SHong Zhang 5084222ddf1SHong Zhang Input Parameters: 5094222ddf1SHong Zhang + A - the first matrix 5104222ddf1SHong Zhang . B - the second matrix 5114222ddf1SHong Zhang . C - the third matrix 5124222ddf1SHong Zhang - n - number of random vectors to be tested 5134222ddf1SHong Zhang 5144222ddf1SHong Zhang Output Parameter: 51520f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 5164222ddf1SHong Zhang 5174222ddf1SHong Zhang Level: intermediate 5184222ddf1SHong Zhang 519*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 5204222ddf1SHong Zhang @*/ 521d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 522d71ae5a4SJacob Faibussowitsch { 5234222ddf1SHong Zhang PetscFunctionBegin; 5249566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_FALSE, flg)); 5253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5264222ddf1SHong Zhang } 5274222ddf1SHong Zhang 528447fed29SStefano Zampini /*@ 529447fed29SStefano Zampini MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t 530447fed29SStefano Zampini 531c3339decSBarry Smith Collective 532447fed29SStefano Zampini 533447fed29SStefano Zampini Input Parameters: 534447fed29SStefano Zampini + A - the first matrix 535447fed29SStefano Zampini . B - the second matrix 536447fed29SStefano Zampini . C - the third matrix 537447fed29SStefano Zampini - n - number of random vectors to be tested 538447fed29SStefano Zampini 539447fed29SStefano Zampini Output Parameter: 54020f4b53cSBarry Smith . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 541447fed29SStefano Zampini 542447fed29SStefano Zampini Level: intermediate 543447fed29SStefano Zampini 544*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 545447fed29SStefano Zampini @*/ 546d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRARtMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 547d71ae5a4SJacob Faibussowitsch { 548447fed29SStefano Zampini PetscFunctionBegin; 5499566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_TRUE, flg)); 5503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5514222ddf1SHong Zhang } 5524222ddf1SHong Zhang 5534222ddf1SHong Zhang /*@ 55420f4b53cSBarry Smith MatIsLinear - Check if a shell matrix `A` is a linear operator. 55586919fd6SHong Zhang 556c3339decSBarry Smith Collective 55786919fd6SHong Zhang 55886919fd6SHong Zhang Input Parameters: 55986919fd6SHong Zhang + A - the shell matrix 56086919fd6SHong Zhang - n - number of random vectors to be tested 56186919fd6SHong Zhang 56286919fd6SHong Zhang Output Parameter: 56320f4b53cSBarry Smith . flg - `PETSC_TRUE` if the shell matrix is linear; `PETSC_FALSE` otherwise. 56486919fd6SHong Zhang 56586919fd6SHong Zhang Level: intermediate 56620f4b53cSBarry Smith 567*52cd20c4SPierre Jolivet .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 56886919fd6SHong Zhang @*/ 569d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsLinear(Mat A, PetscInt n, PetscBool *flg) 570d71ae5a4SJacob Faibussowitsch { 57186919fd6SHong Zhang Vec x, y, s1, s2; 57286919fd6SHong Zhang PetscRandom rctx; 57386919fd6SHong Zhang PetscScalar a; 57486919fd6SHong Zhang PetscInt k; 57586919fd6SHong Zhang PetscReal norm, normA; 57686919fd6SHong Zhang MPI_Comm comm; 57786919fd6SHong Zhang PetscMPIInt rank; 57886919fd6SHong Zhang 57986919fd6SHong Zhang PetscFunctionBegin; 58086919fd6SHong Zhang PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5819566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 5829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 58386919fd6SHong Zhang 5849566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rctx)); 5859566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 5869566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &s1)); 5879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 5889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s1, &s2)); 58986919fd6SHong Zhang 59086919fd6SHong Zhang *flg = PETSC_TRUE; 59186919fd6SHong Zhang for (k = 0; k < n; k++) { 5929566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 5939566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 59448a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscRandomGetValue(rctx, &a)); 5959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm)); 59686919fd6SHong Zhang 59786919fd6SHong Zhang /* s2 = a*A*x + A*y */ 5989566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s2)); /* s2 = A*y */ 5999566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, s1)); /* s1 = A*x */ 6009566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, a, s1)); /* s2 = a s1 + s2 */ 60186919fd6SHong Zhang 60286919fd6SHong Zhang /* s1 = A * (a x + y) */ 6039566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, a, x)); /* y = a x + y */ 6049566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s1)); 6059566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_INFINITY, &normA)); 60686919fd6SHong Zhang 6079566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, -1.0, s1)); /* s2 = - s1 + s2 */ 6089566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &norm)); 6091b8dac88SHong Zhang if (norm / normA > 100. * PETSC_MACHINE_EPSILON) { 61086919fd6SHong Zhang *flg = PETSC_FALSE; 6119566063dSJacob 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))); 61286919fd6SHong Zhang break; 61386919fd6SHong Zhang } 61486919fd6SHong Zhang } 6159566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 6169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 6179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 6189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 6199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 6203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62186919fd6SHong Zhang } 622