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