xref: /petsc/src/mat/tests/ex139.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
1 
2 const char help[] = "Test MatCreateLocalRef()\n\n";
3 
4 #include <petscmat.h>
5 
6 static PetscErrorCode GetLocalRef(Mat A,IS isrow,IS iscol,Mat *B)
7 {
8   IS             istmp;
9 
10   PetscFunctionBegin;
11   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Extracting LocalRef with isrow:\n"));
12   PetscCall(ISOnComm(isrow,PETSC_COMM_WORLD,PETSC_COPY_VALUES,&istmp));
13   PetscCall(ISView(istmp,PETSC_VIEWER_STDOUT_WORLD));
14   PetscCall(ISDestroy(&istmp));
15   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Extracting LocalRef with iscol (only rank=0 shown):\n"));
16   PetscCall(ISOnComm(iscol,PETSC_COMM_WORLD,PETSC_COPY_VALUES,&istmp));
17   PetscCall(ISView(istmp,PETSC_VIEWER_STDOUT_WORLD));
18   PetscCall(ISDestroy(&istmp));
19 
20   PetscCall(MatCreateLocalRef(A,isrow,iscol,B));
21   PetscFunctionReturn(0);
22 }
23 
24 int main(int argc,char *argv[])
25 {
26   MPI_Comm               comm;
27   Mat                    J,B;
28   PetscInt               i,j,k,l,rstart,rend,m,n,top_bs,row_bs,col_bs,nlocblocks,*idx,nrowblocks,ncolblocks,*ridx,*cidx,*icol,*irow;
29   PetscScalar            *vals;
30   ISLocalToGlobalMapping brmap;
31   IS                     is0,is1;
32   PetscBool              diag,blocked;
33 
34   PetscCall(PetscInitialize(&argc,&argv,0,help));
35   comm = PETSC_COMM_WORLD;
36 
37   PetscOptionsBegin(comm,NULL,"LocalRef Test Options",NULL);
38   {
39     top_bs = 2; row_bs = 2; col_bs = 2; diag = PETSC_FALSE; blocked = PETSC_FALSE;
40     PetscCall(PetscOptionsInt("-top_bs","Block size of top-level matrix",0,top_bs,&top_bs,NULL));
41     PetscCall(PetscOptionsInt("-row_bs","Block size of row map",0,row_bs,&row_bs,NULL));
42     PetscCall(PetscOptionsInt("-col_bs","Block size of col map",0,col_bs,&col_bs,NULL));
43     PetscCall(PetscOptionsBool("-diag","Extract a diagonal black",0,diag,&diag,NULL));
44     PetscCall(PetscOptionsBool("-blocked","Use block insertion",0,blocked,&blocked,NULL));
45   }
46   PetscOptionsEnd();
47 
48   PetscCall(MatCreate(comm,&J));
49   PetscCall(MatSetSizes(J,6,6,PETSC_DETERMINE,PETSC_DETERMINE));
50   PetscCall(MatSetBlockSize(J,top_bs));
51   PetscCall(MatSetFromOptions(J));
52   PetscCall(MatSeqBAIJSetPreallocation(J,top_bs,PETSC_DECIDE,0));
53   PetscCall(MatMPIBAIJSetPreallocation(J,top_bs,PETSC_DECIDE,0,PETSC_DECIDE,0));
54   PetscCall(MatSetUp(J));
55   PetscCall(MatGetSize(J,&m,&n));
56   PetscCall(MatGetOwnershipRange(J,&rstart,&rend));
57 
58   nlocblocks = (rend-rstart)/top_bs + 2;
59   PetscCall(PetscMalloc1(nlocblocks,&idx));
60   for (i=0; i<nlocblocks; i++) {
61     idx[i] = (rstart/top_bs + i - 1 + m/top_bs) % (m/top_bs);
62   }
63   PetscCall(ISLocalToGlobalMappingCreate(comm,top_bs,nlocblocks,idx,PETSC_OWN_POINTER,&brmap));
64   PetscCall(PetscPrintf(comm,"Block ISLocalToGlobalMapping:\n"));
65   PetscCall(ISLocalToGlobalMappingView(brmap,PETSC_VIEWER_STDOUT_WORLD));
66 
67   PetscCall(MatSetLocalToGlobalMapping(J,brmap,brmap));
68   PetscCall(ISLocalToGlobalMappingDestroy(&brmap));
69 
70   /* Create index sets for local submatrix */
71   nrowblocks = (rend-rstart)/row_bs;
72   ncolblocks = (rend-rstart)/col_bs;
73   PetscCall(PetscMalloc2(nrowblocks,&ridx,ncolblocks,&cidx));
74   for (i=0; i<nrowblocks; i++) ridx[i] = i + ((i > nrowblocks/2) ^ !rstart);
75   for (i=0; i<ncolblocks; i++) cidx[i] = i + ((i < ncolblocks/2) ^ !!rstart);
76   PetscCall(ISCreateBlock(PETSC_COMM_SELF,row_bs,nrowblocks,ridx,PETSC_COPY_VALUES,&is0));
77   PetscCall(ISCreateBlock(PETSC_COMM_SELF,col_bs,ncolblocks,cidx,PETSC_COPY_VALUES,&is1));
78   PetscCall(PetscFree2(ridx,cidx));
79 
80   if (diag) {
81     PetscCall(ISDestroy(&is1));
82     PetscCall(PetscObjectReference((PetscObject)is0));
83     is1        = is0;
84     ncolblocks = nrowblocks;
85   }
86 
87   PetscCall(GetLocalRef(J,is0,is1,&B));
88 
89   PetscCall(MatZeroEntries(J));
90 
91   PetscCall(PetscMalloc3(row_bs,&irow,col_bs,&icol,row_bs*col_bs,&vals));
92   for (i=0; i<nrowblocks; i++) {
93     for (j=0; j<ncolblocks; j++) {
94       for (k=0; k<row_bs; k++) {
95         for (l=0; l<col_bs; l++) {
96           vals[k*col_bs+l] = i * 100000 + j * 1000 + k * 10 + l;
97         }
98         irow[k] = i*row_bs+k;
99       }
100       for (l=0; l<col_bs; l++) icol[l] = j*col_bs+l;
101       if (blocked) {
102         PetscCall(MatSetValuesBlockedLocal(B,1,&i,1,&j,vals,ADD_VALUES));
103       } else {
104         PetscCall(MatSetValuesLocal(B,row_bs,irow,col_bs,icol,vals,ADD_VALUES));
105       }
106     }
107   }
108   PetscCall(PetscFree3(irow,icol,vals));
109 
110   PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
111   PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
112 
113   PetscCall(MatView(J,PETSC_VIEWER_STDOUT_WORLD));
114 
115   PetscCall(ISDestroy(&is0));
116   PetscCall(ISDestroy(&is1));
117   PetscCall(MatDestroy(&B));
118   PetscCall(MatDestroy(&J));
119   PetscCall(PetscFinalize());
120   return 0;
121 }
122 
123 /*TEST
124 
125    test:
126       nsize: 2
127       filter: grep -v "type: mpi"
128       args: -blocked {{0 1}} -mat_type {{aij baij}}
129 
130 TEST*/
131