1c4762a1bSJed Brown static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Example: 4c4762a1bSJed Brown mpiexec -n <np> ./ex104 -mat_type elemental 5c4762a1bSJed Brown */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petscmat.h> 8c4762a1bSJed Brown 99371c9d4SSatish Balay int main(int argc, char **argv) { 10c4762a1bSJed Brown Mat A, B, C, D; 11c4762a1bSJed Brown PetscInt i, M = 10, N = 5, j, nrows, ncols, am, an, rstart, rend; 12c4762a1bSJed Brown PetscRandom r; 13c20d7725SJed Brown PetscBool equal, Aiselemental; 14c4762a1bSJed Brown PetscReal fill = 1.0; 15c4762a1bSJed Brown IS isrows, iscols; 16c4762a1bSJed Brown const PetscInt *rows, *cols; 17c4762a1bSJed Brown PetscScalar *v, rval; 18c4762a1bSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 19c4762a1bSJed Brown PetscBool Test_MatMatMult = PETSC_TRUE; 20c4762a1bSJed Brown #else 21c4762a1bSJed Brown PetscBool Test_MatMatMult = PETSC_FALSE; 22c4762a1bSJed Brown #endif 23c4762a1bSJed Brown PetscMPIInt size; 24c4762a1bSJed Brown 25327415f7SBarry Smith PetscFunctionBeginUser; 269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 28c4762a1bSJed Brown 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL)); 319566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 339566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 359566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 369566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r)); 379566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Set local matrix entries */ 409566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A, &isrows, &iscols)); 419566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows, &nrows)); 429566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows, &rows)); 439566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols, &ncols)); 449566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols, &cols)); 459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows * ncols, &v)); 46c4762a1bSJed Brown for (i = 0; i < nrows; i++) { 47c4762a1bSJed Brown for (j = 0; j < ncols; j++) { 489566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rval)); 49c4762a1bSJed Brown v[i * ncols + j] = rval; 50c4762a1bSJed Brown } 51c4762a1bSJed Brown } 529566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, nrows, rows, ncols, cols, v, INSERT_VALUES)); 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 559566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows, &rows)); 569566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols, &cols)); 579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 599566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 60c4762a1bSJed Brown 619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATELEMENTAL, &Aiselemental)); 62c20d7725SJed Brown 63c20d7725SJed Brown /* Test MatCreateTranspose() and MatTranspose() */ 649566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A, &C)); 659566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B)); /* B = A^T */ 669566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C, B, 10, &equal)); 6728b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A^T*x != (x^T*A)^T"); 689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 69c20d7725SJed Brown 709566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B)); 71c20d7725SJed Brown if (!Aiselemental) { 729566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_INPLACE_MATRIX, &B)); 739566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C, B, 10, &equal)); 7428b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "C*x != B*x"); 75c20d7725SJed Brown } 769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 77c20d7725SJed Brown 78c20d7725SJed Brown /* Test B = C*A for matrix type transpose and seqdense */ 79c20d7725SJed Brown if (size == 1 && !Aiselemental) { 809566063dSJacob Faibussowitsch PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &B)); 819566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(C, A, B, 10, &equal)); 8228b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B != C*A for matrix type transpose and seqdense"); 839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 84c20d7725SJed Brown } 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Test MatMatMult() */ 88c4762a1bSJed Brown if (Test_MatMatMult) { 89c4762a1bSJed Brown #if !defined(PETSC_HAVE_ELEMENTAL) 90be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "This test requires ELEMENTAL"); 91c4762a1bSJed Brown #endif 929566063dSJacob Faibussowitsch PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B)); /* B = A^T */ 939566063dSJacob Faibussowitsch PetscCall(MatMatMult(B, A, MAT_INITIAL_MATRIX, fill, &C)); /* C = B*A = A^T*A */ 949566063dSJacob Faibussowitsch PetscCall(MatMatMult(B, A, MAT_REUSE_MATRIX, fill, &C)); 959566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(B, A, C, 10, &equal)); 9628b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B*A*x != C*x"); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 999566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &D)); 100c20d7725SJed Brown 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* Test MatTransposeMatMult() */ 107c20d7725SJed Brown if (!Aiselemental) { 1089566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, A, MAT_INITIAL_MATRIX, fill, &D)); /* D = A^T*A */ 1099566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, A, MAT_REUSE_MATRIX, fill, &D)); 1109566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMultEqual(A, A, D, 10, &equal)); 11128b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x"); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 1149566063dSJacob Faibussowitsch PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C)); 1159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* Test D*x = A^T*C*A*x, where C is in AIJ format */ 1199566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 1209566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 121c4762a1bSJed Brown if (size == 1) { 1229566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, am, am)); 123c4762a1bSJed Brown } else { 1249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, am, am, PETSC_DECIDE, PETSC_DECIDE)); 125c4762a1bSJed Brown } 1269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 1279566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 1289566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 129c4762a1bSJed Brown v[0] = 1.0; 130*48a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(C, 1, &i, 1, &i, v, INSERT_VALUES)); 1319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 1329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* B = C*A, D = A^T*B */ 1359566063dSJacob Faibussowitsch PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, 1.0, &B)); 1369566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, fill, &D)); 1379566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMultEqual(A, B, D, 10, &equal)); 13828b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*B*x"); 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Test MatMatTransposeMult() */ 146c20d7725SJed Brown if (!Aiselemental) { 147c4762a1bSJed Brown PetscReal diff, scale; 148c4762a1bSJed Brown PetscInt am, an, aM, aN; 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 1519566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &aM, &aN)); 1529566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), PETSC_DECIDE, an, aM + 10, aN, NULL, &B)); 1539566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B, NULL)); 1549566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(A, B, MAT_INITIAL_MATRIX, fill, &D)); /* D = A*A^T */ 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* Test MatDuplicate for matrix product */ 1579566063dSJacob Faibussowitsch PetscCall(MatDuplicate(D, MAT_COPY_VALUES, &C)); 158c4762a1bSJed Brown 1599566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(A, B, MAT_REUSE_MATRIX, fill, &D)); 1609566063dSJacob Faibussowitsch PetscCall(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN)); 161c4762a1bSJed Brown 1629566063dSJacob Faibussowitsch PetscCall(MatNorm(C, NORM_FROBENIUS, &diff)); 1639566063dSJacob Faibussowitsch PetscCall(MatNorm(D, NORM_FROBENIUS, &scale)); 164e00437b9SBarry Smith PetscCheck(diff <= PETSC_SMALL * scale, PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX"); 1659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultEqual(A, B, D, 10, &equal)); 16828b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "D*x != A^T*A*x"); 1699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 1709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 1739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 1759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 176b122ec5aSJacob Faibussowitsch return 0; 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /*TEST 180c4762a1bSJed Brown 181c4762a1bSJed Brown test: 182c4762a1bSJed Brown output_file: output/ex104.out 183c4762a1bSJed Brown 184c4762a1bSJed Brown test: 185c4762a1bSJed Brown suffix: 2 186c4762a1bSJed Brown nsize: 2 187c4762a1bSJed Brown output_file: output/ex104.out 188c4762a1bSJed Brown 189c4762a1bSJed Brown test: 190c4762a1bSJed Brown suffix: 3 191c4762a1bSJed Brown nsize: 4 192c4762a1bSJed Brown output_file: output/ex104.out 193c4762a1bSJed Brown args: -M 23 -N 31 194c4762a1bSJed Brown 195c4762a1bSJed Brown test: 196c4762a1bSJed Brown suffix: 4 197c4762a1bSJed Brown nsize: 4 198c4762a1bSJed Brown output_file: output/ex104.out 199c4762a1bSJed Brown args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic 200c4762a1bSJed Brown 201c4762a1bSJed Brown test: 202c4762a1bSJed Brown suffix: 5 203c4762a1bSJed Brown nsize: 4 204c4762a1bSJed Brown output_file: output/ex104.out 205c4762a1bSJed Brown args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv 206c4762a1bSJed Brown 207c20d7725SJed Brown test: 208c20d7725SJed Brown suffix: 6 209c20d7725SJed Brown args: -mat_type elemental 210c20d7725SJed Brown requires: elemental 211c20d7725SJed Brown output_file: output/ex104.out 212c20d7725SJed Brown 213c20d7725SJed Brown test: 214c20d7725SJed Brown suffix: 7 215c20d7725SJed Brown nsize: 2 216c20d7725SJed Brown args: -mat_type elemental 217c20d7725SJed Brown requires: elemental 218c20d7725SJed Brown output_file: output/ex104.out 219c20d7725SJed Brown 220c4762a1bSJed Brown TEST*/ 221