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 59371c9d4SSatish Balay PetscErrorCode TestInitialMatrix(void) { 6c4762a1bSJed Brown const PetscInt nr = 2, nc = 3, nk = 10; 7c4762a1bSJed Brown PetscInt n, N, m, M; 8c4762a1bSJed Brown const PetscInt arow[2 * 3] = {2, 2, 2, 3, 3, 3}; 9c4762a1bSJed Brown const PetscInt acol[2 * 3] = {3, 2, 4, 3, 2, 4}; 10c4762a1bSJed Brown Mat A, Atranspose, B, C; 11c4762a1bSJed Brown Mat subs[2 * 3], **block; 12c4762a1bSJed Brown Vec x, y, Ax, ATy; 13c4762a1bSJed Brown PetscInt i, j; 14c20d7725SJed Brown PetscScalar dot1, dot2, zero = 0.0, one = 1.0, *valsB, *valsC; 15c4762a1bSJed Brown PetscReal norm; 16c4762a1bSJed Brown PetscRandom rctx; 17c20d7725SJed Brown PetscBool equal; 18c4762a1bSJed Brown 19c4762a1bSJed Brown PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 21c4762a1bSJed 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 */ 229566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rctx, zero, one)); 239566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 24*48a46eb9SPierre Jolivet for (i = 0; i < (nr * nc); i++) PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, arow[i], acol[i], NULL, &subs[i])); 259566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, nr, NULL, nc, NULL, subs, &A)); 269566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &y)); 289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &ATy)); 299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(y, &Ax)); 309566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A, rctx)); 319566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atranspose)); 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 349566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(A, NULL, NULL, &block)); 35c4762a1bSJed Brown for (i = 0; i < nr; i++) { 36*48a46eb9SPierre Jolivet for (j = 0; j < nc; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 399566063dSJacob Faibussowitsch PetscCall(MatView(Atranspose, PETSC_VIEWER_STDOUT_WORLD)); 409566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(Atranspose, NULL, NULL, &block)); 41c4762a1bSJed Brown for (i = 0; i < nc; i++) { 42*48a46eb9SPierre Jolivet for (j = 0; j < nr; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Check <Ax, y> = <x, A^Ty> */ 46c4762a1bSJed Brown for (i = 0; i < 10; i++) { 479566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rctx)); 489566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y, rctx)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, Ax)); 519566063dSJacob Faibussowitsch PetscCall(VecDot(Ax, y, &dot1)); 529566063dSJacob Faibussowitsch PetscCall(MatMult(Atranspose, y, ATy)); 539566063dSJacob Faibussowitsch PetscCall(VecDot(ATy, x, &dot2)); 54c4762a1bSJed Brown 559566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "<Ax, y> = %g\n", (double)PetscRealPart(dot1))); 569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "<x, A^Ty> = %g\n", (double)PetscRealPart(dot2))); 57c4762a1bSJed Brown } 589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, acol[0] + acol[nr] + acol[2 * nr], nk, NULL, &B)); 639566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B, rctx)); 649566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 659566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 669566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &equal)); 6728b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in C != A*B"); 68c20d7725SJed Brown 699566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 709566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 71c4762a1bSJed Brown for (i = 0; i < nk; i++) { 729566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(B, i, &valsB)); 739566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, N, valsB, &x)); 749566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &Ax)); 759566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, Ax)); 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 779566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(C, i, &valsC)); 789566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, m, M, valsC, &y)); 799566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, Ax)); 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ax)); 819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 829566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(C, &valsC)); 839566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(B, &valsB)); 84c4762a1bSJed Brown } 859566063dSJacob Faibussowitsch PetscCall(MatNorm(C, NORM_INFINITY, &norm)); 8608401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatMatMult(): %g", (double)norm); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 89c4762a1bSJed Brown 90*48a46eb9SPierre Jolivet for (i = 0; i < (nr * nc); i++) PetscCall(MatDestroy(&subs[i])); 919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atranspose)); 939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ATy)); 949566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 95c4762a1bSJed Brown PetscFunctionReturn(0); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 989371c9d4SSatish Balay PetscErrorCode TestReuseMatrix(void) { 99c4762a1bSJed Brown const PetscInt n = 2; 100c4762a1bSJed Brown Mat A; 101c4762a1bSJed Brown Mat subs[2 * 2], **block; 102c4762a1bSJed Brown PetscInt i, j; 103c4762a1bSJed Brown PetscRandom rctx; 104c4762a1bSJed Brown PetscScalar zero = 0.0, one = 1.0; 105c4762a1bSJed Brown 106c4762a1bSJed Brown PetscFunctionBegin; 1079566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 1089566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rctx, zero, one)); 1099566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rctx)); 110*48a46eb9SPierre Jolivet for (i = 0; i < (n * n); i++) PetscCall(MatCreateSeqDense(PETSC_COMM_WORLD, n, n, NULL, &subs[i])); 1119566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, n, NULL, n, NULL, subs, &A)); 1129566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A, rctx)); 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1159566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(A, NULL, NULL, &block)); 116c4762a1bSJed Brown for (i = 0; i < n; i++) { 117*48a46eb9SPierre Jolivet for (j = 0; j < n; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 118c4762a1bSJed Brown } 1199566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INPLACE_MATRIX, &A)); 1209566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1219566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(A, NULL, NULL, &block)); 122c4762a1bSJed Brown for (i = 0; i < n; i++) { 123*48a46eb9SPierre Jolivet for (j = 0; j < n; j++) PetscCall(MatView(block[i][j], PETSC_VIEWER_STDOUT_WORLD)); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126*48a46eb9SPierre Jolivet for (i = 0; i < (n * n); i++) PetscCall(MatDestroy(&subs[i])); 1279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1289566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 129c4762a1bSJed Brown PetscFunctionReturn(0); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 1329371c9d4SSatish Balay int main(int argc, char **args) { 133327415f7SBarry Smith PetscFunctionBeginUser; 1349566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 1359566063dSJacob Faibussowitsch PetscCall(TestInitialMatrix()); 1369566063dSJacob Faibussowitsch PetscCall(TestReuseMatrix()); 1379566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 138b122ec5aSJacob Faibussowitsch return 0; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /*TEST 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown args: -malloc_dump 145c4762a1bSJed Brown 146c4762a1bSJed Brown TEST*/ 147