1544b287aSSatish Balay static char help[] = "Test MatCreateSubMatrices\n\n";
2544b287aSSatish Balay
3544b287aSSatish Balay #include <petscis.h>
4544b287aSSatish Balay #include <petscmat.h>
5544b287aSSatish Balay
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8544b287aSSatish Balay Mat A, *submats, *submats2;
9544b287aSSatish Balay IS *irow, *icol;
10544b287aSSatish Balay PetscInt i, n;
11544b287aSSatish Balay PetscMPIInt rank;
12544b287aSSatish Balay PetscViewer matfd, rowfd, colfd;
13544b287aSSatish Balay PetscBool same;
14544b287aSSatish Balay char matfile[PETSC_MAX_PATH_LEN], rowfile[PETSC_MAX_PATH_LEN], colfile[PETSC_MAX_PATH_LEN];
15544b287aSSatish Balay char rankstr[16] = {0};
16544b287aSSatish Balay
17327415f7SBarry Smith PetscFunctionBeginUser;
18c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
20544b287aSSatish Balay
219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-A", matfile, sizeof(matfile), NULL));
229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-row", rowfile, sizeof(rowfile), NULL));
239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-col", colfile, sizeof(colfile), NULL));
24544b287aSSatish Balay
25544b287aSSatish Balay /* Each rank has its own files for row/col ISes */
269566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(rankstr, 16, "-%d", rank));
279566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(rowfile, rankstr, PETSC_MAX_PATH_LEN));
289566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(colfile, rankstr, PETSC_MAX_PATH_LEN));
29544b287aSSatish Balay
309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, matfile, FILE_MODE_READ, &matfd));
319566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, rowfile, FILE_MODE_READ, &rowfd));
329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, colfile, FILE_MODE_READ, &colfd));
33544b287aSSatish Balay
349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
359566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
369566063dSJacob Faibussowitsch PetscCall(MatLoad(A, matfd));
37544b287aSSatish Balay
38544b287aSSatish Balay /* We stored the number of ISes at the beginning of rowfd */
399566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(rowfd, &n, 1, NULL, PETSC_INT));
409566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &irow, n, &icol));
41544b287aSSatish Balay for (i = 0; i < n; i++) {
429566063dSJacob Faibussowitsch PetscCall(ISCreate(PETSC_COMM_SELF, &irow[i]));
439566063dSJacob Faibussowitsch PetscCall(ISCreate(PETSC_COMM_SELF, &icol[i]));
449566063dSJacob Faibussowitsch PetscCall(ISLoad(irow[i], rowfd));
459566063dSJacob Faibussowitsch PetscCall(ISLoad(icol[i], colfd));
46544b287aSSatish Balay }
47544b287aSSatish Balay
489566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&matfd));
499566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&rowfd));
509566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&colfd));
51544b287aSSatish Balay
52544b287aSSatish Balay /* Create submats for the first time */
539566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, n, irow, icol, MAT_INITIAL_MATRIX, &submats));
54544b287aSSatish Balay
55544b287aSSatish Balay /* Dup submats to submats2 for later comparison */
569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &submats2));
5748a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(MatDuplicate(submats[i], MAT_COPY_VALUES, &submats2[i]));
58544b287aSSatish Balay
59544b287aSSatish Balay /* Create submats again */
609566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, n, irow, icol, MAT_REUSE_MATRIX, &submats));
61544b287aSSatish Balay
62544b287aSSatish Balay /* Compare submats and submats2 */
63544b287aSSatish Balay for (i = 0; i < n; i++) {
649566063dSJacob Faibussowitsch PetscCall(MatEqual(submats[i], submats2[i], &same));
6528b400f6SJacob Faibussowitsch PetscCheck(same, PETSC_COMM_SELF, PETSC_ERR_PLIB, "submatrix %" PetscInt_FMT " is not same", i);
66544b287aSSatish Balay }
67544b287aSSatish Balay
689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
69544b287aSSatish Balay for (i = 0; i < n; i++) {
709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow[i]));
719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol[i]));
72544b287aSSatish Balay }
739566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(n, &submats));
749566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(n, &submats2));
759566063dSJacob Faibussowitsch PetscCall(PetscFree2(irow, icol));
769566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
77b122ec5aSJacob Faibussowitsch return 0;
78544b287aSSatish Balay }
79544b287aSSatish Balay
80544b287aSSatish Balay /*TEST
81544b287aSSatish Balay
82544b287aSSatish Balay test:
83544b287aSSatish Balay suffix: 1
84544b287aSSatish Balay nsize: 2
85dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
86544b287aSSatish Balay args: -mat_type {{aij baij}} -A ${DATAFILESPATH}/matrices/CreateSubMatrices/A -row ${DATAFILESPATH}/matrices/CreateSubMatrices/row -col ${DATAFILESPATH}/matrices/CreateSubMatrices/col
87*3886731fSPierre Jolivet output_file: output/empty.out
88544b287aSSatish Balay
89544b287aSSatish Balay TEST*/
90