1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 59371c9d4SSatish Balay int main(int argc, char **argv) { 6c20d7725SJed Brown Mat A, B, C, D, AT; 7c4762a1bSJed Brown PetscInt i, M, N, Istart, Iend, n = 7, j, J, Ii, m = 8, am, an; 8c4762a1bSJed Brown PetscScalar v; 9c4762a1bSJed Brown PetscRandom r; 103fbca975SStefano Zampini PetscBool equal = PETSC_FALSE, flg; 11c20d7725SJed Brown PetscReal fill = 1.0, norm; 12c4762a1bSJed Brown PetscMPIInt size; 13c4762a1bSJed Brown 14327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL)); 19c4762a1bSJed Brown 209566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r)); 219566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Create a aij matrix A */ 24c4762a1bSJed Brown M = N = m * n; 259566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N)); 279566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 289566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 299566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, 5, NULL, 5, NULL)); 309566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL)); 31c4762a1bSJed Brown 329566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &Istart, &Iend)); 33c4762a1bSJed Brown am = Iend - Istart; 34c4762a1bSJed Brown for (Ii = Istart; Ii < Iend; Ii++) { 359371c9d4SSatish Balay v = -1.0; 369371c9d4SSatish Balay i = Ii / n; 379371c9d4SSatish Balay j = Ii - i * n; 389371c9d4SSatish Balay if (i > 0) { 399371c9d4SSatish Balay J = Ii - n; 409371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 419371c9d4SSatish Balay } 429371c9d4SSatish Balay if (i < m - 1) { 439371c9d4SSatish Balay J = Ii + n; 449371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 459371c9d4SSatish Balay } 469371c9d4SSatish Balay if (j > 0) { 479371c9d4SSatish Balay J = Ii - 1; 489371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 499371c9d4SSatish Balay } 509371c9d4SSatish Balay if (j < n - 1) { 519371c9d4SSatish Balay J = Ii + 1; 529371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 539371c9d4SSatish Balay } 549371c9d4SSatish Balay v = 4.0; 559371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 56c4762a1bSJed Brown } 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create a dense matrix B */ 619566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &am, &an)); 629566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 639566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, an, PETSC_DECIDE, PETSC_DECIDE, M)); 649566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATDENSE)); 659566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(B, NULL)); 669566063dSJacob Faibussowitsch PetscCall(MatMPIDenseSetPreallocation(B, NULL)); 679566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 689566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B, r)); 699566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 70c4762a1bSJed Brown 71c20d7725SJed Brown /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */ 729566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 739566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 749566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, am, PETSC_DECIDE, PETSC_DECIDE, N)); 759566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 769566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 779566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(C)); 789566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, fill, &C)); 799566063dSJacob Faibussowitsch PetscCall(MatNorm(C, NORM_INFINITY, &norm)); 809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 81c20d7725SJed Brown 82c4762a1bSJed Brown /* Test C = A*B (aij*dense) */ 839566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, fill, &C)); 849566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, fill, &C)); 85c4762a1bSJed Brown 86c20d7725SJed Brown /* Test developer API */ 879566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A, B, NULL, &D)); 889566063dSJacob Faibussowitsch PetscCall(MatProductSetType(D, MATPRODUCT_AB)); 899566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(D, "default")); 909566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(D, fill)); 919566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(D)); 929566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(D)); 93*48a46eb9SPierre Jolivet for (i = 0; i < 2; i++) PetscCall(MatProductNumeric(D)); 949566063dSJacob Faibussowitsch PetscCall(MatEqual(C, D, &equal)); 9528b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "C != D"); 969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 97c4762a1bSJed Brown 98c20d7725SJed Brown /* Test D = AT*B (transpose(aij)*dense) */ 999566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A, &AT)); 1009566063dSJacob Faibussowitsch PetscCall(MatMatMult(AT, B, MAT_INITIAL_MATRIX, fill, &D)); 1019566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(AT, B, D, 10, &equal)); 10228b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != AT*B (transpose(aij)*dense)"); 1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AT)); 105c20d7725SJed Brown 106c4762a1bSJed Brown /* Test D = C*A (dense*aij) */ 1079566063dSJacob Faibussowitsch PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &D)); 1089566063dSJacob Faibussowitsch PetscCall(MatMatMult(C, A, MAT_REUSE_MATRIX, fill, &D)); 1099566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(C, A, D, 10, &equal)); 11028b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != C*A (dense*aij)"); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Test D = A*C (aij*dense) */ 1149566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, C, MAT_INITIAL_MATRIX, fill, &D)); 1159566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, C, MAT_REUSE_MATRIX, fill, &D)); 1169566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, C, D, 10, &equal)); 11728b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != A*C (aij*dense)"); 1189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Test D = B*C (dense*dense) */ 1219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 122c4762a1bSJed Brown if (size == 1) { 1239566063dSJacob Faibussowitsch PetscCall(MatMatMult(B, C, MAT_INITIAL_MATRIX, fill, &D)); 1249566063dSJacob Faibussowitsch PetscCall(MatMatMult(B, C, MAT_REUSE_MATRIX, fill, &D)); 1259566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(B, C, D, 10, &equal)); 12628b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != B*C (dense*dense)"); 1279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c20d7725SJed Brown /* Test D = B*C^T (dense*dense) */ 1319566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(B, C, MAT_INITIAL_MATRIX, fill, &D)); 1329566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(B, C, MAT_REUSE_MATRIX, fill, &D)); 1339566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultEqual(B, C, D, 10, &equal)); 13428b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "D != B*C^T (dense*dense)"); 1359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 136c20d7725SJed Brown 1373fbca975SStefano Zampini /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */ 1383fbca975SStefano Zampini flg = PETSC_FALSE; 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_userAPI", &flg)); 1403fbca975SStefano Zampini if (flg) { 1413fbca975SStefano Zampini /* user driver */ 1429566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, C, MAT_REUSE_MATRIX, fill, &B)); 1433fbca975SStefano Zampini } else { 1443fbca975SStefano Zampini /* clear internal data structures related with previous products to avoid circular references */ 1459566063dSJacob Faibussowitsch PetscCall(MatProductClear(A)); 1469566063dSJacob Faibussowitsch PetscCall(MatProductClear(B)); 1479566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 1489566063dSJacob Faibussowitsch PetscCall(MatProductCreateWithMat(A, C, NULL, B)); 1499566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B, MATPRODUCT_AB)); 1509566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 1519566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 1529566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 1533fbca975SStefano Zampini } 1543fbca975SStefano Zampini 1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1589566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 159b122ec5aSJacob Faibussowitsch return 0; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /*TEST 163c4762a1bSJed Brown 164c4762a1bSJed Brown test: 165c4762a1bSJed Brown args: -M 10 -N 10 166c4762a1bSJed Brown output_file: output/ex109.out 167c4762a1bSJed Brown 168c4762a1bSJed Brown test: 169c4762a1bSJed Brown suffix: 2 170c4762a1bSJed Brown nsize: 3 171c4762a1bSJed Brown output_file: output/ex109.out 172c4762a1bSJed Brown 173c20d7725SJed Brown test: 174c20d7725SJed Brown suffix: 3 175c20d7725SJed Brown nsize: 2 176c20d7725SJed Brown args: -matmattransmult_mpidense_mpidense_via cyclic 177c20d7725SJed Brown output_file: output/ex109.out 178c20d7725SJed Brown 1793fbca975SStefano Zampini test: 1803fbca975SStefano Zampini suffix: 4 1813fbca975SStefano Zampini args: -test_userAPI 1823fbca975SStefano Zampini output_file: output/ex109.out 1833fbca975SStefano Zampini 1843fbca975SStefano Zampini test: 1853fbca975SStefano Zampini suffix: 5 1863fbca975SStefano Zampini nsize: 3 1873fbca975SStefano Zampini args: -test_userAPI 1883fbca975SStefano Zampini output_file: output/ex109.out 1893fbca975SStefano Zampini 190c4762a1bSJed Brown TEST*/ 191