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; 11580c7c76SPierre Jolivet const char *sops[] = {"MatMult", "MatMultAdd", "MatMultTranspose", "MatMultTransposeAdd", "MatMultHermitianTranspose", "MatMultHermitianTransposeAdd"}; 12447fed29SStefano Zampini const char *sop; 13447fed29SStefano Zampini 14447fed29SStefano Zampini PetscFunctionBegin; 15447fed29SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 16447fed29SStefano Zampini PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 17447fed29SStefano Zampini PetscCheckSameComm(A, 1, B, 2); 18447fed29SStefano Zampini PetscValidLogicalCollectiveInt(A, n, 3); 19dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 4); 20e573050aSPierre Jolivet PetscValidLogicalCollectiveInt(A, t, 5); 21447fed29SStefano Zampini PetscValidLogicalCollectiveBool(A, add, 6); 229566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 239566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 24aed4548fSBarry 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); 25e573050aSPierre Jolivet sop = sops[(add ? 1 : 0) + 2 * t]; /* t = 0 => no transpose, t = 1 => transpose, t = 2 => Hermitian transpose */ 269566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)A), &rctx)); 279566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 28447fed29SStefano Zampini if (t) { 299566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &s1, &Ax)); 309566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &s2, &Bx)); 31447fed29SStefano Zampini } else { 329566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Ax, &s1)); 339566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &Bx, &s2)); 34447fed29SStefano Zampini } 35447fed29SStefano Zampini if (add) { 369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s1, &Ay)); 379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s2, &By)); 38447fed29SStefano Zampini } 39447fed29SStefano Zampini 40447fed29SStefano Zampini *flg = PETSC_TRUE; 41447fed29SStefano Zampini for (k = 0; k < n; k++) { 429566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ax, rctx)); 439566063dSJacob Faibussowitsch PetscCall(VecCopy(Ax, Bx)); 44447fed29SStefano Zampini if (add) { 459566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Ay, rctx)); 469566063dSJacob Faibussowitsch PetscCall(VecCopy(Ay, By)); 47447fed29SStefano Zampini } 48e573050aSPierre Jolivet if (t == 1) { 49447fed29SStefano Zampini if (add) { 509566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A, Ax, Ay, s1)); 519566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(B, Bx, By, s2)); 52447fed29SStefano Zampini } else { 539566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Ax, s1)); 549566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Bx, s2)); 55447fed29SStefano Zampini } 56e573050aSPierre Jolivet } else if (t == 2) { 57e573050aSPierre Jolivet if (add) { 589566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(A, Ax, Ay, s1)); 599566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTransposeAdd(B, Bx, By, s2)); 60e573050aSPierre Jolivet } else { 619566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(A, Ax, s1)); 629566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(B, Bx, s2)); 63e573050aSPierre Jolivet } 64447fed29SStefano Zampini } else { 65447fed29SStefano Zampini if (add) { 669566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A, Ax, Ay, s1)); 679566063dSJacob Faibussowitsch PetscCall(MatMultAdd(B, Bx, By, s2)); 68447fed29SStefano Zampini } else { 699566063dSJacob Faibussowitsch PetscCall(MatMult(A, Ax, s1)); 709566063dSJacob Faibussowitsch PetscCall(MatMult(B, Bx, s2)); 71447fed29SStefano Zampini } 72447fed29SStefano Zampini } 739566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 74447fed29SStefano Zampini if (r2 < tol) { 759566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_INFINITY, &r1)); 76447fed29SStefano Zampini } else { 779566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, none, s1)); 789566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 79447fed29SStefano Zampini r1 /= r2; 80447fed29SStefano Zampini } 81447fed29SStefano Zampini if (r1 > tol) { 82447fed29SStefano Zampini *flg = PETSC_FALSE; 839566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s() %g\n", k, sop, (double)r1)); 84447fed29SStefano Zampini break; 85447fed29SStefano Zampini } 86447fed29SStefano Zampini } 879566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ay)); 919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&By)); 929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 94447fed29SStefano Zampini PetscFunctionReturn(0); 95447fed29SStefano Zampini } 96447fed29SStefano Zampini 97d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg, PetscBool At, PetscBool Bt) 98d71ae5a4SJacob Faibussowitsch { 99447fed29SStefano Zampini Vec Ax, Bx, Cx, s1, s2, s3; 100447fed29SStefano Zampini PetscRandom rctx; 101447fed29SStefano Zampini PetscReal r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON; 102447fed29SStefano Zampini PetscInt am, an, bm, bn, cm, cn, k; 103447fed29SStefano Zampini PetscScalar none = -1.0; 104580c7c76SPierre Jolivet const char *sops[] = {"MatMatMult", "MatTransposeMatMult", "MatMatTransposeMult", "MatTransposeMatTransposeMult"}; 105447fed29SStefano Zampini const char *sop; 106447fed29SStefano Zampini 107447fed29SStefano Zampini PetscFunctionBegin; 108447fed29SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 109447fed29SStefano Zampini PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 110447fed29SStefano Zampini PetscCheckSameComm(A, 1, B, 2); 111447fed29SStefano Zampini PetscValidHeaderSpecific(C, MAT_CLASSID, 3); 112447fed29SStefano Zampini PetscCheckSameComm(A, 1, C, 3); 113447fed29SStefano Zampini PetscValidLogicalCollectiveInt(A, n, 4); 114dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 5); 115447fed29SStefano Zampini PetscValidLogicalCollectiveBool(A, At, 6); 116447fed29SStefano Zampini PetscValidLogicalCollectiveBool(B, Bt, 7); 1179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 1189566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 1199566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &cm, &cn)); 1209371c9d4SSatish Balay if (At) { 1219371c9d4SSatish Balay PetscInt tt = an; 1229371c9d4SSatish Balay an = am; 1239371c9d4SSatish Balay am = tt; 1249371c9d4SSatish Balay }; 1259371c9d4SSatish Balay if (Bt) { 1269371c9d4SSatish Balay PetscInt tt = bn; 1279371c9d4SSatish Balay bn = bm; 1289371c9d4SSatish Balay bm = tt; 1299371c9d4SSatish Balay }; 130aed4548fSBarry 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); 131447fed29SStefano Zampini 132447fed29SStefano Zampini sop = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)]; 1339566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)C), &rctx)); 1349566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 135447fed29SStefano Zampini if (Bt) { 1369566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &s1, &Bx)); 137447fed29SStefano Zampini } else { 1389566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &Bx, &s1)); 139447fed29SStefano Zampini } 140447fed29SStefano Zampini if (At) { 1419566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &s2, &Ax)); 142447fed29SStefano Zampini } else { 1439566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Ax, &s2)); 144447fed29SStefano Zampini } 1459566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &Cx, &s3)); 146447fed29SStefano Zampini 147447fed29SStefano Zampini *flg = PETSC_TRUE; 148447fed29SStefano Zampini for (k = 0; k < n; k++) { 1499566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Bx, rctx)); 150447fed29SStefano Zampini if (Bt) { 1519566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, Bx, s1)); 152447fed29SStefano Zampini } else { 1539566063dSJacob Faibussowitsch PetscCall(MatMult(B, Bx, s1)); 154447fed29SStefano Zampini } 1559566063dSJacob Faibussowitsch PetscCall(VecCopy(s1, Ax)); 156447fed29SStefano Zampini if (At) { 1579566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, Ax, s2)); 158447fed29SStefano Zampini } else { 1599566063dSJacob Faibussowitsch PetscCall(MatMult(A, Ax, s2)); 160447fed29SStefano Zampini } 1619566063dSJacob Faibussowitsch PetscCall(VecCopy(Bx, Cx)); 1629566063dSJacob Faibussowitsch PetscCall(MatMult(C, Cx, s3)); 163447fed29SStefano Zampini 1649566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 165447fed29SStefano Zampini if (r2 < tol) { 1669566063dSJacob Faibussowitsch PetscCall(VecNorm(s3, NORM_INFINITY, &r1)); 167447fed29SStefano Zampini } else { 1689566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, none, s3)); 1699566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 170447fed29SStefano Zampini r1 /= r2; 171447fed29SStefano Zampini } 172447fed29SStefano Zampini if (r1 > tol) { 173447fed29SStefano Zampini *flg = PETSC_FALSE; 1749566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s %g\n", k, sop, (double)r1)); 175447fed29SStefano Zampini break; 176447fed29SStefano Zampini } 177447fed29SStefano Zampini } 1789566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 1819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Cx)); 1829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s3)); 185447fed29SStefano Zampini PetscFunctionReturn(0); 186447fed29SStefano Zampini } 187447fed29SStefano Zampini 18886a22c91SHong Zhang /*@ 18986a22c91SHong Zhang MatMultEqual - Compares matrix-vector products of two matrices. 19086a22c91SHong Zhang 191*c3339decSBarry Smith Collective 19286a22c91SHong Zhang 19386a22c91SHong Zhang Input Parameters: 19486a22c91SHong Zhang + A - the first matrix 1954222ddf1SHong Zhang . B - the second matrix 19686a22c91SHong Zhang - n - number of random vectors to be tested 19786a22c91SHong Zhang 19886a22c91SHong Zhang Output Parameter: 19986a22c91SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 20086a22c91SHong Zhang 20186a22c91SHong Zhang Level: intermediate 20286a22c91SHong Zhang 20386a22c91SHong Zhang @*/ 204d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 205d71ae5a4SJacob Faibussowitsch { 20686a22c91SHong Zhang PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_FALSE)); 20886a22c91SHong Zhang PetscFunctionReturn(0); 20986a22c91SHong Zhang } 21086a22c91SHong Zhang 21186a22c91SHong Zhang /*@ 21286a22c91SHong Zhang MatMultAddEqual - Compares matrix-vector products of two matrices. 21386a22c91SHong Zhang 214*c3339decSBarry 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: 22286a22c91SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 22386a22c91SHong Zhang 22486a22c91SHong Zhang Level: intermediate 22586a22c91SHong Zhang 22686a22c91SHong Zhang @*/ 227d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 228d71ae5a4SJacob Faibussowitsch { 22986a22c91SHong Zhang PetscFunctionBegin; 2309566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 0, PETSC_TRUE)); 23163879571SHong Zhang PetscFunctionReturn(0); 23263879571SHong Zhang } 23363879571SHong Zhang 23463879571SHong Zhang /*@ 23563879571SHong Zhang MatMultTransposeEqual - Compares matrix-vector products of two matrices. 23663879571SHong Zhang 237*c3339decSBarry Smith Collective 23863879571SHong Zhang 23963879571SHong Zhang Input Parameters: 24063879571SHong Zhang + A - the first matrix 2414222ddf1SHong Zhang . B - the second matrix 24263879571SHong Zhang - n - number of random vectors to be tested 24363879571SHong Zhang 24463879571SHong Zhang Output Parameter: 24563879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 24663879571SHong Zhang 24763879571SHong Zhang Level: intermediate 24863879571SHong Zhang 24963879571SHong Zhang @*/ 250d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 251d71ae5a4SJacob Faibussowitsch { 25263879571SHong Zhang PetscFunctionBegin; 2539566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_FALSE)); 25463879571SHong Zhang PetscFunctionReturn(0); 25563879571SHong Zhang } 25663879571SHong Zhang 25763879571SHong Zhang /*@ 25863879571SHong Zhang MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 25963879571SHong Zhang 260*c3339decSBarry Smith Collective 26163879571SHong Zhang 26263879571SHong Zhang Input Parameters: 26363879571SHong Zhang + A - the first matrix 2644222ddf1SHong Zhang . B - the second matrix 26563879571SHong Zhang - n - number of random vectors to be tested 26663879571SHong Zhang 26763879571SHong Zhang Output Parameter: 26863879571SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 26963879571SHong Zhang 27063879571SHong Zhang Level: intermediate 27163879571SHong Zhang 27263879571SHong Zhang @*/ 273d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 274d71ae5a4SJacob Faibussowitsch { 27563879571SHong Zhang PetscFunctionBegin; 2769566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 1, PETSC_TRUE)); 277e573050aSPierre Jolivet PetscFunctionReturn(0); 278e573050aSPierre Jolivet } 279e573050aSPierre Jolivet 280e573050aSPierre Jolivet /*@ 281e573050aSPierre Jolivet MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices. 282e573050aSPierre Jolivet 283*c3339decSBarry Smith Collective 284e573050aSPierre Jolivet 285e573050aSPierre Jolivet Input Parameters: 286e573050aSPierre Jolivet + A - the first matrix 287e573050aSPierre Jolivet . B - the second matrix 288e573050aSPierre Jolivet - n - number of random vectors to be tested 289e573050aSPierre Jolivet 290e573050aSPierre Jolivet Output Parameter: 291e573050aSPierre Jolivet . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 292e573050aSPierre Jolivet 293e573050aSPierre Jolivet Level: intermediate 294e573050aSPierre Jolivet 295e573050aSPierre Jolivet @*/ 296d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 297d71ae5a4SJacob Faibussowitsch { 298e573050aSPierre Jolivet PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 2, PETSC_FALSE)); 300e573050aSPierre Jolivet PetscFunctionReturn(0); 301e573050aSPierre Jolivet } 302e573050aSPierre Jolivet 303e573050aSPierre Jolivet /*@ 304e573050aSPierre Jolivet MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices. 305e573050aSPierre Jolivet 306*c3339decSBarry Smith Collective 307e573050aSPierre Jolivet 308e573050aSPierre Jolivet Input Parameters: 309e573050aSPierre Jolivet + A - the first matrix 310e573050aSPierre Jolivet . B - the second matrix 311e573050aSPierre Jolivet - n - number of random vectors to be tested 312e573050aSPierre Jolivet 313e573050aSPierre Jolivet Output Parameter: 314e573050aSPierre Jolivet . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 315e573050aSPierre Jolivet 316e573050aSPierre Jolivet Level: intermediate 317e573050aSPierre Jolivet 318e573050aSPierre Jolivet @*/ 319d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultHermitianTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 320d71ae5a4SJacob Faibussowitsch { 321e573050aSPierre Jolivet PetscFunctionBegin; 3229566063dSJacob Faibussowitsch PetscCall(MatMultEqual_Private(A, B, n, flg, 2, PETSC_TRUE)); 32386a22c91SHong Zhang PetscFunctionReturn(0); 32486a22c91SHong Zhang } 325a52676ecSHong Zhang 326a52676ecSHong Zhang /*@ 327a52676ecSHong Zhang MatMatMultEqual - Test A*B*x = C*x for n random vector x 328a52676ecSHong Zhang 329*c3339decSBarry Smith Collective 330a52676ecSHong Zhang 331a52676ecSHong Zhang Input Parameters: 332a52676ecSHong Zhang + A - the first matrix 333c2916339SPierre Jolivet . B - the second matrix 334c2916339SPierre Jolivet . C - the third matrix 335a52676ecSHong Zhang - n - number of random vectors to be tested 336a52676ecSHong Zhang 337a52676ecSHong Zhang Output Parameter: 338f0fc11ceSJed Brown . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 339a52676ecSHong Zhang 340a52676ecSHong Zhang Level: intermediate 341a52676ecSHong Zhang 342a52676ecSHong Zhang @*/ 343d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 344d71ae5a4SJacob Faibussowitsch { 345a52676ecSHong Zhang PetscFunctionBegin; 3469566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_FALSE)); 347a52676ecSHong Zhang PetscFunctionReturn(0); 348a52676ecSHong Zhang } 349a52676ecSHong Zhang 350a52676ecSHong Zhang /*@ 351a52676ecSHong Zhang MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x 352a52676ecSHong Zhang 353*c3339decSBarry Smith Collective 354a52676ecSHong Zhang 355a52676ecSHong Zhang Input Parameters: 356a52676ecSHong Zhang + A - the first matrix 3574222ddf1SHong Zhang . B - the second matrix 3584222ddf1SHong Zhang . C - the third matrix 359a52676ecSHong Zhang - n - number of random vectors to be tested 360a52676ecSHong Zhang 361a52676ecSHong Zhang Output Parameter: 362a52676ecSHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 363a52676ecSHong Zhang 364a52676ecSHong Zhang Level: intermediate 365a52676ecSHong Zhang 366a52676ecSHong Zhang @*/ 367d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTransposeMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 368d71ae5a4SJacob Faibussowitsch { 369a52676ecSHong Zhang PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_TRUE, PETSC_FALSE)); 371a52676ecSHong Zhang PetscFunctionReturn(0); 372a52676ecSHong Zhang } 37386919fd6SHong Zhang 37426546c1bSToby Isaac /*@ 37526546c1bSToby Isaac MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x 37626546c1bSToby Isaac 377*c3339decSBarry Smith Collective 37826546c1bSToby Isaac 37926546c1bSToby Isaac Input Parameters: 38026546c1bSToby Isaac + A - the first matrix 3814222ddf1SHong Zhang . B - the second matrix 3824222ddf1SHong Zhang . C - the third matrix 38326546c1bSToby Isaac - n - number of random vectors to be tested 38426546c1bSToby Isaac 38526546c1bSToby Isaac Output Parameter: 38626546c1bSToby Isaac . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 38726546c1bSToby Isaac 38826546c1bSToby Isaac Level: intermediate 38926546c1bSToby Isaac 39026546c1bSToby Isaac @*/ 391d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMatTransposeMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 392d71ae5a4SJacob Faibussowitsch { 393447fed29SStefano Zampini PetscFunctionBegin; 3949566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_TRUE)); 395447fed29SStefano Zampini PetscFunctionReturn(0); 396447fed29SStefano Zampini } 397447fed29SStefano Zampini 398d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProjMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool rart, PetscBool *flg) 399d71ae5a4SJacob Faibussowitsch { 400447fed29SStefano Zampini Vec x, v1, v2, v3, v4, Cx, Bx; 401447fed29SStefano Zampini PetscReal norm_abs, norm_rel, tol = PETSC_SQRT_MACHINE_EPSILON; 402447fed29SStefano Zampini PetscInt i, am, an, bm, bn, cm, cn; 403447fed29SStefano Zampini PetscRandom rdm; 404cc48ffa7SToby Isaac PetscScalar none = -1.0; 405cc48ffa7SToby Isaac 406cc48ffa7SToby Isaac PetscFunctionBegin; 4079566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 4089566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B, &bm, &bn)); 4099371c9d4SSatish Balay if (rart) { 4109371c9d4SSatish Balay PetscInt t = bm; 4119371c9d4SSatish Balay bm = bn; 4129371c9d4SSatish Balay bn = t; 4139371c9d4SSatish Balay } 4149566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C, &cm, &cn)); 415aed4548fSBarry 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); 416cc48ffa7SToby Isaac 417447fed29SStefano Zampini /* Create left vector of A: v2 */ 4189566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Bx, &v2)); 419447fed29SStefano Zampini 420447fed29SStefano Zampini /* Create right vectors of B: x, v3, v4 */ 421447fed29SStefano Zampini if (rart) { 4229566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &v1, &x)); 423447fed29SStefano Zampini } else { 4249566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &x, &v1)); 425447fed29SStefano Zampini } 4269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &v3)); 427447fed29SStefano Zampini 4289566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &Cx, &v4)); 4299566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 4309566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 431cc48ffa7SToby Isaac 432cc48ffa7SToby Isaac *flg = PETSC_TRUE; 433447fed29SStefano Zampini for (i = 0; i < n; i++) { 4349566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 4359566063dSJacob Faibussowitsch PetscCall(VecCopy(x, Cx)); 4369566063dSJacob Faibussowitsch PetscCall(MatMult(C, Cx, v4)); /* v4 = C*x */ 437447fed29SStefano Zampini if (rart) { 4389566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, x, v1)); 439cc48ffa7SToby Isaac } else { 4409566063dSJacob Faibussowitsch PetscCall(MatMult(B, x, v1)); 441cc48ffa7SToby Isaac } 4429566063dSJacob Faibussowitsch PetscCall(VecCopy(v1, Bx)); 4439566063dSJacob Faibussowitsch PetscCall(MatMult(A, Bx, v2)); /* v2 = A*B*x */ 4449566063dSJacob Faibussowitsch PetscCall(VecCopy(v2, v1)); 445447fed29SStefano Zampini if (rart) { 4469566063dSJacob Faibussowitsch PetscCall(MatMult(B, v1, v3)); /* v3 = R*A*R^t*x */ 447447fed29SStefano Zampini } else { 4489566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(B, v1, v3)); /* v3 = Bt*A*B*x */ 449447fed29SStefano Zampini } 4509566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_abs)); 4519566063dSJacob Faibussowitsch PetscCall(VecAXPY(v4, none, v3)); 4529566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_rel)); 453447fed29SStefano Zampini 454447fed29SStefano Zampini if (norm_abs > tol) norm_rel /= norm_abs; 455447fed29SStefano Zampini if (norm_rel > tol) { 456cc48ffa7SToby Isaac *flg = PETSC_FALSE; 4579566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th Mat%sMult() %g\n", i, rart ? "RARt" : "PtAP", (double)norm_rel)); 458cc48ffa7SToby Isaac break; 459cc48ffa7SToby Isaac } 460cc48ffa7SToby Isaac } 461447fed29SStefano Zampini 4629566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 4639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 4649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Bx)); 4659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Cx)); 4669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v1)); 4679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v2)); 4689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v3)); 4699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v4)); 470cc48ffa7SToby Isaac PetscFunctionReturn(0); 471cc48ffa7SToby Isaac } 472cc48ffa7SToby Isaac 47386919fd6SHong Zhang /*@ 4744222ddf1SHong Zhang MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B 4754222ddf1SHong Zhang 476*c3339decSBarry Smith Collective 4774222ddf1SHong Zhang 4784222ddf1SHong Zhang Input Parameters: 4794222ddf1SHong Zhang + A - the first matrix 4804222ddf1SHong Zhang . B - the second matrix 4814222ddf1SHong Zhang . C - the third matrix 4824222ddf1SHong Zhang - n - number of random vectors to be tested 4834222ddf1SHong Zhang 4844222ddf1SHong Zhang Output Parameter: 4854222ddf1SHong Zhang . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 4864222ddf1SHong Zhang 4874222ddf1SHong Zhang Level: intermediate 4884222ddf1SHong Zhang 4894222ddf1SHong Zhang @*/ 490d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 491d71ae5a4SJacob Faibussowitsch { 4924222ddf1SHong Zhang PetscFunctionBegin; 4939566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_FALSE, flg)); 494447fed29SStefano Zampini PetscFunctionReturn(0); 4954222ddf1SHong Zhang } 4964222ddf1SHong Zhang 497447fed29SStefano Zampini /*@ 498447fed29SStefano Zampini MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t 499447fed29SStefano Zampini 500*c3339decSBarry Smith Collective 501447fed29SStefano Zampini 502447fed29SStefano Zampini Input Parameters: 503447fed29SStefano Zampini + A - the first matrix 504447fed29SStefano Zampini . B - the second matrix 505447fed29SStefano Zampini . C - the third matrix 506447fed29SStefano Zampini - n - number of random vectors to be tested 507447fed29SStefano Zampini 508447fed29SStefano Zampini Output Parameter: 509447fed29SStefano Zampini . flg - PETSC_TRUE if the products are equal; PETSC_FALSE otherwise. 510447fed29SStefano Zampini 511447fed29SStefano Zampini Level: intermediate 512447fed29SStefano Zampini 513447fed29SStefano Zampini @*/ 514d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRARtMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 515d71ae5a4SJacob Faibussowitsch { 516447fed29SStefano Zampini PetscFunctionBegin; 5179566063dSJacob Faibussowitsch PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_TRUE, flg)); 5184222ddf1SHong Zhang PetscFunctionReturn(0); 5194222ddf1SHong Zhang } 5204222ddf1SHong Zhang 5214222ddf1SHong Zhang /*@ 52286919fd6SHong Zhang MatIsLinear - Check if a shell matrix A is a linear operator. 52386919fd6SHong Zhang 524*c3339decSBarry Smith Collective 52586919fd6SHong Zhang 52686919fd6SHong Zhang Input Parameters: 52786919fd6SHong Zhang + A - the shell matrix 52886919fd6SHong Zhang - n - number of random vectors to be tested 52986919fd6SHong Zhang 53086919fd6SHong Zhang Output Parameter: 53186919fd6SHong Zhang . flg - PETSC_TRUE if the shell matrix is linear; PETSC_FALSE otherwise. 53286919fd6SHong Zhang 53386919fd6SHong Zhang Level: intermediate 53486919fd6SHong Zhang @*/ 535d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsLinear(Mat A, PetscInt n, PetscBool *flg) 536d71ae5a4SJacob Faibussowitsch { 53786919fd6SHong Zhang Vec x, y, s1, s2; 53886919fd6SHong Zhang PetscRandom rctx; 53986919fd6SHong Zhang PetscScalar a; 54086919fd6SHong Zhang PetscInt k; 54186919fd6SHong Zhang PetscReal norm, normA; 54286919fd6SHong Zhang MPI_Comm comm; 54386919fd6SHong Zhang PetscMPIInt rank; 54486919fd6SHong Zhang 54586919fd6SHong Zhang PetscFunctionBegin; 54686919fd6SHong Zhang PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5479566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 5489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 54986919fd6SHong Zhang 5509566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rctx)); 5519566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 5529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &s1)); 5539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 5549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s1, &s2)); 55586919fd6SHong Zhang 55686919fd6SHong Zhang *flg = PETSC_TRUE; 55786919fd6SHong Zhang for (k = 0; k < n; k++) { 5589566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 5599566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 56048a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscRandomGetValue(rctx, &a)); 5619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm)); 56286919fd6SHong Zhang 56386919fd6SHong Zhang /* s2 = a*A*x + A*y */ 5649566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s2)); /* s2 = A*y */ 5659566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, s1)); /* s1 = A*x */ 5669566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, a, s1)); /* s2 = a s1 + s2 */ 56786919fd6SHong Zhang 56886919fd6SHong Zhang /* s1 = A * (a x + y) */ 5699566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, a, x)); /* y = a x + y */ 5709566063dSJacob Faibussowitsch PetscCall(MatMult(A, y, s1)); 5719566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_INFINITY, &normA)); 57286919fd6SHong Zhang 5739566063dSJacob Faibussowitsch PetscCall(VecAXPY(s2, -1.0, s1)); /* s2 = - s1 + s2 */ 5749566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_INFINITY, &norm)); 5751b8dac88SHong Zhang if (norm / normA > 100. * PETSC_MACHINE_EPSILON) { 57686919fd6SHong Zhang *flg = PETSC_FALSE; 5779566063dSJacob 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))); 57886919fd6SHong Zhang break; 57986919fd6SHong Zhang } 58086919fd6SHong Zhang } 5819566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 5829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 5839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 5849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 5859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 58686919fd6SHong Zhang PetscFunctionReturn(0); 58786919fd6SHong Zhang } 588