Actual source code: ex139.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: const char help[] = "Test MatCreateLocalRef()\n\n";

  4: #include <petscmat.h>

  8: static PetscErrorCode GetLocalRef(Mat A,IS isrow,IS iscol,Mat *B)
  9: {
 11:   IS             istmp;

 14:   PetscPrintf(PETSC_COMM_WORLD,"Extracting LocalRef with isrow:\n");
 15:   ISOnComm(isrow,PETSC_COMM_WORLD,PETSC_COPY_VALUES,&istmp);
 16:   ISView(istmp,PETSC_VIEWER_STDOUT_WORLD);
 17:   ISDestroy(&istmp);
 18:   PetscPrintf(PETSC_COMM_WORLD,"Extracting LocalRef with iscol (only rank=0 shown):\n");
 19:   ISOnComm(iscol,PETSC_COMM_WORLD,PETSC_COPY_VALUES,&istmp);
 20:   ISView(istmp,PETSC_VIEWER_STDOUT_WORLD);
 21:   ISDestroy(&istmp);

 23:   MatCreateLocalRef(A,isrow,iscol,B);
 24:   return(0);
 25: }

 29: int main(int argc,char *argv[])
 30: {
 31:   PetscErrorCode         ierr;
 32:   MPI_Comm               comm;
 33:   Mat                    J,B;
 34:   PetscInt               i,j,k,l,rstart,rend,m,n,top_bs,row_bs,col_bs,nlocblocks,*idx,nrowblocks,ncolblocks,*ridx,*cidx,*icol,*irow;
 35:   PetscScalar            *vals;
 36:   ISLocalToGlobalMapping brmap;
 37:   IS                     is0,is1;
 38:   PetscBool              diag,blocked;

 40:   PetscInitialize(&argc,&argv,0,help);
 41:   comm = PETSC_COMM_WORLD;

 43:   PetscOptionsBegin(comm,NULL,"LocalRef Test Options",NULL);
 44:   {
 45:     top_bs = 2; row_bs = 2; col_bs = 2; diag = PETSC_FALSE; blocked = PETSC_FALSE;
 46:     PetscOptionsInt("-top_bs","Block size of top-level matrix",0,top_bs,&top_bs,NULL);
 47:     PetscOptionsInt("-row_bs","Block size of row map",0,row_bs,&row_bs,NULL);
 48:     PetscOptionsInt("-col_bs","Block size of col map",0,col_bs,&col_bs,NULL);
 49:     PetscOptionsBool("-diag","Extract a diagonal black",0,diag,&diag,NULL);
 50:     PetscOptionsBool("-blocked","Use block insertion",0,blocked,&blocked,NULL);
 51:   }
 52:   PetscOptionsEnd();

 54:   MatCreate(comm,&J);
 55:   MatSetSizes(J,6,6,PETSC_DETERMINE,PETSC_DETERMINE);
 56:   MatSetBlockSize(J,top_bs);
 57:   MatSetFromOptions(J);
 58:   MatSeqBAIJSetPreallocation(J,top_bs,PETSC_DECIDE,0);
 59:   MatMPIBAIJSetPreallocation(J,top_bs,PETSC_DECIDE,0,PETSC_DECIDE,0);
 60:   MatSetUp(J);
 61:   MatGetSize(J,&m,&n);
 62:   MatGetOwnershipRange(J,&rstart,&rend);

 64:   nlocblocks = (rend-rstart)/top_bs + 2;
 65:   PetscMalloc1(nlocblocks,&idx);
 66:   for (i=0; i<nlocblocks; i++) {
 67:     idx[i] = (rstart/top_bs + i - 1 + m/top_bs) % (m/top_bs);
 68:   }
 69:   ISLocalToGlobalMappingCreate(comm,top_bs,nlocblocks,idx,PETSC_OWN_POINTER,&brmap);
 70:   PetscPrintf(comm,"Block ISLocalToGlobalMapping:\n");
 71:   ISLocalToGlobalMappingView(brmap,PETSC_VIEWER_STDOUT_WORLD);

 73:   MatSetLocalToGlobalMapping(J,brmap,brmap);
 74:   ISLocalToGlobalMappingDestroy(&brmap);

 76:   /* Create index sets for local submatrix */
 77:   nrowblocks = (rend-rstart)/row_bs;
 78:   ncolblocks = (rend-rstart)/col_bs;
 79:   PetscMalloc2(nrowblocks,&ridx,ncolblocks,&cidx);
 80:   for (i=0; i<nrowblocks; i++) ridx[i] = i + ((i > nrowblocks/2) ^ !rstart);
 81:   for (i=0; i<ncolblocks; i++) cidx[i] = i + ((i < ncolblocks/2) ^ !!rstart);
 82:   ISCreateBlock(PETSC_COMM_SELF,row_bs,nrowblocks,ridx,PETSC_COPY_VALUES,&is0);
 83:   ISCreateBlock(PETSC_COMM_SELF,col_bs,ncolblocks,cidx,PETSC_COPY_VALUES,&is1);
 84:   PetscFree2(ridx,cidx);

 86:   if (diag) {
 87:     ISDestroy(&is1);
 88:     PetscObjectReference((PetscObject)is0);
 89:     is1        = is0;
 90:     ncolblocks = nrowblocks;
 91:   }

 93:   GetLocalRef(J,is0,is1,&B);

 95:   MatZeroEntries(J);

 97:   PetscMalloc3(row_bs,&irow,col_bs,&icol,row_bs*col_bs,&vals);
 98:   for (i=0; i<nrowblocks; i++) {
 99:     for (j=0; j<ncolblocks; j++) {
100:       for (k=0; k<row_bs; k++) {
101:         for (l=0; l<col_bs; l++) {
102:           vals[k*col_bs+l] = i * 100000 + j * 1000 + k * 10 + l;
103:         }
104:         irow[k] = i*row_bs+k;
105:       }
106:       for (l=0; l<col_bs; l++) icol[l] = j*col_bs+l;
107:       if (blocked) {
108:         MatSetValuesBlockedLocal(B,1,&i,1,&j,vals,ADD_VALUES);
109:       } else {
110:         MatSetValuesLocal(B,row_bs,irow,col_bs,icol,vals,ADD_VALUES);
111:       }
112:     }
113:   }
114:   PetscFree3(irow,icol,vals);

116:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
117:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);

119:   MatView(J,PETSC_VIEWER_STDOUT_WORLD);

121:   ISDestroy(&is0);
122:   ISDestroy(&is1);
123:   MatDestroy(&B);
124:   MatDestroy(&J);
125:   PetscFinalize();
126:   return 0;
127: }