1c4762a1bSJed Brown static char help[] = "Tests the use of MatTranspose_Nest and MatMatMult_Nest_Dense\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch PetscErrorCode TestInitialMatrix(void) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown const PetscInt nr = 2, nc = 3, nk = 10; 8c4762a1bSJed Brown PetscInt n, N, m, M; 9c4762a1bSJed Brown const PetscInt arow[2 * 3] = {2, 2, 2, 3, 3, 3}; 10c4762a1bSJed Brown const PetscInt acol[2 * 3] = {3, 2, 4, 3, 2, 4}; 11c4762a1bSJed Brown Mat A, Atranspose, B, C; 12c4762a1bSJed Brown Mat subs[2 * 3], **block; 13c4762a1bSJed Brown Vec x, y, Ax, ATy; 14c4762a1bSJed Brown PetscInt i, j; 15c20d7725SJed Brown PetscScalar dot1, dot2, zero = 0.0, one = 1.0, *valsB, *valsC; 16c4762a1bSJed Brown PetscReal norm; 17c4762a1bSJed Brown PetscRandom rctx; 18c20d7725SJed Brown PetscBool equal; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 22c4762a1bSJed Brown /* Force the random numbers to have imaginary part 0 so printed results are the same for --with-scalar-type=real or --with-scalar-type=complex */ 239566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rctx, zero, one)); 249566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 2548a46eb9SPierre Jolivet for (i = 0; i < (nr * nc); i++) PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, arow[i], acol[i], NULL, &subs[i])); 269566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, nr, NULL, nc, NULL, subs, &A)); 279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 289566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &y)); 299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &ATy)); 309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Ax)); 319566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A, rctx)); 329566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atranspose)); 33c4762a1bSJed Brown 349566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 359566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(A, NULL, NULL, &block)); 36c4762a1bSJed Brown for (i = 0; i < nr; i++) { 3748a46eb9SPierre Jolivet for (j = 0; j < nc; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(MatView(Atranspose, PETSC_VIEWER_STDOUT_WORLD)); 419566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(Atranspose, NULL, NULL, &block)); 42c4762a1bSJed Brown for (i = 0; i < nc; i++) { 4348a46eb9SPierre Jolivet for (j = 0; j < nr; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Check <Ax, y> = <x, A^Ty> */ 47c4762a1bSJed Brown for (i = 0; i < 10; i++) { 489566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 499566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, Ax)); 529566063dSJacob Faibussowitsch PetscCall(VecDot(Ax, y, &dot1)); 539566063dSJacob Faibussowitsch PetscCall(MatMult(Atranspose, y, ATy)); 549566063dSJacob Faibussowitsch PetscCall(VecDot(ATy, x, &dot2)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "<Ax, y> = %g\n", (double)PetscRealPart(dot1))); 579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "<x, A^Ty> = %g\n", (double)PetscRealPart(dot2))); 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, acol[0] + acol[nr] + acol[2 * nr], nk, NULL, &B)); 649566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B, rctx)); 65fb842aefSJose E. Roman PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C)); 66fb842aefSJose E. Roman PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DETERMINE, &C)); 679566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &equal)); 6828b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C != A*B"); 69c20d7725SJed Brown 709566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 719566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 72c4762a1bSJed Brown for (i = 0; i < nk; i++) { 739566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(B, i, &valsB)); 749566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, valsB, &x)); 759566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &Ax)); 769566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, Ax)); 779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 789566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(C, i, &valsC)); 799566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, m, M, valsC, &y)); 809566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, Ax)); 819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 839566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(C, &valsC)); 849566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(B, &valsB)); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(MatNorm(C, NORM_INFINITY, &norm)); 8708401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatMatMult(): %g", (double)norm); 889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 90c4762a1bSJed Brown 9148a46eb9SPierre Jolivet for (i = 0; i < (nr * nc); i++) PetscCall(MatDestroy(&subs[i])); 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atranspose)); 949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ATy)); 959566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99d71ae5a4SJacob Faibussowitsch PetscErrorCode TestReuseMatrix(void) 100d71ae5a4SJacob Faibussowitsch { 101c4762a1bSJed Brown const PetscInt n = 2; 102c4762a1bSJed Brown Mat A; 103c4762a1bSJed Brown Mat subs[2 * 2], **block; 104c4762a1bSJed Brown PetscInt i, j; 105c4762a1bSJed Brown PetscRandom rctx; 106c4762a1bSJed Brown PetscScalar zero = 0.0, one = 1.0; 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscFunctionBegin; 1099566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 1109566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rctx, zero, one)); 1119566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 11248a46eb9SPierre Jolivet for (i = 0; i < (n * n); i++) PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, n, n, NULL, &subs[i])); 1139566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, n, NULL, n, NULL, subs, &A)); 1149566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A, rctx)); 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1179566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(A, NULL, NULL, &block)); 118c4762a1bSJed Brown for (i = 0; i < n; i++) { 11948a46eb9SPierre Jolivet for (j = 0; j < n; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 120c4762a1bSJed Brown } 1219566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INPLACE_MATRIX, &A)); 1229566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1239566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(A, NULL, NULL, &block)); 124c4762a1bSJed Brown for (i = 0; i < n; i++) { 12548a46eb9SPierre Jolivet for (j = 0; j < n; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 12848a46eb9SPierre Jolivet for (i = 0; i < (n * n); i++) PetscCall(MatDestroy(&subs[i])); 1299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1309566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 135d71ae5a4SJacob Faibussowitsch { 136327415f7SBarry Smith PetscFunctionBeginUser; 137*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 1389566063dSJacob Faibussowitsch PetscCall(TestInitialMatrix()); 1399566063dSJacob Faibussowitsch PetscCall(TestReuseMatrix()); 1409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 141b122ec5aSJacob Faibussowitsch return 0; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown /*TEST 145c4762a1bSJed Brown 146c4762a1bSJed Brown test: 147c4762a1bSJed Brown args: -malloc_dump 148c4762a1bSJed Brown 149c4762a1bSJed Brown TEST*/ 150