1*f1af5d2fSBarry Smith /*$Id: mpiadj.c,v 1.31 1999/10/24 14:02:36 bsmith Exp bsmith $*/ 2b97cf49bSBarry Smith 3b97cf49bSBarry Smith /* 4b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 5b97cf49bSBarry Smith */ 6b97cf49bSBarry Smith #include "sys.h" 70752156aSBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h" 8b97cf49bSBarry Smith 9b97cf49bSBarry Smith #undef __FUNC__ 100752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj_ASCII" 110752156aSBarry Smith extern int MatView_MPIAdj_ASCII(Mat A,Viewer viewer) 12b97cf49bSBarry Smith { 130752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 14b97cf49bSBarry Smith int ierr, i,j, m = a->m, format; 15b97cf49bSBarry Smith char *outputname; 16b97cf49bSBarry Smith 17433994e6SBarry Smith PetscFunctionBegin; 1877ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 19d132466eSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 20b97cf49bSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 213a40ed3dSBarry Smith PetscFunctionReturn(0); 220752156aSBarry Smith } else { 236831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 24b97cf49bSBarry Smith for ( i=0; i<m; i++ ) { 256831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr); 26b97cf49bSBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 276831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr); 280752156aSBarry Smith } 296831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 30b97cf49bSBarry Smith } 316831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 32b97cf49bSBarry Smith } 336831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 343a40ed3dSBarry Smith PetscFunctionReturn(0); 35b97cf49bSBarry Smith } 36b97cf49bSBarry Smith 37b97cf49bSBarry Smith #undef __FUNC__ 380752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj" 39e1311b90SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer) 40b97cf49bSBarry Smith { 41b97cf49bSBarry Smith int ierr; 426831982aSBarry Smith PetscTruth isascii; 43b97cf49bSBarry Smith 44433994e6SBarry Smith PetscFunctionBegin; 456831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 460f5bd95cSBarry Smith if (isascii) { 473a40ed3dSBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 485cd90555SBarry Smith } else { 490f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 50b97cf49bSBarry Smith } 513a40ed3dSBarry Smith PetscFunctionReturn(0); 52b97cf49bSBarry Smith } 53b97cf49bSBarry Smith 54b97cf49bSBarry Smith #undef __FUNC__ 550752156aSBarry Smith #define __FUNC__ "MatDestroy_MPIAdj" 5694d884c6SBarry Smith int MatDestroy_MPIAdj(Mat mat) 57b97cf49bSBarry Smith { 5894d884c6SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) mat->data; 5994d884c6SBarry Smith int ierr; 6094d884c6SBarry Smith 6194d884c6SBarry Smith PetscFunctionBegin; 6294d884c6SBarry Smith 6394d884c6SBarry Smith if (mat->mapping) { 6494d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 6594d884c6SBarry Smith } 6694d884c6SBarry Smith if (mat->bmapping) { 6794d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 6894d884c6SBarry Smith } 69c7fcc2eaSBarry Smith if (mat->rmap) { 70c7fcc2eaSBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 71c7fcc2eaSBarry Smith } 72c7fcc2eaSBarry Smith if (mat->cmap) { 73c7fcc2eaSBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 74c7fcc2eaSBarry Smith } 75b97cf49bSBarry Smith 76aa482453SBarry Smith #if defined(PETSC_USE_LOG) 7794d884c6SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz); 78b97cf49bSBarry Smith #endif 79606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 80606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 81606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 82606d414cSSatish Balay ierr = PetscFree(a->rowners);CHKERRQ(ierr); 83606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 84b97cf49bSBarry Smith 8594d884c6SBarry Smith PLogObjectDestroy(mat); 8694d884c6SBarry Smith PetscHeaderDestroy(mat); 873a40ed3dSBarry Smith PetscFunctionReturn(0); 88b97cf49bSBarry Smith } 89b97cf49bSBarry Smith 90b97cf49bSBarry Smith 91b97cf49bSBarry Smith #undef __FUNC__ 920752156aSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj" 930752156aSBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 94b97cf49bSBarry Smith { 950752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 96b97cf49bSBarry Smith 97433994e6SBarry Smith PetscFunctionBegin; 98b97cf49bSBarry Smith if (op == MAT_STRUCTURALLY_SYMMETRIC) { 99b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 100b97cf49bSBarry Smith } else { 101981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 102b97cf49bSBarry Smith } 1033a40ed3dSBarry Smith PetscFunctionReturn(0); 104b97cf49bSBarry Smith } 105b97cf49bSBarry Smith 106b97cf49bSBarry Smith 107b97cf49bSBarry Smith /* 108b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 109b97cf49bSBarry Smith */ 110b97cf49bSBarry Smith 111b97cf49bSBarry Smith #undef __FUNC__ 1120752156aSBarry Smith #define __FUNC__ "MatMarkDiag_MPIAdj" 1130752156aSBarry Smith int MatMarkDiag_MPIAdj(Mat A) 114b97cf49bSBarry Smith { 1150752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 116b97cf49bSBarry Smith int i,j, *diag, m = a->m; 117b97cf49bSBarry Smith 118433994e6SBarry Smith PetscFunctionBegin; 119b97cf49bSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 120b97cf49bSBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 121b97cf49bSBarry Smith for ( i=0; i<a->m; i++ ) { 122b97cf49bSBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 123b97cf49bSBarry Smith if (a->j[j] == i) { 124b97cf49bSBarry Smith diag[i] = j; 125b97cf49bSBarry Smith break; 126b97cf49bSBarry Smith } 127b97cf49bSBarry Smith } 128b97cf49bSBarry Smith } 129b97cf49bSBarry Smith a->diag = diag; 1303a40ed3dSBarry Smith PetscFunctionReturn(0); 131b97cf49bSBarry Smith } 132b97cf49bSBarry Smith 133b97cf49bSBarry Smith #undef __FUNC__ 1340752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj" 1350752156aSBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n) 136b97cf49bSBarry Smith { 137433994e6SBarry Smith PetscFunctionBegin; 1380752156aSBarry Smith if (m) *m = A->M; 1390752156aSBarry Smith if (n) *n = A->N; 1403a40ed3dSBarry Smith PetscFunctionReturn(0); 141b97cf49bSBarry Smith } 142b97cf49bSBarry Smith 143b97cf49bSBarry Smith #undef __FUNC__ 1440752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj" 1450752156aSBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n) 146b97cf49bSBarry Smith { 1470752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 148433994e6SBarry Smith PetscFunctionBegin; 1490752156aSBarry Smith if (m) *m = a->m; 1500752156aSBarry Smith if (n) *n = A->N; 1513a40ed3dSBarry Smith PetscFunctionReturn(0); 152b97cf49bSBarry Smith } 153b97cf49bSBarry Smith 154b97cf49bSBarry Smith #undef __FUNC__ 1550752156aSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj" 1560752156aSBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n) 157b97cf49bSBarry Smith { 1580752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 159433994e6SBarry Smith PetscFunctionBegin; 1600752156aSBarry Smith *m = a->rstart; *n = a->rend; 1613a40ed3dSBarry Smith PetscFunctionReturn(0); 162b97cf49bSBarry Smith } 163b97cf49bSBarry Smith 164b97cf49bSBarry Smith #undef __FUNC__ 1650752156aSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj" 1660752156aSBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 167b97cf49bSBarry Smith { 1680752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 169b97cf49bSBarry Smith int *itmp; 170b97cf49bSBarry Smith 171433994e6SBarry Smith PetscFunctionBegin; 1720752156aSBarry Smith row -= a->rstart; 1730752156aSBarry Smith 174a8c6a408SBarry Smith if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 175b97cf49bSBarry Smith 176b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 177b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 178b97cf49bSBarry Smith if (idx) { 179b97cf49bSBarry Smith itmp = a->j + a->i[row]; 180b97cf49bSBarry Smith if (*nz) { 181b97cf49bSBarry Smith *idx = itmp; 182b97cf49bSBarry Smith } 183b97cf49bSBarry Smith else *idx = 0; 184b97cf49bSBarry Smith } 1853a40ed3dSBarry Smith PetscFunctionReturn(0); 186b97cf49bSBarry Smith } 187b97cf49bSBarry Smith 188b97cf49bSBarry Smith #undef __FUNC__ 1890752156aSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj" 1900752156aSBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 191b97cf49bSBarry Smith { 192433994e6SBarry Smith PetscFunctionBegin; 1933a40ed3dSBarry Smith PetscFunctionReturn(0); 194b97cf49bSBarry Smith } 195b97cf49bSBarry Smith 196b97cf49bSBarry Smith #undef __FUNC__ 1970752156aSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj" 1980752156aSBarry Smith int MatGetBlockSize_MPIAdj(Mat A, int *bs) 199b97cf49bSBarry Smith { 200433994e6SBarry Smith PetscFunctionBegin; 201b97cf49bSBarry Smith *bs = 1; 2023a40ed3dSBarry Smith PetscFunctionReturn(0); 203b97cf49bSBarry Smith } 204b97cf49bSBarry Smith 205b97cf49bSBarry Smith 206b97cf49bSBarry Smith #undef __FUNC__ 2070752156aSBarry Smith #define __FUNC__ "MatEqual_MPIAdj" 2080752156aSBarry Smith int MatEqual_MPIAdj(Mat A,Mat B, PetscTruth* flg) 209b97cf49bSBarry Smith { 2100752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data; 2110f5bd95cSBarry Smith int ierr; 2120f5bd95cSBarry Smith PetscTruth flag; 213b97cf49bSBarry Smith 214433994e6SBarry Smith PetscFunctionBegin; 215a8c6a408SBarry Smith if (B->type != MATMPIADJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 216b97cf49bSBarry Smith 217b97cf49bSBarry Smith /* If the matrix dimensions are not equal, or no of nonzeros */ 2180752156aSBarry Smith if ((a->m != b->m ) ||( a->nz != b->nz)) { 2190f5bd95cSBarry Smith flag = PETSC_FALSE; 220b97cf49bSBarry Smith } 221b97cf49bSBarry Smith 222b97cf49bSBarry Smith /* if the a->i are the same */ 2230f5bd95cSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int),&flag);CHKERRQ(ierr); 224b97cf49bSBarry Smith 225b97cf49bSBarry Smith /* if a->j are the same */ 2260f5bd95cSBarry Smith ierr = PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int),&flag);CHKERRQ(ierr); 227b97cf49bSBarry Smith 228ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 2293a40ed3dSBarry Smith PetscFunctionReturn(0); 230b97cf49bSBarry Smith } 231b97cf49bSBarry Smith 232b97cf49bSBarry Smith 233b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2340b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2350b3b1487SBarry Smith MatGetRow_MPIAdj, 2360b3b1487SBarry Smith MatRestoreRow_MPIAdj, 237b97cf49bSBarry Smith 0, 238b97cf49bSBarry Smith 0, 239b97cf49bSBarry Smith 0, 240b97cf49bSBarry Smith 0, 2410b3b1487SBarry Smith 0, 2420b3b1487SBarry Smith 0, 2430b3b1487SBarry Smith 0, 2440b3b1487SBarry Smith 0, 2450b3b1487SBarry Smith 0, 2460b3b1487SBarry Smith 0, 2470b3b1487SBarry Smith 0, 2480b3b1487SBarry Smith 0, 2490b3b1487SBarry Smith 0, 2500b3b1487SBarry Smith MatEqual_MPIAdj, 2510b3b1487SBarry Smith 0, 2520b3b1487SBarry Smith 0, 2530b3b1487SBarry Smith 0, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith MatSetOption_MPIAdj, 2580b3b1487SBarry Smith 0, 2590b3b1487SBarry Smith 0, 2600b3b1487SBarry Smith 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith MatGetSize_MPIAdj, 2650b3b1487SBarry Smith MatGetLocalSize_MPIAdj, 2660b3b1487SBarry Smith MatGetOwnershipRange_MPIAdj, 2670b3b1487SBarry Smith 0, 2680b3b1487SBarry Smith 0, 2690b3b1487SBarry Smith 0, 2700b3b1487SBarry Smith 0, 2710b3b1487SBarry Smith 0, 2720b3b1487SBarry Smith 0, 2730b3b1487SBarry Smith 0, 2740b3b1487SBarry Smith 0, 2750b3b1487SBarry Smith 0, 2760b3b1487SBarry Smith 0, 2770b3b1487SBarry Smith 0, 2780b3b1487SBarry Smith 0, 2790b3b1487SBarry Smith 0, 2800b3b1487SBarry Smith 0, 2810b3b1487SBarry Smith 0, 2820b3b1487SBarry Smith 0, 2830b3b1487SBarry Smith 0, 2840b3b1487SBarry Smith 0, 285b97cf49bSBarry Smith 0, 2860752156aSBarry Smith MatGetBlockSize_MPIAdj, 287b97cf49bSBarry Smith 0, 288b97cf49bSBarry Smith 0, 289b97cf49bSBarry Smith 0, 290b97cf49bSBarry Smith 0, 291b97cf49bSBarry Smith 0, 2920752156aSBarry Smith 0, 2930752156aSBarry Smith 0, 2940b3b1487SBarry Smith 0, 2950b3b1487SBarry Smith 0, 2960b3b1487SBarry Smith 0, 2970b3b1487SBarry Smith 0, 2980b3b1487SBarry Smith 0, 2990b3b1487SBarry Smith MatGetMaps_Petsc}; 300b97cf49bSBarry Smith 301b97cf49bSBarry Smith 302b97cf49bSBarry Smith #undef __FUNC__ 3030752156aSBarry Smith #define __FUNC__ "MatCreateMPIAdj" 304b97cf49bSBarry Smith /*@C 3050752156aSBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 306b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 307b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 308b97cf49bSBarry Smith 309ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 310ef5ee4d1SLois Curfman McInnes 311b97cf49bSBarry Smith Input Parameters: 312ef5ee4d1SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 3130752156aSBarry Smith . m - number of local rows 314b97cf49bSBarry Smith . n - number of columns 315b97cf49bSBarry Smith . i - the indices into j for the start of each row 316ef5ee4d1SLois Curfman McInnes - j - the column indices for each row (sorted for each row). 317ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 318b97cf49bSBarry Smith 319b97cf49bSBarry Smith Output Parameter: 320b97cf49bSBarry Smith . A - the matrix 321b97cf49bSBarry Smith 32215091d37SBarry Smith Level: intermediate 32315091d37SBarry Smith 3244bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 3254bc6d8c0SBarry Smith MatSetValues(). 3264bc6d8c0SBarry Smith You must NOT free the ii and jj arrays yourself. PETSc will free them 327b97cf49bSBarry Smith when the matrix is destroyed. 328b97cf49bSBarry Smith 329ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 330b97cf49bSBarry Smith 33191e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering() 332b97cf49bSBarry Smith @*/ 3330752156aSBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j, Mat *A) 334b97cf49bSBarry Smith { 335b97cf49bSBarry Smith Mat B; 3360752156aSBarry Smith Mat_MPIAdj *b; 337*f1af5d2fSBarry Smith int ii,ierr, size,rank; 338*f1af5d2fSBarry Smith PetscTruth flg; 339b97cf49bSBarry Smith 340433994e6SBarry Smith PetscFunctionBegin; 341d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 342d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 343b97cf49bSBarry Smith 344b97cf49bSBarry Smith *A = 0; 3453f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,"Mat",comm,MatDestroy,MatView); 346b97cf49bSBarry Smith PLogObjectCreate(B); 3470752156aSBarry Smith B->data = (void *) (b = PetscNew(Mat_MPIAdj));CHKPTRQ(b); 348549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 349549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 350e1311b90SBarry Smith B->ops->destroy = MatDestroy_MPIAdj; 351e1311b90SBarry Smith B->ops->view = MatView_MPIAdj; 352b97cf49bSBarry Smith B->factor = 0; 353b97cf49bSBarry Smith B->lupivotthreshold = 1.0; 354b97cf49bSBarry Smith B->mapping = 0; 355b97cf49bSBarry Smith B->assembled = PETSC_FALSE; 356b97cf49bSBarry Smith 3570752156aSBarry Smith b->m = m; B->m = m; 358ca161407SBarry Smith ierr = MPI_Allreduce(&m,&B->M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 3590752156aSBarry Smith B->n = n; B->N = n; 3600752156aSBarry Smith 361c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 362c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 363488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,B->M,&B->rmap);CHKERRQ(ierr); 364c7fcc2eaSBarry Smith /* we don't know the "local columns" so just use the row information :-( */ 365488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,B->M,&B->cmap);CHKERRQ(ierr); 366c7fcc2eaSBarry Smith 3670752156aSBarry Smith b->rowners = (int *) PetscMalloc((size+1)*sizeof(int));CHKPTRQ(b->rowners); 3680752156aSBarry Smith PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 369ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 3700752156aSBarry Smith b->rowners[0] = 0; 3710752156aSBarry Smith for ( ii=2; ii<=size; ii++ ) { 3720752156aSBarry Smith b->rowners[ii] += b->rowners[ii-1]; 3730752156aSBarry Smith } 3740752156aSBarry Smith b->rstart = b->rowners[rank]; 3750752156aSBarry Smith b->rend = b->rowners[rank+1]; 376b97cf49bSBarry Smith 377b97cf49bSBarry Smith b->j = j; 378b97cf49bSBarry Smith b->i = i; 379b97cf49bSBarry Smith 380b97cf49bSBarry Smith b->nz = i[m]; 381b97cf49bSBarry Smith b->diag = 0; 382b97cf49bSBarry Smith b->symmetric = PETSC_FALSE; 383b97cf49bSBarry Smith 384b97cf49bSBarry Smith *A = B; 385b97cf49bSBarry Smith 386b97cf49bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 387b97cf49bSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 388b97cf49bSBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 389b97cf49bSBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3903a40ed3dSBarry Smith PetscFunctionReturn(0); 391b97cf49bSBarry Smith } 392b97cf49bSBarry Smith 393b97cf49bSBarry Smith 394b97cf49bSBarry Smith 395433994e6SBarry Smith 396433994e6SBarry Smith 397433994e6SBarry Smith 398433994e6SBarry Smith 399433994e6SBarry Smith 400433994e6SBarry Smith 401