1 2 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3 4 static PetscErrorCode MatMultEqual_Private(Mat A, Mat B, PetscInt n, PetscBool *flg, PetscInt t, PetscInt add) 5 { 6 Vec Ax = NULL, Bx = NULL, s1 = NULL, s2 = NULL, Ay = NULL, By = NULL; 7 PetscRandom rctx; 8 PetscReal r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON; 9 PetscInt am, an, bm, bn, k; 10 PetscScalar none = -1.0; 11 #if defined(PETSC_USE_INFO) 12 const char *sops[] = {"MatMult", "MatMultAdd", "MatMultAdd (update)", "MatMultTranspose", "MatMultTransposeAdd", "MatMultTransposeAdd (update)", "MatMultHermitianTranspose", "MatMultHermitianTransposeAdd", "MatMultHermitianTransposeAdd (update)"}; 13 const char *sop; 14 #endif 15 16 PetscFunctionBegin; 17 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 18 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 19 PetscCheckSameComm(A, 1, B, 2); 20 PetscValidLogicalCollectiveInt(A, n, 3); 21 PetscValidBoolPointer(flg, 4); 22 PetscValidLogicalCollectiveInt(A, t, 5); 23 PetscValidLogicalCollectiveInt(A, add, 6); 24 PetscCall(MatGetLocalSize(A, &am, &an)); 25 PetscCall(MatGetLocalSize(B, &bm, &bn)); 26 PetscCheck(am == bm && an == bn, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat A,Mat B: local dim %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, am, bm, an, bn); 27 #if defined(PETSC_USE_INFO) 28 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 */ 29 #endif 30 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)A), &rctx)); 31 PetscCall(PetscRandomSetFromOptions(rctx)); 32 if (t) { 33 PetscCall(MatCreateVecs(A, &s1, &Ax)); 34 PetscCall(MatCreateVecs(B, &s2, &Bx)); 35 } else { 36 PetscCall(MatCreateVecs(A, &Ax, &s1)); 37 PetscCall(MatCreateVecs(B, &Bx, &s2)); 38 } 39 if (add) { 40 PetscCall(VecDuplicate(s1, &Ay)); 41 PetscCall(VecDuplicate(s2, &By)); 42 } 43 44 *flg = PETSC_TRUE; 45 for (k = 0; k < n; k++) { 46 Vec Aadd = NULL, Badd = NULL; 47 48 PetscCall(VecSetRandom(Ax, rctx)); 49 PetscCall(VecCopy(Ax, Bx)); 50 if (add) { 51 PetscCall(VecSetRandom(Ay, rctx)); 52 PetscCall(VecCopy(Ay, By)); 53 Aadd = Ay; 54 Badd = By; 55 if (add == 2) { 56 PetscCall(VecCopy(Ay, s1)); 57 PetscCall(VecCopy(By, s2)); 58 Aadd = s1; 59 Badd = s2; 60 } 61 } 62 if (t == 1) { 63 if (add) { 64 PetscCall(MatMultTransposeAdd(A, Ax, Aadd, s1)); 65 PetscCall(MatMultTransposeAdd(B, Bx, Badd, s2)); 66 } else { 67 PetscCall(MatMultTranspose(A, Ax, s1)); 68 PetscCall(MatMultTranspose(B, Bx, s2)); 69 } 70 } else if (t == 2) { 71 if (add) { 72 PetscCall(MatMultHermitianTransposeAdd(A, Ax, Aadd, s1)); 73 PetscCall(MatMultHermitianTransposeAdd(B, Bx, Badd, s2)); 74 } else { 75 PetscCall(MatMultHermitianTranspose(A, Ax, s1)); 76 PetscCall(MatMultHermitianTranspose(B, Bx, s2)); 77 } 78 } else { 79 if (add) { 80 PetscCall(MatMultAdd(A, Ax, Aadd, s1)); 81 PetscCall(MatMultAdd(B, Bx, Badd, s2)); 82 } else { 83 PetscCall(MatMult(A, Ax, s1)); 84 PetscCall(MatMult(B, Bx, s2)); 85 } 86 } 87 PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 88 if (r2 < tol) { 89 PetscCall(VecNorm(s1, NORM_INFINITY, &r1)); 90 } else { 91 PetscCall(VecAXPY(s2, none, s1)); 92 PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 93 r1 /= r2; 94 } 95 if (r1 > tol) { 96 *flg = PETSC_FALSE; 97 PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s() %g\n", k, sop, (double)r1)); 98 break; 99 } 100 } 101 PetscCall(PetscRandomDestroy(&rctx)); 102 PetscCall(VecDestroy(&Ax)); 103 PetscCall(VecDestroy(&Bx)); 104 PetscCall(VecDestroy(&Ay)); 105 PetscCall(VecDestroy(&By)); 106 PetscCall(VecDestroy(&s1)); 107 PetscCall(VecDestroy(&s2)); 108 PetscFunctionReturn(PETSC_SUCCESS); 109 } 110 111 static PetscErrorCode MatMatMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg, PetscBool At, PetscBool Bt) 112 { 113 Vec Ax, Bx, Cx, s1, s2, s3; 114 PetscRandom rctx; 115 PetscReal r1, r2, tol = PETSC_SQRT_MACHINE_EPSILON; 116 PetscInt am, an, bm, bn, cm, cn, k; 117 PetscScalar none = -1.0; 118 #if defined(PETSC_USE_INFO) 119 const char *sops[] = {"MatMatMult", "MatTransposeMatMult", "MatMatTransposeMult", "MatTransposeMatTransposeMult"}; 120 const char *sop; 121 #endif 122 123 PetscFunctionBegin; 124 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 125 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 126 PetscCheckSameComm(A, 1, B, 2); 127 PetscValidHeaderSpecific(C, MAT_CLASSID, 3); 128 PetscCheckSameComm(A, 1, C, 3); 129 PetscValidLogicalCollectiveInt(A, n, 4); 130 PetscValidBoolPointer(flg, 5); 131 PetscValidLogicalCollectiveBool(A, At, 6); 132 PetscValidLogicalCollectiveBool(B, Bt, 7); 133 PetscCall(MatGetLocalSize(A, &am, &an)); 134 PetscCall(MatGetLocalSize(B, &bm, &bn)); 135 PetscCall(MatGetLocalSize(C, &cm, &cn)); 136 if (At) { 137 PetscInt tt = an; 138 an = am; 139 am = tt; 140 } 141 if (Bt) { 142 PetscInt tt = bn; 143 bn = bm; 144 bm = tt; 145 } 146 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); 147 148 #if defined(PETSC_USE_INFO) 149 sop = sops[(At ? 1 : 0) + 2 * (Bt ? 1 : 0)]; 150 #endif 151 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)C), &rctx)); 152 PetscCall(PetscRandomSetFromOptions(rctx)); 153 if (Bt) { 154 PetscCall(MatCreateVecs(B, &s1, &Bx)); 155 } else { 156 PetscCall(MatCreateVecs(B, &Bx, &s1)); 157 } 158 if (At) { 159 PetscCall(MatCreateVecs(A, &s2, &Ax)); 160 } else { 161 PetscCall(MatCreateVecs(A, &Ax, &s2)); 162 } 163 PetscCall(MatCreateVecs(C, &Cx, &s3)); 164 165 *flg = PETSC_TRUE; 166 for (k = 0; k < n; k++) { 167 PetscCall(VecSetRandom(Bx, rctx)); 168 if (Bt) { 169 PetscCall(MatMultTranspose(B, Bx, s1)); 170 } else { 171 PetscCall(MatMult(B, Bx, s1)); 172 } 173 PetscCall(VecCopy(s1, Ax)); 174 if (At) { 175 PetscCall(MatMultTranspose(A, Ax, s2)); 176 } else { 177 PetscCall(MatMult(A, Ax, s2)); 178 } 179 PetscCall(VecCopy(Bx, Cx)); 180 PetscCall(MatMult(C, Cx, s3)); 181 182 PetscCall(VecNorm(s2, NORM_INFINITY, &r2)); 183 if (r2 < tol) { 184 PetscCall(VecNorm(s3, NORM_INFINITY, &r1)); 185 } else { 186 PetscCall(VecAXPY(s2, none, s3)); 187 PetscCall(VecNorm(s2, NORM_INFINITY, &r1)); 188 r1 /= r2; 189 } 190 if (r1 > tol) { 191 *flg = PETSC_FALSE; 192 PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th %s %g\n", k, sop, (double)r1)); 193 break; 194 } 195 } 196 PetscCall(PetscRandomDestroy(&rctx)); 197 PetscCall(VecDestroy(&Ax)); 198 PetscCall(VecDestroy(&Bx)); 199 PetscCall(VecDestroy(&Cx)); 200 PetscCall(VecDestroy(&s1)); 201 PetscCall(VecDestroy(&s2)); 202 PetscCall(VecDestroy(&s3)); 203 PetscFunctionReturn(PETSC_SUCCESS); 204 } 205 206 /*@ 207 MatMultEqual - Compares matrix-vector products of two matrices. 208 209 Collective 210 211 Input Parameters: 212 + A - the first matrix 213 . B - the second matrix 214 - n - number of random vectors to be tested 215 216 Output Parameter: 217 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 218 219 Level: intermediate 220 221 .seealso: `Mat`, `MatMultAddEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()`, `MatIsLinear()` 222 @*/ 223 PetscErrorCode MatMultEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 224 { 225 PetscFunctionBegin; 226 PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 0)); 227 PetscFunctionReturn(PETSC_SUCCESS); 228 } 229 230 /*@ 231 MatMultAddEqual - Compares matrix-vector product plus vector add of two matrices. 232 233 Collective 234 235 Input Parameters: 236 + A - the first matrix 237 . B - the second matrix 238 - n - number of random vectors to be tested 239 240 Output Parameter: 241 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 242 243 Level: intermediate 244 245 .seealso: `Mat`, `MatMultEqual()`, `MatMultTransposeEqual()`, `MatMultTransposeAddEqual()` 246 @*/ 247 PetscErrorCode MatMultAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 248 { 249 PetscFunctionBegin; 250 PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 1)); 251 PetscCall(MatMultEqual_Private(A, B, n, flg, 0, 2)); 252 PetscFunctionReturn(PETSC_SUCCESS); 253 } 254 255 /*@ 256 MatMultTransposeEqual - Compares matrix-vector products of two matrices. 257 258 Collective 259 260 Input Parameters: 261 + A - the first matrix 262 . B - the second matrix 263 - n - number of random vectors to be tested 264 265 Output Parameter: 266 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 267 268 Level: intermediate 269 270 .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeAddEqual()` 271 @*/ 272 PetscErrorCode MatMultTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 273 { 274 PetscFunctionBegin; 275 PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 0)); 276 PetscFunctionReturn(PETSC_SUCCESS); 277 } 278 279 /*@ 280 MatMultTransposeAddEqual - Compares matrix-vector products of two matrices. 281 282 Collective 283 284 Input Parameters: 285 + A - the first matrix 286 . B - the second matrix 287 - n - number of random vectors to be tested 288 289 Output Parameter: 290 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 291 292 Level: intermediate 293 294 .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 295 @*/ 296 PetscErrorCode MatMultTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 297 { 298 PetscFunctionBegin; 299 PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 1)); 300 PetscCall(MatMultEqual_Private(A, B, n, flg, 1, 2)); 301 PetscFunctionReturn(PETSC_SUCCESS); 302 } 303 304 /*@ 305 MatMultHermitianTransposeEqual - Compares matrix-vector products of two matrices. 306 307 Collective 308 309 Input Parameters: 310 + A - the first matrix 311 . B - the second matrix 312 - n - number of random vectors to be tested 313 314 Output Parameter: 315 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 316 317 Level: intermediate 318 319 .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 320 @*/ 321 PetscErrorCode MatMultHermitianTransposeEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 322 { 323 PetscFunctionBegin; 324 PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 0)); 325 PetscFunctionReturn(PETSC_SUCCESS); 326 } 327 328 /*@ 329 MatMultHermitianTransposeAddEqual - Compares matrix-vector products of two matrices. 330 331 Collective 332 333 Input Parameters: 334 + A - the first matrix 335 . B - the second matrix 336 - n - number of random vectors to be tested 337 338 Output Parameter: 339 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 340 341 Level: intermediate 342 343 .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 344 @*/ 345 PetscErrorCode MatMultHermitianTransposeAddEqual(Mat A, Mat B, PetscInt n, PetscBool *flg) 346 { 347 PetscFunctionBegin; 348 PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 1)); 349 PetscCall(MatMultEqual_Private(A, B, n, flg, 2, 2)); 350 PetscFunctionReturn(PETSC_SUCCESS); 351 } 352 353 /*@ 354 MatMatMultEqual - Test A*B*x = C*x for n random vector x 355 356 Collective 357 358 Input Parameters: 359 + A - the first matrix 360 . B - the second matrix 361 . C - the third matrix 362 - n - number of random vectors to be tested 363 364 Output Parameter: 365 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 366 367 Level: intermediate 368 369 .seealso: `Mat`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 370 @*/ 371 PetscErrorCode MatMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 372 { 373 PetscFunctionBegin; 374 PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_FALSE)); 375 PetscFunctionReturn(PETSC_SUCCESS); 376 } 377 378 /*@ 379 MatTransposeMatMultEqual - Test A^T*B*x = C*x for n random vector x 380 381 Collective 382 383 Input Parameters: 384 + A - the first matrix 385 . B - the second matrix 386 . C - the third matrix 387 - n - number of random vectors to be tested 388 389 Output Parameter: 390 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 391 392 Level: intermediate 393 394 .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 395 @*/ 396 PetscErrorCode MatTransposeMatMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 397 { 398 PetscFunctionBegin; 399 PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_TRUE, PETSC_FALSE)); 400 PetscFunctionReturn(PETSC_SUCCESS); 401 } 402 403 /*@ 404 MatMatTransposeMultEqual - Test A*B^T*x = C*x for n random vector x 405 406 Collective 407 408 Input Parameters: 409 + A - the first matrix 410 . B - the second matrix 411 . C - the third matrix 412 - n - number of random vectors to be tested 413 414 Output Parameter: 415 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 416 417 Level: intermediate 418 419 .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 420 @*/ 421 PetscErrorCode MatMatTransposeMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 422 { 423 PetscFunctionBegin; 424 PetscCall(MatMatMultEqual_Private(A, B, C, n, flg, PETSC_FALSE, PETSC_TRUE)); 425 PetscFunctionReturn(PETSC_SUCCESS); 426 } 427 428 static PetscErrorCode MatProjMultEqual_Private(Mat A, Mat B, Mat C, PetscInt n, PetscBool rart, PetscBool *flg) 429 { 430 Vec x, v1, v2, v3, v4, Cx, Bx; 431 PetscReal norm_abs, norm_rel, tol = PETSC_SQRT_MACHINE_EPSILON; 432 PetscInt i, am, an, bm, bn, cm, cn; 433 PetscRandom rdm; 434 PetscScalar none = -1.0; 435 436 PetscFunctionBegin; 437 PetscCall(MatGetLocalSize(A, &am, &an)); 438 PetscCall(MatGetLocalSize(B, &bm, &bn)); 439 if (rart) { 440 PetscInt t = bm; 441 bm = bn; 442 bn = t; 443 } 444 PetscCall(MatGetLocalSize(C, &cm, &cn)); 445 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); 446 447 /* Create left vector of A: v2 */ 448 PetscCall(MatCreateVecs(A, &Bx, &v2)); 449 450 /* Create right vectors of B: x, v3, v4 */ 451 if (rart) { 452 PetscCall(MatCreateVecs(B, &v1, &x)); 453 } else { 454 PetscCall(MatCreateVecs(B, &x, &v1)); 455 } 456 PetscCall(VecDuplicate(x, &v3)); 457 458 PetscCall(MatCreateVecs(C, &Cx, &v4)); 459 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 460 PetscCall(PetscRandomSetFromOptions(rdm)); 461 462 *flg = PETSC_TRUE; 463 for (i = 0; i < n; i++) { 464 PetscCall(VecSetRandom(x, rdm)); 465 PetscCall(VecCopy(x, Cx)); 466 PetscCall(MatMult(C, Cx, v4)); /* v4 = C*x */ 467 if (rart) { 468 PetscCall(MatMultTranspose(B, x, v1)); 469 } else { 470 PetscCall(MatMult(B, x, v1)); 471 } 472 PetscCall(VecCopy(v1, Bx)); 473 PetscCall(MatMult(A, Bx, v2)); /* v2 = A*B*x */ 474 PetscCall(VecCopy(v2, v1)); 475 if (rart) { 476 PetscCall(MatMult(B, v1, v3)); /* v3 = R*A*R^t*x */ 477 } else { 478 PetscCall(MatMultTranspose(B, v1, v3)); /* v3 = Bt*A*B*x */ 479 } 480 PetscCall(VecNorm(v4, NORM_2, &norm_abs)); 481 PetscCall(VecAXPY(v4, none, v3)); 482 PetscCall(VecNorm(v4, NORM_2, &norm_rel)); 483 484 if (norm_abs > tol) norm_rel /= norm_abs; 485 if (norm_rel > tol) { 486 *flg = PETSC_FALSE; 487 PetscCall(PetscInfo(A, "Error: %" PetscInt_FMT "-th Mat%sMult() %g\n", i, rart ? "RARt" : "PtAP", (double)norm_rel)); 488 break; 489 } 490 } 491 492 PetscCall(PetscRandomDestroy(&rdm)); 493 PetscCall(VecDestroy(&x)); 494 PetscCall(VecDestroy(&Bx)); 495 PetscCall(VecDestroy(&Cx)); 496 PetscCall(VecDestroy(&v1)); 497 PetscCall(VecDestroy(&v2)); 498 PetscCall(VecDestroy(&v3)); 499 PetscCall(VecDestroy(&v4)); 500 PetscFunctionReturn(PETSC_SUCCESS); 501 } 502 503 /*@ 504 MatPtAPMultEqual - Compares matrix-vector products of C = Bt*A*B 505 506 Collective 507 508 Input Parameters: 509 + A - the first matrix 510 . B - the second matrix 511 . C - the third matrix 512 - n - number of random vectors to be tested 513 514 Output Parameter: 515 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 516 517 Level: intermediate 518 519 .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 520 @*/ 521 PetscErrorCode MatPtAPMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 522 { 523 PetscFunctionBegin; 524 PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_FALSE, flg)); 525 PetscFunctionReturn(PETSC_SUCCESS); 526 } 527 528 /*@ 529 MatRARtMultEqual - Compares matrix-vector products of C = B*A*B^t 530 531 Collective 532 533 Input Parameters: 534 + A - the first matrix 535 . B - the second matrix 536 . C - the third matrix 537 - n - number of random vectors to be tested 538 539 Output Parameter: 540 . flg - `PETSC_TRUE` if the products are equal; `PETSC_FALSE` otherwise. 541 542 Level: intermediate 543 544 .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 545 @*/ 546 PetscErrorCode MatRARtMultEqual(Mat A, Mat B, Mat C, PetscInt n, PetscBool *flg) 547 { 548 PetscFunctionBegin; 549 PetscCall(MatProjMultEqual_Private(A, B, C, n, PETSC_TRUE, flg)); 550 PetscFunctionReturn(PETSC_SUCCESS); 551 } 552 553 /*@ 554 MatIsLinear - Check if a shell matrix `A` is a linear operator. 555 556 Collective 557 558 Input Parameters: 559 + A - the shell matrix 560 - n - number of random vectors to be tested 561 562 Output Parameter: 563 . flg - `PETSC_TRUE` if the shell matrix is linear; `PETSC_FALSE` otherwise. 564 565 Level: intermediate 566 567 .seealso: `Mat`, `MatMatMultEqual()`, `MatMultEqual()`, `MatMultAddEqual()`, `MatMultTransposeEqual()` 568 @*/ 569 PetscErrorCode MatIsLinear(Mat A, PetscInt n, PetscBool *flg) 570 { 571 Vec x, y, s1, s2; 572 PetscRandom rctx; 573 PetscScalar a; 574 PetscInt k; 575 PetscReal norm, normA; 576 MPI_Comm comm; 577 PetscMPIInt rank; 578 579 PetscFunctionBegin; 580 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 581 PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 582 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 583 584 PetscCall(PetscRandomCreate(comm, &rctx)); 585 PetscCall(PetscRandomSetFromOptions(rctx)); 586 PetscCall(MatCreateVecs(A, &x, &s1)); 587 PetscCall(VecDuplicate(x, &y)); 588 PetscCall(VecDuplicate(s1, &s2)); 589 590 *flg = PETSC_TRUE; 591 for (k = 0; k < n; k++) { 592 PetscCall(VecSetRandom(x, rctx)); 593 PetscCall(VecSetRandom(y, rctx)); 594 if (rank == 0) PetscCall(PetscRandomGetValue(rctx, &a)); 595 PetscCallMPI(MPI_Bcast(&a, 1, MPIU_SCALAR, 0, comm)); 596 597 /* s2 = a*A*x + A*y */ 598 PetscCall(MatMult(A, y, s2)); /* s2 = A*y */ 599 PetscCall(MatMult(A, x, s1)); /* s1 = A*x */ 600 PetscCall(VecAXPY(s2, a, s1)); /* s2 = a s1 + s2 */ 601 602 /* s1 = A * (a x + y) */ 603 PetscCall(VecAXPY(y, a, x)); /* y = a x + y */ 604 PetscCall(MatMult(A, y, s1)); 605 PetscCall(VecNorm(s1, NORM_INFINITY, &normA)); 606 607 PetscCall(VecAXPY(s2, -1.0, s1)); /* s2 = - s1 + s2 */ 608 PetscCall(VecNorm(s2, NORM_INFINITY, &norm)); 609 if (norm / normA > 100. * PETSC_MACHINE_EPSILON) { 610 *flg = PETSC_FALSE; 611 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))); 612 break; 613 } 614 } 615 PetscCall(PetscRandomDestroy(&rctx)); 616 PetscCall(VecDestroy(&x)); 617 PetscCall(VecDestroy(&y)); 618 PetscCall(VecDestroy(&s1)); 619 PetscCall(VecDestroy(&s2)); 620 PetscFunctionReturn(PETSC_SUCCESS); 621 } 622