xref: /petsc/src/mat/tests/ex212.c (revision 73b6653f4f4dddc0b51913ad0a8248d4f409633b)
1*73b6653fSLawrence Mitchell static char help[] = "Test MatCreateSubMatrix with -mat_type nest and block sizes.\n";
2*73b6653fSLawrence Mitchell 
3*73b6653fSLawrence Mitchell #include <petscmat.h>
4*73b6653fSLawrence Mitchell 
5*73b6653fSLawrence Mitchell int main(int argc, char **argv)
6*73b6653fSLawrence Mitchell {
7*73b6653fSLawrence Mitchell   Mat                    A, B, C, mats[6];
8*73b6653fSLawrence Mitchell   IS                     rows[2];
9*73b6653fSLawrence Mitchell   ISLocalToGlobalMapping cmap, rmap;
10*73b6653fSLawrence Mitchell   const PetscInt         indices[3] = {0, 1, 2};
11*73b6653fSLawrence Mitchell   PetscInt               i;
12*73b6653fSLawrence Mitchell   PetscMPIInt            size;
13*73b6653fSLawrence Mitchell   PetscErrorCode         ierr;
14*73b6653fSLawrence Mitchell 
15*73b6653fSLawrence Mitchell   ierr = PetscInitialize(&argc,&argv,NULL,help); if (ierr) return ierr;
16*73b6653fSLawrence Mitchell   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
17*73b6653fSLawrence Mitchell   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Only coded for one process");
18*73b6653fSLawrence Mitchell   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
19*73b6653fSLawrence Mitchell   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
20*73b6653fSLawrence Mitchell   ierr = MatSetBlockSizes(A,2,1);CHKERRQ(ierr);
21*73b6653fSLawrence Mitchell   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
22*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,2,1,indices,PETSC_COPY_VALUES,&rmap);CHKERRQ(ierr);
23*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,1,indices,PETSC_COPY_VALUES,&cmap);CHKERRQ(ierr);
24*73b6653fSLawrence Mitchell   ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);
25*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
26*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
27*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,2,indices,PETSC_COPY_VALUES,&rmap);CHKERRQ(ierr);
28*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,indices,PETSC_COPY_VALUES,&cmap);CHKERRQ(ierr);
29*73b6653fSLawrence Mitchell   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
30*73b6653fSLawrence Mitchell   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,2,3);CHKERRQ(ierr);
31*73b6653fSLawrence Mitchell   ierr = MatSetBlockSizes(B,1,1);CHKERRQ(ierr);
32*73b6653fSLawrence Mitchell   ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
33*73b6653fSLawrence Mitchell   ierr = MatSetLocalToGlobalMapping(B,rmap,cmap);CHKERRQ(ierr);
34*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
35*73b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
36*73b6653fSLawrence Mitchell   ierr = MatSetUp(A);CHKERRQ(ierr);
37*73b6653fSLawrence Mitchell   ierr = MatSetUp(B);CHKERRQ(ierr);
38*73b6653fSLawrence Mitchell   mats[0] = A;
39*73b6653fSLawrence Mitchell   mats[1] = B;
40*73b6653fSLawrence Mitchell   mats[2] = A;
41*73b6653fSLawrence Mitchell   mats[3] = NULL;
42*73b6653fSLawrence Mitchell   mats[4] = B;
43*73b6653fSLawrence Mitchell   mats[5] = A;
44*73b6653fSLawrence Mitchell   ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,3,NULL,mats,&C);CHKERRQ(ierr);
45*73b6653fSLawrence Mitchell   ierr = MatSetUp(C);CHKERRQ(ierr);
46*73b6653fSLawrence Mitchell   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47*73b6653fSLawrence Mitchell   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48*73b6653fSLawrence Mitchell   ierr = MatView(C,NULL);CHKERRQ(ierr);
49*73b6653fSLawrence Mitchell   ierr = MatNestGetISs(C,rows,NULL);CHKERRQ(ierr);
50*73b6653fSLawrence Mitchell   for (i=0; i<2; i++) {
51*73b6653fSLawrence Mitchell     Mat      submat;
52*73b6653fSLawrence Mitchell     IS       cols[3];
53*73b6653fSLawrence Mitchell     PetscInt j;
54*73b6653fSLawrence Mitchell     ierr = ISView(rows[i],NULL);CHKERRQ(ierr);
55*73b6653fSLawrence Mitchell     ierr = MatCreateSubMatrix(C,rows[i],NULL,MAT_INITIAL_MATRIX,&submat);CHKERRQ(ierr);
56*73b6653fSLawrence Mitchell     ierr = MatView(submat,NULL);CHKERRQ(ierr);
57*73b6653fSLawrence Mitchell     ierr = MatNestGetISs(submat,NULL,cols);CHKERRQ(ierr);
58*73b6653fSLawrence Mitchell     for (j=0; j<3; j++) {
59*73b6653fSLawrence Mitchell       ierr = ISView(cols[j],NULL);CHKERRQ(ierr);
60*73b6653fSLawrence Mitchell     }
61*73b6653fSLawrence Mitchell     ierr = MatDestroy(&submat);CHKERRQ(ierr);
62*73b6653fSLawrence Mitchell   }
63*73b6653fSLawrence Mitchell   ierr = MatDestroy(&C);CHKERRQ(ierr);
64*73b6653fSLawrence Mitchell   ierr = MatDestroy(&B);CHKERRQ(ierr);
65*73b6653fSLawrence Mitchell   ierr = MatDestroy(&A);CHKERRQ(ierr);
66*73b6653fSLawrence Mitchell   ierr = PetscFinalize();
67*73b6653fSLawrence Mitchell   return ierr;
68*73b6653fSLawrence Mitchell }
69*73b6653fSLawrence Mitchell 
70*73b6653fSLawrence Mitchell /*TEST
71*73b6653fSLawrence Mitchell 
72*73b6653fSLawrence Mitchell     test:
73*73b6653fSLawrence Mitchell 
74*73b6653fSLawrence Mitchell TEST*/
75