1544b287aSSatish Balay static char help[] = "Test MatCreateSubMatrices\n\n"; 2544b287aSSatish Balay 3544b287aSSatish Balay #include <petscis.h> 4544b287aSSatish Balay #include <petscmat.h> 5544b287aSSatish Balay 69371c9d4SSatish Balay int main(int argc, char **args) { 7544b287aSSatish Balay Mat A, *submats, *submats2; 8544b287aSSatish Balay IS *irow, *icol; 9544b287aSSatish Balay PetscInt i, n; 10544b287aSSatish Balay PetscMPIInt rank; 11544b287aSSatish Balay PetscViewer matfd, rowfd, colfd; 12544b287aSSatish Balay PetscBool same; 13544b287aSSatish Balay char matfile[PETSC_MAX_PATH_LEN], rowfile[PETSC_MAX_PATH_LEN], colfile[PETSC_MAX_PATH_LEN]; 14544b287aSSatish Balay char rankstr[16] = {0}; 15544b287aSSatish Balay 16327415f7SBarry Smith PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 19544b287aSSatish Balay 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-A", matfile, sizeof(matfile), NULL)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-row", rowfile, sizeof(rowfile), NULL)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-col", colfile, sizeof(colfile), NULL)); 23544b287aSSatish Balay 24544b287aSSatish Balay /* Each rank has its own files for row/col ISes */ 259566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(rankstr, 16, "-%d", rank)); 269566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(rowfile, rankstr, PETSC_MAX_PATH_LEN)); 279566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(colfile, rankstr, PETSC_MAX_PATH_LEN)); 28544b287aSSatish Balay 299566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, matfile, FILE_MODE_READ, &matfd)); 309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, rowfile, FILE_MODE_READ, &rowfd)); 319566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, colfile, FILE_MODE_READ, &colfd)); 32544b287aSSatish Balay 339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 359566063dSJacob Faibussowitsch PetscCall(MatLoad(A, matfd)); 36544b287aSSatish Balay 37544b287aSSatish Balay /* We stored the number of ISes at the beginning of rowfd */ 389566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(rowfd, &n, 1, NULL, PETSC_INT)); 399566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &irow, n, &icol)); 40544b287aSSatish Balay for (i = 0; i < n; i++) { 419566063dSJacob Faibussowitsch PetscCall(ISCreate(PETSC_COMM_SELF, &irow[i])); 429566063dSJacob Faibussowitsch PetscCall(ISCreate(PETSC_COMM_SELF, &icol[i])); 439566063dSJacob Faibussowitsch PetscCall(ISLoad(irow[i], rowfd)); 449566063dSJacob Faibussowitsch PetscCall(ISLoad(icol[i], colfd)); 45544b287aSSatish Balay } 46544b287aSSatish Balay 479566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&matfd)); 489566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&rowfd)); 499566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&colfd)); 50544b287aSSatish Balay 51544b287aSSatish Balay /* Create submats for the first time */ 529566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, n, irow, icol, MAT_INITIAL_MATRIX, &submats)); 53544b287aSSatish Balay 54544b287aSSatish Balay /* Dup submats to submats2 for later comparison */ 559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &submats2)); 56*48a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(MatDuplicate(submats[i], MAT_COPY_VALUES, &submats2[i])); 57544b287aSSatish Balay 58544b287aSSatish Balay /* Create submats again */ 599566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, n, irow, icol, MAT_REUSE_MATRIX, &submats)); 60544b287aSSatish Balay 61544b287aSSatish Balay /* Compare submats and submats2 */ 62544b287aSSatish Balay for (i = 0; i < n; i++) { 639566063dSJacob Faibussowitsch PetscCall(MatEqual(submats[i], submats2[i], &same)); 6428b400f6SJacob Faibussowitsch PetscCheck(same, PETSC_COMM_SELF, PETSC_ERR_PLIB, "submatrix %" PetscInt_FMT " is not same", i); 65544b287aSSatish Balay } 66544b287aSSatish Balay 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 68544b287aSSatish Balay for (i = 0; i < n; i++) { 699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow[i])); 709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol[i])); 71544b287aSSatish Balay } 729566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(n, &submats)); 739566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(n, &submats2)); 749566063dSJacob Faibussowitsch PetscCall(PetscFree2(irow, icol)); 759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 76b122ec5aSJacob Faibussowitsch return 0; 77544b287aSSatish Balay } 78544b287aSSatish Balay 79544b287aSSatish Balay /*TEST 80544b287aSSatish Balay 81544b287aSSatish Balay test: 82544b287aSSatish Balay suffix: 1 83544b287aSSatish Balay nsize: 2 84dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 85544b287aSSatish Balay args: -mat_type {{aij baij}} -A ${DATAFILESPATH}/matrices/CreateSubMatrices/A -row ${DATAFILESPATH}/matrices/CreateSubMatrices/row -col ${DATAFILESPATH}/matrices/CreateSubMatrices/col 86544b287aSSatish Balay 87544b287aSSatish Balay TEST*/ 88