1b97cf49bSBarry Smith #ifdef PETSC_RCS_HEADER 2*6831982aSBarry Smith static char vcid[] = "$Id: mpiadj.c,v 1.29 1999/10/13 20:37:34 bsmith Exp bsmith $"; 3b97cf49bSBarry Smith #endif 4b97cf49bSBarry Smith 5b97cf49bSBarry Smith /* 6b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 7b97cf49bSBarry Smith */ 8b97cf49bSBarry Smith #include "sys.h" 90752156aSBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h" 10b97cf49bSBarry Smith 11b97cf49bSBarry Smith #undef __FUNC__ 120752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj_ASCII" 130752156aSBarry Smith extern int MatView_MPIAdj_ASCII(Mat A,Viewer viewer) 14b97cf49bSBarry Smith { 150752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 16b97cf49bSBarry Smith int ierr, i,j, m = a->m, format; 17b97cf49bSBarry Smith char *outputname; 18b97cf49bSBarry Smith 19433994e6SBarry Smith PetscFunctionBegin; 2077ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 21d132466eSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 22b97cf49bSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 233a40ed3dSBarry Smith PetscFunctionReturn(0); 240752156aSBarry Smith } else { 25*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 26b97cf49bSBarry Smith for ( i=0; i<m; i++ ) { 27*6831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr); 28b97cf49bSBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 29*6831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr); 300752156aSBarry Smith } 31*6831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 32b97cf49bSBarry Smith } 33*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 34b97cf49bSBarry Smith } 35*6831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 363a40ed3dSBarry Smith PetscFunctionReturn(0); 37b97cf49bSBarry Smith } 38b97cf49bSBarry Smith 39b97cf49bSBarry Smith #undef __FUNC__ 400752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj" 41e1311b90SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer) 42b97cf49bSBarry Smith { 43b97cf49bSBarry Smith int ierr; 44*6831982aSBarry Smith PetscTruth isascii; 45b97cf49bSBarry Smith 46433994e6SBarry Smith PetscFunctionBegin; 47*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 480f5bd95cSBarry Smith if (isascii) { 493a40ed3dSBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 505cd90555SBarry Smith } else { 510f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 52b97cf49bSBarry Smith } 533a40ed3dSBarry Smith PetscFunctionReturn(0); 54b97cf49bSBarry Smith } 55b97cf49bSBarry Smith 56b97cf49bSBarry Smith #undef __FUNC__ 570752156aSBarry Smith #define __FUNC__ "MatDestroy_MPIAdj" 5894d884c6SBarry Smith int MatDestroy_MPIAdj(Mat mat) 59b97cf49bSBarry Smith { 6094d884c6SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) mat->data; 6194d884c6SBarry Smith int ierr; 6294d884c6SBarry Smith 6394d884c6SBarry Smith PetscFunctionBegin; 6494d884c6SBarry Smith 6594d884c6SBarry Smith if (mat->mapping) { 6694d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 6794d884c6SBarry Smith } 6894d884c6SBarry Smith if (mat->bmapping) { 6994d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 7094d884c6SBarry Smith } 71c7fcc2eaSBarry Smith if (mat->rmap) { 72c7fcc2eaSBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 73c7fcc2eaSBarry Smith } 74c7fcc2eaSBarry Smith if (mat->cmap) { 75c7fcc2eaSBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 76c7fcc2eaSBarry Smith } 77b97cf49bSBarry Smith 78aa482453SBarry Smith #if defined(PETSC_USE_LOG) 7994d884c6SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz); 80b97cf49bSBarry Smith #endif 81606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 82606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 83606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 84606d414cSSatish Balay ierr = PetscFree(a->rowners);CHKERRQ(ierr); 85606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 86b97cf49bSBarry Smith 8794d884c6SBarry Smith PLogObjectDestroy(mat); 8894d884c6SBarry Smith PetscHeaderDestroy(mat); 893a40ed3dSBarry Smith PetscFunctionReturn(0); 90b97cf49bSBarry Smith } 91b97cf49bSBarry Smith 92b97cf49bSBarry Smith 93b97cf49bSBarry Smith #undef __FUNC__ 940752156aSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj" 950752156aSBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 96b97cf49bSBarry Smith { 970752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 98b97cf49bSBarry Smith 99433994e6SBarry Smith PetscFunctionBegin; 100b97cf49bSBarry Smith if (op == MAT_STRUCTURALLY_SYMMETRIC) { 101b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 102b97cf49bSBarry Smith } else { 103981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 104b97cf49bSBarry Smith } 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 106b97cf49bSBarry Smith } 107b97cf49bSBarry Smith 108b97cf49bSBarry Smith 109b97cf49bSBarry Smith /* 110b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 111b97cf49bSBarry Smith */ 112b97cf49bSBarry Smith 113b97cf49bSBarry Smith #undef __FUNC__ 1140752156aSBarry Smith #define __FUNC__ "MatMarkDiag_MPIAdj" 1150752156aSBarry Smith int MatMarkDiag_MPIAdj(Mat A) 116b97cf49bSBarry Smith { 1170752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 118b97cf49bSBarry Smith int i,j, *diag, m = a->m; 119b97cf49bSBarry Smith 120433994e6SBarry Smith PetscFunctionBegin; 121b97cf49bSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 122b97cf49bSBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 123b97cf49bSBarry Smith for ( i=0; i<a->m; i++ ) { 124b97cf49bSBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 125b97cf49bSBarry Smith if (a->j[j] == i) { 126b97cf49bSBarry Smith diag[i] = j; 127b97cf49bSBarry Smith break; 128b97cf49bSBarry Smith } 129b97cf49bSBarry Smith } 130b97cf49bSBarry Smith } 131b97cf49bSBarry Smith a->diag = diag; 1323a40ed3dSBarry Smith PetscFunctionReturn(0); 133b97cf49bSBarry Smith } 134b97cf49bSBarry Smith 135b97cf49bSBarry Smith #undef __FUNC__ 1360752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj" 1370752156aSBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n) 138b97cf49bSBarry Smith { 139433994e6SBarry Smith PetscFunctionBegin; 1400752156aSBarry Smith if (m) *m = A->M; 1410752156aSBarry Smith if (n) *n = A->N; 1423a40ed3dSBarry Smith PetscFunctionReturn(0); 143b97cf49bSBarry Smith } 144b97cf49bSBarry Smith 145b97cf49bSBarry Smith #undef __FUNC__ 1460752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj" 1470752156aSBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n) 148b97cf49bSBarry Smith { 1490752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 150433994e6SBarry Smith PetscFunctionBegin; 1510752156aSBarry Smith if (m) *m = a->m; 1520752156aSBarry Smith if (n) *n = A->N; 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 154b97cf49bSBarry Smith } 155b97cf49bSBarry Smith 156b97cf49bSBarry Smith #undef __FUNC__ 1570752156aSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj" 1580752156aSBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n) 159b97cf49bSBarry Smith { 1600752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 161433994e6SBarry Smith PetscFunctionBegin; 1620752156aSBarry Smith *m = a->rstart; *n = a->rend; 1633a40ed3dSBarry Smith PetscFunctionReturn(0); 164b97cf49bSBarry Smith } 165b97cf49bSBarry Smith 166b97cf49bSBarry Smith #undef __FUNC__ 1670752156aSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj" 1680752156aSBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 169b97cf49bSBarry Smith { 1700752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 171b97cf49bSBarry Smith int *itmp; 172b97cf49bSBarry Smith 173433994e6SBarry Smith PetscFunctionBegin; 1740752156aSBarry Smith row -= a->rstart; 1750752156aSBarry Smith 176a8c6a408SBarry Smith if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 177b97cf49bSBarry Smith 178b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 179b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 180b97cf49bSBarry Smith if (idx) { 181b97cf49bSBarry Smith itmp = a->j + a->i[row]; 182b97cf49bSBarry Smith if (*nz) { 183b97cf49bSBarry Smith *idx = itmp; 184b97cf49bSBarry Smith } 185b97cf49bSBarry Smith else *idx = 0; 186b97cf49bSBarry Smith } 1873a40ed3dSBarry Smith PetscFunctionReturn(0); 188b97cf49bSBarry Smith } 189b97cf49bSBarry Smith 190b97cf49bSBarry Smith #undef __FUNC__ 1910752156aSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj" 1920752156aSBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 193b97cf49bSBarry Smith { 194433994e6SBarry Smith PetscFunctionBegin; 1953a40ed3dSBarry Smith PetscFunctionReturn(0); 196b97cf49bSBarry Smith } 197b97cf49bSBarry Smith 198b97cf49bSBarry Smith #undef __FUNC__ 1990752156aSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj" 2000752156aSBarry Smith int MatGetBlockSize_MPIAdj(Mat A, int *bs) 201b97cf49bSBarry Smith { 202433994e6SBarry Smith PetscFunctionBegin; 203b97cf49bSBarry Smith *bs = 1; 2043a40ed3dSBarry Smith PetscFunctionReturn(0); 205b97cf49bSBarry Smith } 206b97cf49bSBarry Smith 207b97cf49bSBarry Smith 208b97cf49bSBarry Smith #undef __FUNC__ 2090752156aSBarry Smith #define __FUNC__ "MatEqual_MPIAdj" 2100752156aSBarry Smith int MatEqual_MPIAdj(Mat A,Mat B, PetscTruth* flg) 211b97cf49bSBarry Smith { 2120752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data; 2130f5bd95cSBarry Smith int ierr; 2140f5bd95cSBarry Smith PetscTruth flag; 215b97cf49bSBarry Smith 216433994e6SBarry Smith PetscFunctionBegin; 217a8c6a408SBarry Smith if (B->type != MATMPIADJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 218b97cf49bSBarry Smith 219b97cf49bSBarry Smith /* If the matrix dimensions are not equal, or no of nonzeros */ 2200752156aSBarry Smith if ((a->m != b->m ) ||( a->nz != b->nz)) { 2210f5bd95cSBarry Smith flag = PETSC_FALSE; 222b97cf49bSBarry Smith } 223b97cf49bSBarry Smith 224b97cf49bSBarry Smith /* if the a->i are the same */ 2250f5bd95cSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int),&flag);CHKERRQ(ierr); 226b97cf49bSBarry Smith 227b97cf49bSBarry Smith /* if a->j are the same */ 2280f5bd95cSBarry Smith ierr = PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int),&flag);CHKERRQ(ierr); 229b97cf49bSBarry Smith 230ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 2313a40ed3dSBarry Smith PetscFunctionReturn(0); 232b97cf49bSBarry Smith } 233b97cf49bSBarry Smith 234b97cf49bSBarry Smith 235b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2360b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2370b3b1487SBarry Smith MatGetRow_MPIAdj, 2380b3b1487SBarry Smith MatRestoreRow_MPIAdj, 239b97cf49bSBarry Smith 0, 240b97cf49bSBarry Smith 0, 241b97cf49bSBarry Smith 0, 242b97cf49bSBarry 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 0, 2510b3b1487SBarry Smith 0, 2520b3b1487SBarry Smith MatEqual_MPIAdj, 2530b3b1487SBarry Smith 0, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 2580b3b1487SBarry Smith 0, 2590b3b1487SBarry Smith MatSetOption_MPIAdj, 2600b3b1487SBarry Smith 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 2660b3b1487SBarry Smith MatGetSize_MPIAdj, 2670b3b1487SBarry Smith MatGetLocalSize_MPIAdj, 2680b3b1487SBarry Smith MatGetOwnershipRange_MPIAdj, 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, 2850b3b1487SBarry Smith 0, 2860b3b1487SBarry Smith 0, 287b97cf49bSBarry Smith 0, 2880752156aSBarry Smith MatGetBlockSize_MPIAdj, 289b97cf49bSBarry Smith 0, 290b97cf49bSBarry Smith 0, 291b97cf49bSBarry Smith 0, 292b97cf49bSBarry Smith 0, 293b97cf49bSBarry Smith 0, 2940752156aSBarry Smith 0, 2950752156aSBarry Smith 0, 2960b3b1487SBarry Smith 0, 2970b3b1487SBarry Smith 0, 2980b3b1487SBarry Smith 0, 2990b3b1487SBarry Smith 0, 3000b3b1487SBarry Smith 0, 3010b3b1487SBarry Smith MatGetMaps_Petsc}; 302b97cf49bSBarry Smith 303b97cf49bSBarry Smith 304b97cf49bSBarry Smith #undef __FUNC__ 3050752156aSBarry Smith #define __FUNC__ "MatCreateMPIAdj" 306b97cf49bSBarry Smith /*@C 3070752156aSBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 308b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 309b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 310b97cf49bSBarry Smith 311ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 312ef5ee4d1SLois Curfman McInnes 313b97cf49bSBarry Smith Input Parameters: 314ef5ee4d1SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 3150752156aSBarry Smith . m - number of local rows 316b97cf49bSBarry Smith . n - number of columns 317b97cf49bSBarry Smith . i - the indices into j for the start of each row 318ef5ee4d1SLois Curfman McInnes - j - the column indices for each row (sorted for each row). 319ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 320b97cf49bSBarry Smith 321b97cf49bSBarry Smith Output Parameter: 322b97cf49bSBarry Smith . A - the matrix 323b97cf49bSBarry Smith 32415091d37SBarry Smith Level: intermediate 32515091d37SBarry Smith 3264bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 3274bc6d8c0SBarry Smith MatSetValues(). 3284bc6d8c0SBarry Smith You must NOT free the ii and jj arrays yourself. PETSc will free them 329b97cf49bSBarry Smith when the matrix is destroyed. 330b97cf49bSBarry Smith 331ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 332b97cf49bSBarry Smith 33391e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering() 334b97cf49bSBarry Smith @*/ 3350752156aSBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j, Mat *A) 336b97cf49bSBarry Smith { 337b97cf49bSBarry Smith Mat B; 3380752156aSBarry Smith Mat_MPIAdj *b; 3390752156aSBarry Smith int ii,ierr, flg,size,rank; 340b97cf49bSBarry Smith 341433994e6SBarry Smith PetscFunctionBegin; 342d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 343d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 344b97cf49bSBarry Smith 345b97cf49bSBarry Smith *A = 0; 3463f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,"Mat",comm,MatDestroy,MatView); 347b97cf49bSBarry Smith PLogObjectCreate(B); 3480752156aSBarry Smith B->data = (void *) (b = PetscNew(Mat_MPIAdj));CHKPTRQ(b); 349549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 350549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 351e1311b90SBarry Smith B->ops->destroy = MatDestroy_MPIAdj; 352e1311b90SBarry Smith B->ops->view = MatView_MPIAdj; 353b97cf49bSBarry Smith B->factor = 0; 354b97cf49bSBarry Smith B->lupivotthreshold = 1.0; 355b97cf49bSBarry Smith B->mapping = 0; 356b97cf49bSBarry Smith B->assembled = PETSC_FALSE; 357b97cf49bSBarry Smith 3580752156aSBarry Smith b->m = m; B->m = m; 359ca161407SBarry Smith ierr = MPI_Allreduce(&m,&B->M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 3600752156aSBarry Smith B->n = n; B->N = n; 3610752156aSBarry Smith 362c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 363c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 364488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,B->M,&B->rmap);CHKERRQ(ierr); 365c7fcc2eaSBarry Smith /* we don't know the "local columns" so just use the row information :-( */ 366488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,B->M,&B->cmap);CHKERRQ(ierr); 367c7fcc2eaSBarry Smith 3680752156aSBarry Smith b->rowners = (int *) PetscMalloc((size+1)*sizeof(int));CHKPTRQ(b->rowners); 3690752156aSBarry Smith PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 370ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 3710752156aSBarry Smith b->rowners[0] = 0; 3720752156aSBarry Smith for ( ii=2; ii<=size; ii++ ) { 3730752156aSBarry Smith b->rowners[ii] += b->rowners[ii-1]; 3740752156aSBarry Smith } 3750752156aSBarry Smith b->rstart = b->rowners[rank]; 3760752156aSBarry Smith b->rend = b->rowners[rank+1]; 377b97cf49bSBarry Smith 378b97cf49bSBarry Smith b->j = j; 379b97cf49bSBarry Smith b->i = i; 380b97cf49bSBarry Smith 381b97cf49bSBarry Smith b->nz = i[m]; 382b97cf49bSBarry Smith b->diag = 0; 383b97cf49bSBarry Smith b->symmetric = PETSC_FALSE; 384b97cf49bSBarry Smith 385b97cf49bSBarry Smith *A = B; 386b97cf49bSBarry Smith 387b97cf49bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 388b97cf49bSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 389b97cf49bSBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 390b97cf49bSBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3913a40ed3dSBarry Smith PetscFunctionReturn(0); 392b97cf49bSBarry Smith } 393b97cf49bSBarry Smith 394b97cf49bSBarry Smith 395b97cf49bSBarry Smith 396433994e6SBarry Smith 397433994e6SBarry Smith 398433994e6SBarry Smith 399433994e6SBarry Smith 400433994e6SBarry Smith 401433994e6SBarry Smith 402