1c4762a1bSJed Brown static char help[] = "Tests sequential and parallel MatMatMatMult() and MatPtAP(). Modified from ex96.c \n\ 2c4762a1bSJed Brown -Mx <xg>, where <xg> = number of coarse grid points in the x-direction\n\ 3c4762a1bSJed Brown -My <yg>, where <yg> = number of coarse grid points in the y-direction\n\ 4c4762a1bSJed Brown -Mz <zg>, where <zg> = number of coarse grid points in the z-direction\n\ 5c4762a1bSJed Brown -Npx <npx>, where <npx> = number of processors in the x-direction\n\ 6c4762a1bSJed Brown -Npy <npy>, where <npy> = number of processors in the y-direction\n\ 7c4762a1bSJed Brown -Npz <npz>, where <npz> = number of processors in the z-direction\n\n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Example of usage: mpiexec -n 3 ./ex41 -Mx 10 -My 10 -Mz 10 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscdm.h> 14c4762a1bSJed Brown #include <petscdmda.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* User-defined application contexts */ 17c4762a1bSJed Brown typedef struct { 18c4762a1bSJed Brown PetscInt mx, my, mz; /* number grid points in x, y and z direction */ 19c4762a1bSJed Brown Vec localX, localF; /* local vectors with ghost region */ 20c4762a1bSJed Brown DM da; 21c4762a1bSJed Brown Vec x, b, r; /* global vectors */ 22c4762a1bSJed Brown Mat J; /* Jacobian on grid */ 23c4762a1bSJed Brown } GridCtx; 24c4762a1bSJed Brown typedef struct { 25c4762a1bSJed Brown GridCtx fine; 26c4762a1bSJed Brown GridCtx coarse; 27c4762a1bSJed Brown PetscInt ratio; 28c4762a1bSJed Brown Mat Ii; /* interpolation from coarse to fine */ 29c4762a1bSJed Brown } AppCtx; 30c4762a1bSJed Brown 31c4762a1bSJed Brown #define COARSE_LEVEL 0 32c4762a1bSJed Brown #define FINE_LEVEL 1 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown Mm_ratio - ration of grid lines between fine and coarse grids. 36c4762a1bSJed Brown */ 37*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 38*d71ae5a4SJacob Faibussowitsch { 39c4762a1bSJed Brown AppCtx user; 40c4762a1bSJed Brown PetscMPIInt size, rank; 41c4762a1bSJed Brown PetscInt m, n, M, N, i, nrows; 42c4762a1bSJed Brown PetscScalar one = 1.0; 43c4762a1bSJed Brown PetscReal fill = 2.0; 44c20d7725SJed Brown Mat A, P, R, C, PtAP, D; 45c4762a1bSJed Brown PetscScalar *array; 46c4762a1bSJed Brown PetscRandom rdm; 47c4762a1bSJed Brown PetscBool Test_3D = PETSC_FALSE, flg; 48c4762a1bSJed Brown const PetscInt *ia, *ja; 49c4762a1bSJed Brown 50327415f7SBarry Smith PetscFunctionBeginUser; 519566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Get size of fine grids and coarse grids */ 56c4762a1bSJed Brown user.ratio = 2; 579371c9d4SSatish Balay user.coarse.mx = 4; 589371c9d4SSatish Balay user.coarse.my = 4; 599371c9d4SSatish Balay user.coarse.mz = 4; 60c4762a1bSJed Brown 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-Mx", &user.coarse.mx, NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-My", &user.coarse.my, NULL)); 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-Mz", &user.coarse.mz, NULL)); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ratio", &user.ratio, NULL)); 65c4762a1bSJed Brown if (user.coarse.mz) Test_3D = PETSC_TRUE; 66c4762a1bSJed Brown 67c4762a1bSJed Brown user.fine.mx = user.ratio * (user.coarse.mx - 1) + 1; 68c4762a1bSJed Brown user.fine.my = user.ratio * (user.coarse.my - 1) + 1; 69c4762a1bSJed Brown user.fine.mz = user.ratio * (user.coarse.mz - 1) + 1; 70c4762a1bSJed Brown 71dd400576SPatrick Sanan if (rank == 0) { 72c4762a1bSJed Brown if (!Test_3D) { 739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "coarse grids: %" PetscInt_FMT " %" PetscInt_FMT "; fine grids: %" PetscInt_FMT " %" PetscInt_FMT "\n", user.coarse.mx, user.coarse.my, user.fine.mx, user.fine.my)); 74c4762a1bSJed Brown } else { 759371c9d4SSatish Balay PetscCall(PetscPrintf(PETSC_COMM_SELF, "coarse grids: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "; fine grids: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", user.coarse.mx, user.coarse.my, user.coarse.mz, user.fine.mx, 769371c9d4SSatish Balay user.fine.my, user.fine.mz)); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Set up distributed array for fine grid */ 81c4762a1bSJed Brown if (!Test_3D) { 829566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.fine.mx, user.fine.my, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &user.fine.da)); 83c4762a1bSJed Brown } else { 849371c9d4SSatish Balay PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.fine.mx, user.fine.my, user.fine.mz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, NULL, &user.fine.da)); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.fine.da)); 879566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.fine.da)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* Create and set A at fine grids */ 909566063dSJacob Faibussowitsch PetscCall(DMSetMatType(user.fine.da, MATAIJ)); 919566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.fine.da, &A)); 929566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 939566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* set val=one to A (replace with random values!) */ 969566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 979566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 98c4762a1bSJed Brown if (size == 1) { 999566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 100c4762a1bSJed Brown if (flg) { 1019566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(A, &array)); 102c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one; 1039566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(A, &array)); 104c4762a1bSJed Brown } 1059566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 106c4762a1bSJed Brown } else { 107c4762a1bSJed Brown Mat AA, AB; 1089566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(A, &AA, &AB, NULL)); 1099566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 110c4762a1bSJed Brown if (flg) { 1119566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(AA, &array)); 112c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one; 1139566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(AA, &array)); 114c4762a1bSJed Brown } 1159566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(AA, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 1169566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 117c4762a1bSJed Brown if (flg) { 1189566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArray(AB, &array)); 119c4762a1bSJed Brown for (i = 0; i < ia[nrows]; i++) array[i] = one; 1209566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArray(AB, &array)); 121c4762a1bSJed Brown } 1229566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(AB, 0, PETSC_FALSE, PETSC_FALSE, &nrows, &ia, &ja, &flg)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown /* Set up distributed array for coarse grid */ 125c4762a1bSJed Brown if (!Test_3D) { 1269566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.coarse.mx, user.coarse.my, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &user.coarse.da)); 127c4762a1bSJed Brown } else { 1289566063dSJacob Faibussowitsch PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.coarse.mx, user.coarse.my, user.coarse.mz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, NULL, &user.coarse.da)); 129c4762a1bSJed Brown } 1309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.coarse.da)); 1319566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.coarse.da)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 1349566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(user.coarse.da, user.fine.da, &P, NULL)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* Get R = P^T */ 1379566063dSJacob Faibussowitsch PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &R)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* C = R*A*P */ 140c20d7725SJed Brown /* Developer's API */ 1419566063dSJacob Faibussowitsch PetscCall(MatProductCreate(R, A, P, &D)); 1429566063dSJacob Faibussowitsch PetscCall(MatProductSetType(D, MATPRODUCT_ABC)); 1439566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(D)); 1449566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(D)); 1459566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(D)); 1469566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(D)); /* Test reuse symbolic D */ 147c20d7725SJed Brown 148c20d7725SJed Brown /* User's API */ 149c20d7725SJed Brown { /* Test MatMatMatMult_Basic() */ 150c20d7725SJed Brown Mat Adense, Cdense; 1519566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &Adense)); 1529566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, Adense, P, MAT_INITIAL_MATRIX, fill, &Cdense)); 1539566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, Adense, P, MAT_REUSE_MATRIX, fill, &Cdense)); 154c20d7725SJed Brown 1559566063dSJacob Faibussowitsch PetscCall(MatMultEqual(D, Cdense, 10, &flg)); 15628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "D*v != Cdense*v"); 1579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adense)); 1589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense)); 159c20d7725SJed Brown } 160c20d7725SJed Brown 1619566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, A, P, MAT_INITIAL_MATRIX, fill, &C)); 1629566063dSJacob Faibussowitsch PetscCall(MatMatMatMult(R, A, P, MAT_REUSE_MATRIX, fill, &C)); 1639566063dSJacob Faibussowitsch PetscCall(MatProductClear(C)); 164c4762a1bSJed Brown 165c20d7725SJed Brown /* Test D == C */ 1669566063dSJacob Faibussowitsch PetscCall(MatEqual(D, C, &flg)); 16728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "D != C"); 168c20d7725SJed Brown 169c4762a1bSJed Brown /* Test C == PtAP */ 1709566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, fill, &PtAP)); 1719566063dSJacob Faibussowitsch PetscCall(MatPtAP(A, P, MAT_REUSE_MATRIX, fill, &PtAP)); 1729566063dSJacob Faibussowitsch PetscCall(MatEqual(C, PtAP, &flg)); 17328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "C != PtAP"); 1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PtAP)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Clean up */ 1779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1789566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.fine.da)); 1809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.coarse.da)); 1819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 1829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 1839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 1859566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 186b122ec5aSJacob Faibussowitsch return 0; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /*TEST 190c4762a1bSJed Brown 191c4762a1bSJed Brown test: 192c4762a1bSJed Brown 193c4762a1bSJed Brown test: 194c4762a1bSJed Brown suffix: 2 195c4762a1bSJed Brown nsize: 2 196c4762a1bSJed Brown args: -matmatmatmult_via scalable 197c4762a1bSJed Brown 198c4762a1bSJed Brown test: 199c4762a1bSJed Brown suffix: 3 200c4762a1bSJed Brown nsize: 2 201c4762a1bSJed Brown args: -matmatmatmult_via nonscalable 202c4762a1bSJed Brown output_file: output/ex111_1.out 203c4762a1bSJed Brown 204c4762a1bSJed Brown TEST*/ 205