1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests sequential and parallel MatMatMult() and MatPtAP(), MatTransposeMatMult(), sequential MatMatTransposeMult(), MatRARt()\n\ 3c4762a1bSJed Brown Input arguments are:\n\ 4c4762a1bSJed Brown -f0 <input_file> -f1 <input_file> -f2 <input_file> -f3 <input_file> : file to load\n\n"; 5c4762a1bSJed Brown /* Example of usage: 6c4762a1bSJed Brown ./ex94 -f0 <A_binary> -f1 <B_binary> -matmatmult_mat_view ascii::ascii_info -matmatmulttr_mat_view 7c4762a1bSJed Brown mpiexec -n 3 ./ex94 -f0 medium -f1 medium -f2 arco1 -f3 arco1 -matmatmult_mat_view 8c4762a1bSJed Brown */ 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown B = A - B 14c4762a1bSJed Brown norm = norm(B) 15c4762a1bSJed Brown */ 16d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNormDifference(Mat A, Mat B, PetscReal *norm) 17d71ae5a4SJacob Faibussowitsch { 18c4762a1bSJed Brown PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(MatAXPY(B, -1.0, A, DIFFERENT_NONZERO_PATTERN)); 209566063dSJacob Faibussowitsch PetscCall(MatNorm(B, NORM_FROBENIUS, norm)); 213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22c4762a1bSJed Brown } 23c4762a1bSJed Brown 24d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown Mat A, A_save, B, AT, ATT, BT, BTT, P, R, C, C1; 27c4762a1bSJed Brown Vec x, v1, v2, v3, v4; 28c4762a1bSJed Brown PetscViewer viewer; 29c4762a1bSJed Brown PetscMPIInt size, rank; 30c4762a1bSJed Brown PetscInt i, m, n, j, *idxn, M, N, nzp, rstart, rend; 31c4762a1bSJed Brown PetscReal norm, norm_abs, norm_tmp, fill = 4.0; 32c4762a1bSJed Brown PetscRandom rdm; 33c4762a1bSJed Brown char file[4][128]; 34c4762a1bSJed Brown PetscBool flg, preload = PETSC_TRUE; 35c4762a1bSJed Brown PetscScalar *a, rval, alpha, none = -1.0; 36c4762a1bSJed Brown PetscBool Test_MatMatMult = PETSC_TRUE, Test_MatMatTr = PETSC_TRUE, Test_MatPtAP = PETSC_TRUE, Test_MatRARt = PETSC_TRUE, Test_MatMatMatMult = PETSC_TRUE; 37c4762a1bSJed Brown PetscBool Test_MatAXPY = PETSC_FALSE, view = PETSC_FALSE; 38c4762a1bSJed Brown PetscInt pm, pn, pM, pN; 39c4762a1bSJed Brown MatInfo info; 40c4762a1bSJed Brown PetscBool seqaij; 41c4762a1bSJed Brown MatType mattype; 42c4762a1bSJed Brown Mat Cdensetest, Pdense, Cdense, Adense; 43c4762a1bSJed Brown 44327415f7SBarry Smith PetscFunctionBeginUser; 459566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 48c4762a1bSJed Brown 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL)); 509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-matops_view", &view, NULL)); 519566063dSJacob Faibussowitsch if (view) PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Load the matrices A_save and B */ 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f0", file[0], sizeof(file[0]), &flg)); 5528b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate a file name for small matrix A with the -f0 option."); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f1", file[1], sizeof(file[1]), &flg)); 5728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate a file name for small matrix B with the -f1 option."); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f2", file[2], sizeof(file[2]), &flg)); 59c4762a1bSJed Brown if (!flg) { 60c4762a1bSJed Brown preload = PETSC_FALSE; 61c4762a1bSJed Brown } else { 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f3", file[3], sizeof(file[3]), &flg)); 6328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate a file name for test matrix B with the -f3 option."); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscPreLoadBegin(preload, "Load system"); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[2 * PetscPreLoadIt], FILE_MODE_READ, &viewer)); 689566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A_save)); 699566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A_save)); 709566063dSJacob Faibussowitsch PetscCall(MatLoad(A_save, viewer)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[2 * PetscPreLoadIt + 1], FILE_MODE_READ, &viewer)); 749566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 759566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 769566063dSJacob Faibussowitsch PetscCall(MatLoad(B, viewer)); 779566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(MatGetType(B, &mattype)); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, &M, &N)); 82c4762a1bSJed Brown nzp = PetscMax((PetscInt)(0.1 * M), 5); 839566063dSJacob Faibussowitsch PetscCall(PetscMalloc((nzp + 1) * (sizeof(PetscInt) + sizeof(PetscScalar)), &idxn)); 84c4762a1bSJed Brown a = (PetscScalar *)(idxn + nzp); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Create vectors v1 and v2 that are compatible with A_save */ 879566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &v1)); 889566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A_save, &m, NULL)); 899566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v1, m, PETSC_DECIDE)); 909566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(v1)); 919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v1, &v2)); 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 949566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Test MatAXPY() */ 98c4762a1bSJed Brown /*-------------------*/ 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_MatAXPY", &Test_MatAXPY)); 100c4762a1bSJed Brown if (Test_MatAXPY) { 101c4762a1bSJed Brown Mat Btmp; 1029566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A_save, MAT_COPY_VALUES, &A)); 1039566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &Btmp)); 1049566063dSJacob Faibussowitsch PetscCall(MatAXPY(A, -1.0, B, DIFFERENT_NONZERO_PATTERN)); /* A = -B + A_save */ 105c4762a1bSJed Brown 1069566063dSJacob Faibussowitsch PetscCall(MatScale(A, -1.0)); /* A = -A = B - A_save */ 1079566063dSJacob Faibussowitsch PetscCall(MatAXPY(Btmp, -1.0, A, DIFFERENT_NONZERO_PATTERN)); /* Btmp = -A + B = A_save */ 1089566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A_save, Btmp, 10, &flg)); 10928b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatAXPY() is incorrect"); 1109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Btmp)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown Test_MatMatMult = PETSC_FALSE; 114c4762a1bSJed Brown Test_MatMatTr = PETSC_FALSE; 115c4762a1bSJed Brown Test_MatPtAP = PETSC_FALSE; 116c4762a1bSJed Brown Test_MatRARt = PETSC_FALSE; 117c4762a1bSJed Brown Test_MatMatMatMult = PETSC_FALSE; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* 1) Test MatMatMult() */ 121c4762a1bSJed Brown /* ---------------------*/ 122c4762a1bSJed Brown if (Test_MatMatMult) { 1239566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A_save, MAT_COPY_VALUES, &A)); 1249566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A, &AT)); 1259566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(AT, &ATT)); 1269566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(B, &BT)); 1279566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(BT, &BTT)); 128c20d7725SJed Brown 1299566063dSJacob Faibussowitsch PetscCall(MatMatMult(AT, B, MAT_INITIAL_MATRIX, fill, &C)); 1309566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(AT, B, C, 10, &flg)); 13128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for C=AT*B"); 1329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 133c20d7725SJed Brown 1349566063dSJacob Faibussowitsch PetscCall(MatMatMult(ATT, B, MAT_INITIAL_MATRIX, fill, &C)); 1359566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(ATT, B, C, 10, &flg)); 13628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for C=ATT*B"); 1379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 138c20d7725SJed Brown 1399566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, fill, &C)); 1409566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &flg)); 14128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for reuse C=A*B"); 142c20d7725SJed Brown /* ATT has different matrix type as A (although they have same internal data structure), 143c20d7725SJed Brown we cannot call MatProductReplaceMats(ATT,NULL,NULL,C) and MatMatMult(ATT,B,MAT_REUSE_MATRIX,fill,&C) */ 1449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 145c20d7725SJed Brown 1469566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, BTT, MAT_INITIAL_MATRIX, fill, &C)); 1479566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, BTT, C, 10, &flg)); 14828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult() for C=A*BTT"); 1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 150c20d7725SJed Brown 1519566063dSJacob Faibussowitsch PetscCall(MatMatMult(ATT, BTT, MAT_INITIAL_MATRIX, fill, &C)); 1529566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &flg)); 15328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult()"); 1549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 155c20d7725SJed Brown 1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BTT)); 1579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 1589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ATT)); 1599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AT)); 160c20d7725SJed Brown 1619566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, fill, &C)); 1629566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C, "matmatmult_")); /* enable option '-matmatmult_' for matrix C */ 1639566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 166c4762a1bSJed Brown alpha = 1.0; 167c4762a1bSJed Brown for (i = 0; i < 2; i++) { 168c4762a1bSJed Brown alpha -= 0.1; 1699566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha)); 1709566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, fill, &C)); 171c4762a1bSJed Brown } 1729566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A, B, C, 10, &flg)); 17328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult()"); 1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Test MatDuplicate() of C=A*B */ 1779566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &C1)); 1789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 1799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 180c4762a1bSJed Brown } /* if (Test_MatMatMult) */ 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* 2) Test MatTransposeMatMult() and MatMatTransposeMult() */ 183c4762a1bSJed Brown /* ------------------------------------------------------- */ 184c4762a1bSJed Brown if (Test_MatMatTr) { 185c4762a1bSJed Brown /* Create P */ 186c4762a1bSJed Brown PetscInt PN, rstart, rend; 187c4762a1bSJed Brown PN = M / 2; 188c4762a1bSJed Brown nzp = 5; /* num of nonzeros in each row of P */ 1899566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &P)); 1909566063dSJacob Faibussowitsch PetscCall(MatSetSizes(P, PETSC_DECIDE, PETSC_DECIDE, M, PN)); 1919566063dSJacob Faibussowitsch PetscCall(MatSetType(P, mattype)); 1929566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(P, nzp, NULL)); 1939566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(P, nzp, NULL, nzp, NULL)); 1949566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(P, &rstart, &rend)); 19548a46eb9SPierre Jolivet for (i = 0; i < nzp; i++) PetscCall(PetscRandomGetValue(rdm, &a[i])); 196c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 197c4762a1bSJed Brown for (j = 0; j < nzp; j++) { 1989566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 199c4762a1bSJed Brown idxn[j] = (PetscInt)(PetscRealPart(rval) * PN); 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(MatSetValues(P, 1, &i, nzp, idxn, a, ADD_VALUES)); 202c4762a1bSJed Brown } 2039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 2049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* Create R = P^T */ 2079566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &R)); 208c4762a1bSJed Brown 209c4762a1bSJed Brown { /* Test R = P^T, C1 = R*B */ 2109566063dSJacob Faibussowitsch PetscCall(MatMatMult(R, B, MAT_INITIAL_MATRIX, fill, &C1)); 2119566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_REUSE_MATRIX, &R)); 2129566063dSJacob Faibussowitsch PetscCall(MatMatMult(R, B, MAT_REUSE_MATRIX, fill, &C1)); 2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* C = P^T*B */ 2179566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(P, B, MAT_INITIAL_MATRIX, fill, &C)); 2189566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info)); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 2219566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(P, B, MAT_REUSE_MATRIX, fill, &C)); 222c4762a1bSJed Brown if (view) { 2239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "C = P^T * B:\n")); 2249566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 225c4762a1bSJed Brown } 2269566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 227c4762a1bSJed Brown if (view) { 2289566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nC = P^T * B after MatProductClear():\n")); 2299566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* Compare P^T*B and R*B */ 2339566063dSJacob Faibussowitsch PetscCall(MatMatMult(R, B, MAT_INITIAL_MATRIX, fill, &C1)); 2349566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, C1, &norm)); 23508401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatTransposeMatMult(): %g", (double)norm); 2369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* Test MatDuplicate() of C=P^T*B */ 2399566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &C1)); 2409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 2419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 242c4762a1bSJed Brown 243c4762a1bSJed Brown /* C = B*R^T */ 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &seqaij)); 245c4762a1bSJed Brown if (size == 1 && seqaij) { 2469566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(B, R, MAT_INITIAL_MATRIX, fill, &C)); 2479566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C, "matmatmulttr_")); /* enable '-matmatmulttr_' for matrix C */ 2489566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info)); 249c4762a1bSJed Brown 250c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 2519566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMult(B, R, MAT_REUSE_MATRIX, fill, &C)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* Check */ 2549566063dSJacob Faibussowitsch PetscCall(MatMatMult(B, P, MAT_INITIAL_MATRIX, fill, &C1)); 2559566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, C1, &norm)); 25608401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatMatTransposeMult() %g", (double)norm); 2579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 2589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 259c4762a1bSJed Brown } 2609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 262c4762a1bSJed Brown } 263c4762a1bSJed Brown 264c4762a1bSJed Brown /* 3) Test MatPtAP() */ 265c4762a1bSJed Brown /*-------------------*/ 266c4762a1bSJed Brown if (Test_MatPtAP) { 267c4762a1bSJed Brown PetscInt PN; 268c4762a1bSJed Brown Mat Cdup; 269c4762a1bSJed Brown 2709566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A_save, MAT_COPY_VALUES, &A)); 2719566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 2729566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 273c4762a1bSJed Brown 274c4762a1bSJed Brown PN = M / 2; 275*bbea24aaSStefano Zampini nzp = (PetscInt)(0.1 * PN + 1); /* num of nonzeros in each row of P */ 2769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &P)); 2779566063dSJacob Faibussowitsch PetscCall(MatSetSizes(P, PETSC_DECIDE, PETSC_DECIDE, N, PN)); 2789566063dSJacob Faibussowitsch PetscCall(MatSetType(P, mattype)); 2799566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(P, nzp, NULL)); 2809566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(P, nzp, NULL, nzp, NULL)); 28148a46eb9SPierre Jolivet for (i = 0; i < nzp; i++) PetscCall(PetscRandomGetValue(rdm, &a[i])); 2829566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(P, &rstart, &rend)); 283c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 284c4762a1bSJed Brown for (j = 0; j < nzp; j++) { 2859566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval)); 286c4762a1bSJed Brown idxn[j] = (PetscInt)(PetscRealPart(rval) * PN); 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(MatSetValues(P, 1, &i, nzp, idxn, a, ADD_VALUES)); 289c4762a1bSJed Brown } 2909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY)); 2919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY)); 292c4762a1bSJed Brown 2939566063dSJacob Faibussowitsch /* PetscCall(MatView(P,PETSC_VIEWER_STDOUT_WORLD)); */ 2949566063dSJacob Faibussowitsch PetscCall(MatGetSize(P, &pM, &pN)); 2959566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(P, &pm, &pn)); 2969566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, fill, &C)); 297c4762a1bSJed Brown 298c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 299c4762a1bSJed Brown alpha = 1.0; 300c4762a1bSJed Brown for (i = 0; i < 2; i++) { 301c4762a1bSJed Brown alpha -= 0.1; 3029566063dSJacob Faibussowitsch PetscCall(MatScale(A, alpha)); 3039566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, fill, &C)); 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306c4762a1bSJed Brown /* Test PtAP ops with P Dense and A either AIJ or SeqDense (it assumes MatPtAP_XAIJ_XAIJ is fine) */ 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &seqaij)); 308c4762a1bSJed Brown if (seqaij) { 3099566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cdensetest)); 3109566063dSJacob Faibussowitsch PetscCall(MatConvert(P, MATSEQDENSE, MAT_INITIAL_MATRIX, &Pdense)); 311c4762a1bSJed Brown } else { 3129566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATMPIDENSE, MAT_INITIAL_MATRIX, &Cdensetest)); 3139566063dSJacob Faibussowitsch PetscCall(MatConvert(P, MATMPIDENSE, MAT_INITIAL_MATRIX, &Pdense)); 314c4762a1bSJed Brown } 315c4762a1bSJed Brown 316c20d7725SJed Brown /* test with A(AIJ), Pdense -- call MatPtAP_Basic() when np>1 */ 3179566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, Pdense, MAT_INITIAL_MATRIX, fill, &Cdense)); 3189566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, Pdense, MAT_REUSE_MATRIX, fill, &Cdense)); 3199566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(A, Pdense, Cdense, 10, &flg)); 32028b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in MatPtAP with A AIJ and P Dense"); 3219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense)); 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* test with A SeqDense */ 324c4762a1bSJed Brown if (seqaij) { 3259566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &Adense)); 3269566063dSJacob Faibussowitsch PetscCall(MatPtAP(Adense, Pdense, MAT_INITIAL_MATRIX, fill, &Cdense)); 3279566063dSJacob Faibussowitsch PetscCall(MatPtAP(Adense, Pdense, MAT_REUSE_MATRIX, fill, &Cdense)); 3289566063dSJacob Faibussowitsch PetscCall(MatPtAPMultEqual(Adense, Pdense, Cdense, 10, &flg)); 32928b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in MatPtAP with A SeqDense and P SeqDense"); 3309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense)); 3319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adense)); 332c4762a1bSJed Brown } 3339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdensetest)); 3349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Pdense)); 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* Test MatDuplicate() of C=PtAP and MatView(Cdup,...) */ 3379566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &Cdup)); 338c4762a1bSJed Brown if (view) { 3399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nC = P^T * A * P:\n")); 3409566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 341c4762a1bSJed Brown 3429566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 3439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nC = P^T * A * P after MatProductClear():\n")); 3449566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 345c4762a1bSJed Brown 3469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nCdup:\n")); 3479566063dSJacob Faibussowitsch PetscCall(MatView(Cdup, PETSC_VIEWER_STDOUT_WORLD)); 348c4762a1bSJed Brown } 3499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdup)); 350c4762a1bSJed Brown 351c4762a1bSJed Brown if (size > 1 || !seqaij) Test_MatRARt = PETSC_FALSE; 352c4762a1bSJed Brown /* 4) Test MatRARt() */ 353c4762a1bSJed Brown /* ----------------- */ 354c4762a1bSJed Brown if (Test_MatRARt) { 355c20d7725SJed Brown Mat R, RARt, Rdense, RARtdense; 3569566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &R)); 357c20d7725SJed Brown 358c20d7725SJed Brown /* Test MatRARt_Basic(), MatMatMatMult_Basic() */ 3599566063dSJacob Faibussowitsch PetscCall(MatConvert(R, MATDENSE, MAT_INITIAL_MATRIX, &Rdense)); 3609566063dSJacob Faibussowitsch PetscCall(MatRARt(A, Rdense, MAT_INITIAL_MATRIX, 2.0, &RARtdense)); 3619566063dSJacob Faibussowitsch PetscCall(MatRARt(A, Rdense, MAT_REUSE_MATRIX, 2.0, &RARtdense)); 362c20d7725SJed Brown 3639566063dSJacob Faibussowitsch PetscCall(MatConvert(RARtdense, MATAIJ, MAT_INITIAL_MATRIX, &RARt)); 3649566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, RARt, &norm)); 36508401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "|PtAP - RARtdense| = %g", (double)norm); 3669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Rdense)); 3679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RARtdense)); 3689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RARt)); 369c20d7725SJed Brown 370c20d7725SJed Brown /* Test MatRARt() for aij matrices */ 3719566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_INITIAL_MATRIX, 2.0, &RARt)); 3729566063dSJacob Faibussowitsch PetscCall(MatRARt(A, R, MAT_REUSE_MATRIX, 2.0, &RARt)); 3739566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, RARt, &norm)); 37408401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "|PtAP - RARt| = %g", (double)norm); 3759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 3769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RARt)); 377c4762a1bSJed Brown } 378c4762a1bSJed Brown 379c4762a1bSJed Brown if (Test_MatMatMatMult && size == 1) { 380c4762a1bSJed Brown Mat R, RAP; 3819566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &R)); 3829566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, A, P, MAT_INITIAL_MATRIX, 2.0, &RAP)); 3839566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, A, P, MAT_REUSE_MATRIX, 2.0, &RAP)); 3849566063dSJacob Faibussowitsch PetscCall(MatNormDifference(C, RAP, &norm)); 38508401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PtAP != RAP %g", (double)norm); 3869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 3879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RAP)); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* Create vector x that is compatible with P */ 3919566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 3929566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(P, &m, &n)); 3939566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, n, PETSC_DECIDE)); 3949566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 395c4762a1bSJed Brown 3969566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &v3)); 3979566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v3, n, PETSC_DECIDE)); 3989566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(v3)); 3999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v3, &v4)); 400c4762a1bSJed Brown 401c4762a1bSJed Brown norm = 0.0; 402c4762a1bSJed Brown for (i = 0; i < 10; i++) { 4039566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 4049566063dSJacob Faibussowitsch PetscCall(MatMult(P, x, v1)); 4059566063dSJacob Faibussowitsch PetscCall(MatMult(A, v1, v2)); /* v2 = A*P*x */ 406c4762a1bSJed Brown 4079566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(P, v2, v3)); /* v3 = Pt*A*P*x */ 4089566063dSJacob Faibussowitsch PetscCall(MatMult(C, x, v4)); /* v3 = C*x */ 4099566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_abs)); 4109566063dSJacob Faibussowitsch PetscCall(VecAXPY(v4, none, v3)); 4119566063dSJacob Faibussowitsch PetscCall(VecNorm(v4, NORM_2, &norm_tmp)); 412c4762a1bSJed Brown 413c4762a1bSJed Brown norm_tmp /= norm_abs; 414c4762a1bSJed Brown if (norm_tmp > norm) norm = norm_tmp; 415c4762a1bSJed Brown } 41608401ef6SPierre Jolivet PetscCheck(norm < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatPtAP(), |v1 - v2|: %g", (double)norm); 417c4762a1bSJed Brown 4189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 4199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 4209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 4219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v3)); 4229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v4)); 4239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 424c4762a1bSJed Brown } 425c4762a1bSJed Brown 426c4762a1bSJed Brown /* Destroy objects */ 4279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v1)); 4289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v2)); 4299566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 4309566063dSJacob Faibussowitsch PetscCall(PetscFree(idxn)); 431c4762a1bSJed Brown 4329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_save)); 4339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 434c4762a1bSJed Brown 435c4762a1bSJed Brown PetscPreLoadEnd(); 4369566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 437b122ec5aSJacob Faibussowitsch return 0; 438c4762a1bSJed Brown } 439c4762a1bSJed Brown 440c4762a1bSJed Brown /*TEST 441c4762a1bSJed Brown 442c4762a1bSJed Brown test: 443c4762a1bSJed Brown suffix: 2_mattransposematmult_matmatmult 444c4762a1bSJed Brown nsize: 3 445dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 446c20d7725SJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via at*b> ex94_2.tmp 2>&1 447c4762a1bSJed Brown 448c4762a1bSJed Brown test: 449c4762a1bSJed Brown suffix: 2_mattransposematmult_scalable 450c4762a1bSJed Brown nsize: 3 451dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 452c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via scalable> ex94_2.tmp 2>&1 453c4762a1bSJed Brown output_file: output/ex94_1.out 454c4762a1bSJed Brown 455c4762a1bSJed Brown test: 456c4762a1bSJed Brown suffix: axpy_mpiaij 457dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 458c4762a1bSJed Brown nsize: 8 459c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY 460c4762a1bSJed Brown output_file: output/ex94_1.out 461c4762a1bSJed Brown 462c4762a1bSJed Brown test: 463c4762a1bSJed Brown suffix: axpy_mpibaij 464dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 465c4762a1bSJed Brown nsize: 8 466c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type baij 467c4762a1bSJed Brown output_file: output/ex94_1.out 468c4762a1bSJed Brown 469c4762a1bSJed Brown test: 470c4762a1bSJed Brown suffix: axpy_mpisbaij 471dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 472c4762a1bSJed Brown nsize: 8 473c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type sbaij 474c4762a1bSJed Brown output_file: output/ex94_1.out 475c4762a1bSJed Brown 476c4762a1bSJed Brown test: 477c4762a1bSJed Brown suffix: matmatmult 478dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 479c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info 480c4762a1bSJed Brown output_file: output/ex94_1.out 481c4762a1bSJed Brown 482c4762a1bSJed Brown test: 483c4762a1bSJed Brown suffix: matmatmult_2 484dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 485c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -mat_type mpiaij -viewer_binary_skip_info 486c4762a1bSJed Brown output_file: output/ex94_1.out 487c4762a1bSJed Brown 488c4762a1bSJed Brown test: 489c4762a1bSJed Brown suffix: matmatmult_scalable 490c4762a1bSJed Brown nsize: 4 491dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 492c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -matmatmult_via scalable 493c4762a1bSJed Brown output_file: output/ex94_1.out 494c4762a1bSJed Brown 495c4762a1bSJed Brown test: 496c4762a1bSJed Brown suffix: ptap 497c4762a1bSJed Brown nsize: 3 498dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 499c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -matptap_via scalable 500c4762a1bSJed Brown output_file: output/ex94_1.out 501c4762a1bSJed Brown 502c4762a1bSJed Brown test: 503c4762a1bSJed Brown suffix: rap 504c4762a1bSJed Brown nsize: 3 505dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 506c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium 507c4762a1bSJed Brown output_file: output/ex94_1.out 508c4762a1bSJed Brown 509c4762a1bSJed Brown test: 510c4762a1bSJed Brown suffix: scalable0 511dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 512c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info 513c4762a1bSJed Brown output_file: output/ex94_1.out 514c4762a1bSJed Brown 515c4762a1bSJed Brown test: 516c4762a1bSJed Brown suffix: scalable1 517dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 518c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -matptap_via scalable 519c4762a1bSJed Brown output_file: output/ex94_1.out 520c4762a1bSJed Brown 521c4762a1bSJed Brown test: 522c4762a1bSJed Brown suffix: view 523c4762a1bSJed Brown nsize: 2 524dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 525c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/tiny -f1 ${DATAFILESPATH}/matrices/tiny -viewer_binary_skip_info -matops_view 526c4762a1bSJed Brown output_file: output/ex94_2.out 527c4762a1bSJed Brown 528c4762a1bSJed Brown TEST*/ 529