xref: /petsc/src/mat/tests/ex139.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown const char help[] = "Test MatCreateLocalRef()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown static PetscErrorCode GetLocalRef(Mat A,IS isrow,IS iscol,Mat *B)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   IS             istmp;
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   PetscFunctionBegin;
115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Extracting LocalRef with isrow:\n"));
125f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(isrow,PETSC_COMM_WORLD,PETSC_COPY_VALUES,&istmp));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(istmp,PETSC_VIEWER_STDOUT_WORLD));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&istmp));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Extracting LocalRef with iscol (only rank=0 shown):\n"));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(ISOnComm(iscol,PETSC_COMM_WORLD,PETSC_COPY_VALUES,&istmp));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(istmp,PETSC_VIEWER_STDOUT_WORLD));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&istmp));
19c4762a1bSJed Brown 
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateLocalRef(A,isrow,iscol,B));
21c4762a1bSJed Brown   PetscFunctionReturn(0);
22c4762a1bSJed Brown }
23c4762a1bSJed Brown 
24c4762a1bSJed Brown int main(int argc,char *argv[])
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   PetscErrorCode         ierr;
27c4762a1bSJed Brown   MPI_Comm               comm;
28c4762a1bSJed Brown   Mat                    J,B;
29c4762a1bSJed Brown   PetscInt               i,j,k,l,rstart,rend,m,n,top_bs,row_bs,col_bs,nlocblocks,*idx,nrowblocks,ncolblocks,*ridx,*cidx,*icol,*irow;
30c4762a1bSJed Brown   PetscScalar            *vals;
31c4762a1bSJed Brown   ISLocalToGlobalMapping brmap;
32c4762a1bSJed Brown   IS                     is0,is1;
33c4762a1bSJed Brown   PetscBool              diag,blocked;
34c4762a1bSJed Brown 
35*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
36c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"LocalRef Test Options",NULL);CHKERRQ(ierr);
39c4762a1bSJed Brown   {
40c4762a1bSJed Brown     top_bs = 2; row_bs = 2; col_bs = 2; diag = PETSC_FALSE; blocked = PETSC_FALSE;
415f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-top_bs","Block size of top-level matrix",0,top_bs,&top_bs,NULL));
425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-row_bs","Block size of row map",0,row_bs,&row_bs,NULL));
435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsInt("-col_bs","Block size of col map",0,col_bs,&col_bs,NULL));
445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsBool("-diag","Extract a diagonal black",0,diag,&diag,NULL));
455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsBool("-blocked","Use block insertion",0,blocked,&blocked,NULL));
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(comm,&J));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,6,6,PETSC_DETERMINE,PETSC_DETERMINE));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(J,top_bs));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqBAIJSetPreallocation(J,top_bs,PETSC_DECIDE,0));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIBAIJSetPreallocation(J,top_bs,PETSC_DECIDE,0,PETSC_DECIDE,0));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(J));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(J,&m,&n));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(J,&rstart,&rend));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   nlocblocks = (rend-rstart)/top_bs + 2;
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nlocblocks,&idx));
61c4762a1bSJed Brown   for (i=0; i<nlocblocks; i++) {
62c4762a1bSJed Brown     idx[i] = (rstart/top_bs + i - 1 + m/top_bs) % (m/top_bs);
63c4762a1bSJed Brown   }
645f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(comm,top_bs,nlocblocks,idx,PETSC_OWN_POINTER,&brmap));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"Block ISLocalToGlobalMapping:\n"));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingView(brmap,PETSC_VIEWER_STDOUT_WORLD));
67c4762a1bSJed Brown 
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetLocalToGlobalMapping(J,brmap,brmap));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&brmap));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Create index sets for local submatrix */
72c4762a1bSJed Brown   nrowblocks = (rend-rstart)/row_bs;
73c4762a1bSJed Brown   ncolblocks = (rend-rstart)/col_bs;
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nrowblocks,&ridx,ncolblocks,&cidx));
75c4762a1bSJed Brown   for (i=0; i<nrowblocks; i++) ridx[i] = i + ((i > nrowblocks/2) ^ !rstart);
76c4762a1bSJed Brown   for (i=0; i<ncolblocks; i++) cidx[i] = i + ((i < ncolblocks/2) ^ !!rstart);
775f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,row_bs,nrowblocks,ridx,PETSC_COPY_VALUES,&is0));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateBlock(PETSC_COMM_SELF,col_bs,ncolblocks,cidx,PETSC_COPY_VALUES,&is1));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(ridx,cidx));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   if (diag) {
825f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectReference((PetscObject)is0));
84c4762a1bSJed Brown     is1        = is0;
85c4762a1bSJed Brown     ncolblocks = nrowblocks;
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown 
885f80ce2aSJacob Faibussowitsch   CHKERRQ(GetLocalRef(J,is0,is1,&B));
89c4762a1bSJed Brown 
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(J));
91c4762a1bSJed Brown 
925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(row_bs,&irow,col_bs,&icol,row_bs*col_bs,&vals));
93c4762a1bSJed Brown   for (i=0; i<nrowblocks; i++) {
94c4762a1bSJed Brown     for (j=0; j<ncolblocks; j++) {
95c4762a1bSJed Brown       for (k=0; k<row_bs; k++) {
96c4762a1bSJed Brown         for (l=0; l<col_bs; l++) {
97c4762a1bSJed Brown           vals[k*col_bs+l] = i * 100000 + j * 1000 + k * 10 + l;
98c4762a1bSJed Brown         }
99c4762a1bSJed Brown         irow[k] = i*row_bs+k;
100c4762a1bSJed Brown       }
101c4762a1bSJed Brown       for (l=0; l<col_bs; l++) icol[l] = j*col_bs+l;
102c4762a1bSJed Brown       if (blocked) {
1035f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesBlockedLocal(B,1,&i,1,&j,vals,ADD_VALUES));
104c4762a1bSJed Brown       } else {
1055f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesLocal(B,row_bs,irow,col_bs,icol,vals,ADD_VALUES));
106c4762a1bSJed Brown       }
107c4762a1bSJed Brown     }
108c4762a1bSJed Brown   }
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(irow,icol,vals));
110c4762a1bSJed Brown 
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
113c4762a1bSJed Brown 
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(J,PETSC_VIEWER_STDOUT_WORLD));
115c4762a1bSJed Brown 
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is0));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
120*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
121*b122ec5aSJacob Faibussowitsch   return 0;
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown /*TEST
125c4762a1bSJed Brown 
126c4762a1bSJed Brown    test:
127c4762a1bSJed Brown       nsize: 2
128c4762a1bSJed Brown       filter: grep -v "type: mpi"
129c4762a1bSJed Brown       args: -blocked {{0 1}} -mat_type {{aij baij}}
130c4762a1bSJed Brown 
131c4762a1bSJed Brown TEST*/
132