xref: /petsc/src/mat/tests/ex212.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
173b6653fSLawrence Mitchell static char help[] = "Test MatCreateSubMatrix with -mat_type nest and block sizes.\n";
273b6653fSLawrence Mitchell 
373b6653fSLawrence Mitchell #include <petscmat.h>
473b6653fSLawrence Mitchell 
59371c9d4SSatish Balay int main(int argc, char **argv) {
673b6653fSLawrence Mitchell   Mat                    A, B, C, mats[6];
773b6653fSLawrence Mitchell   IS                     rows[2];
873b6653fSLawrence Mitchell   ISLocalToGlobalMapping cmap, rmap;
973b6653fSLawrence Mitchell   const PetscInt         indices[3] = {0, 1, 2};
1073b6653fSLawrence Mitchell   PetscInt               i;
1173b6653fSLawrence Mitchell   PetscMPIInt            size;
1273b6653fSLawrence Mitchell 
13327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
16be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only coded for one process");
179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
199566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(A, 2, 1));
209566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
219566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 2, 1, indices, PETSC_COPY_VALUES, &rmap));
229566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 1, indices, PETSC_COPY_VALUES, &cmap));
239566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap));
249566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&rmap));
259566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&cmap));
269566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 2, indices, PETSC_COPY_VALUES, &rmap));
279566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 3, indices, PETSC_COPY_VALUES, &cmap));
289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, 2, 3));
309566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(B, 1, 1));
319566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATAIJ));
329566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(B, rmap, cmap));
339566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&rmap));
349566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&cmap));
359566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
369566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
3773b6653fSLawrence Mitchell   mats[0] = A;
3873b6653fSLawrence Mitchell   mats[1] = B;
3973b6653fSLawrence Mitchell   mats[2] = A;
4073b6653fSLawrence Mitchell   mats[3] = NULL;
4173b6653fSLawrence Mitchell   mats[4] = B;
4273b6653fSLawrence Mitchell   mats[5] = A;
439566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 3, NULL, mats, &C));
449566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
479566063dSJacob Faibussowitsch   PetscCall(MatView(C, NULL));
489566063dSJacob Faibussowitsch   PetscCall(MatNestGetISs(C, rows, NULL));
4973b6653fSLawrence Mitchell   for (i = 0; i < 2; i++) {
5073b6653fSLawrence Mitchell     Mat      submat;
5173b6653fSLawrence Mitchell     IS       cols[3];
5273b6653fSLawrence Mitchell     PetscInt j;
539566063dSJacob Faibussowitsch     PetscCall(ISView(rows[i], NULL));
549566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix(C, rows[i], NULL, MAT_INITIAL_MATRIX, &submat));
559566063dSJacob Faibussowitsch     PetscCall(MatView(submat, NULL));
569566063dSJacob Faibussowitsch     PetscCall(MatNestGetISs(submat, NULL, cols));
57*48a46eb9SPierre Jolivet     for (j = 0; j < 3; j++) PetscCall(ISView(cols[j], NULL));
589566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&submat));
5973b6653fSLawrence Mitchell   }
609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
639566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
64b122ec5aSJacob Faibussowitsch   return 0;
6573b6653fSLawrence Mitchell }
6673b6653fSLawrence Mitchell 
6773b6653fSLawrence Mitchell /*TEST
6873b6653fSLawrence Mitchell 
6973b6653fSLawrence Mitchell     test:
7073b6653fSLawrence Mitchell 
7173b6653fSLawrence Mitchell TEST*/
72