1b97cf49bSBarry Smith #ifdef PETSC_RCS_HEADER 2*0f5bd95cSBarry Smith static char vcid[] = "$Id: mpiadj.c,v 1.28 1999/10/01 21:21:33 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 FILE *fd; 18b97cf49bSBarry Smith char *outputname; 190752156aSBarry Smith MPI_Comm comm = A->comm; 20b97cf49bSBarry Smith 21433994e6SBarry Smith PetscFunctionBegin; 22b97cf49bSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 2377ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 24d132466eSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 25b97cf49bSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 263a40ed3dSBarry Smith PetscFunctionReturn(0); 270752156aSBarry Smith } else { 28b97cf49bSBarry Smith for ( i=0; i<m; i++ ) { 29d132466eSBarry Smith ierr = PetscSynchronizedFPrintf(comm,fd,"row %d:",i+a->rstart);CHKERRQ(ierr); 30b97cf49bSBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 31d132466eSBarry Smith ierr = PetscSynchronizedFPrintf(comm,fd," %d ",a->j[j]);CHKERRQ(ierr); 320752156aSBarry Smith } 33d132466eSBarry Smith ierr = PetscSynchronizedFPrintf(comm,fd,"\n");CHKERRQ(ierr); 34b97cf49bSBarry Smith } 35b97cf49bSBarry Smith } 36d132466eSBarry Smith ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 373a40ed3dSBarry Smith PetscFunctionReturn(0); 38b97cf49bSBarry Smith } 39b97cf49bSBarry Smith 40b97cf49bSBarry Smith #undef __FUNC__ 410752156aSBarry Smith #define __FUNC__ "MatView_MPIAdj" 42e1311b90SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer) 43b97cf49bSBarry Smith { 44b97cf49bSBarry Smith int ierr; 45*0f5bd95cSBarry Smith int isascii; 46b97cf49bSBarry Smith 47433994e6SBarry Smith PetscFunctionBegin; 48*0f5bd95cSBarry Smith isascii = PetscTypeCompare(viewer,ASCII_VIEWER); 49*0f5bd95cSBarry Smith if (isascii) { 503a40ed3dSBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 515cd90555SBarry Smith } else { 52*0f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 53b97cf49bSBarry Smith } 543a40ed3dSBarry Smith PetscFunctionReturn(0); 55b97cf49bSBarry Smith } 56b97cf49bSBarry Smith 57b97cf49bSBarry Smith #undef __FUNC__ 580752156aSBarry Smith #define __FUNC__ "MatDestroy_MPIAdj" 5994d884c6SBarry Smith int MatDestroy_MPIAdj(Mat mat) 60b97cf49bSBarry Smith { 6194d884c6SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) mat->data; 6294d884c6SBarry Smith int ierr; 6394d884c6SBarry Smith 6494d884c6SBarry Smith PetscFunctionBegin; 6594d884c6SBarry Smith 6694d884c6SBarry Smith if (mat->mapping) { 6794d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 6894d884c6SBarry Smith } 6994d884c6SBarry Smith if (mat->bmapping) { 7094d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 7194d884c6SBarry Smith } 72c7fcc2eaSBarry Smith if (mat->rmap) { 73c7fcc2eaSBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 74c7fcc2eaSBarry Smith } 75c7fcc2eaSBarry Smith if (mat->cmap) { 76c7fcc2eaSBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 77c7fcc2eaSBarry Smith } 78b97cf49bSBarry Smith 79aa482453SBarry Smith #if defined(PETSC_USE_LOG) 8094d884c6SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz); 81b97cf49bSBarry Smith #endif 82606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 83606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 84606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 85606d414cSSatish Balay ierr = PetscFree(a->rowners);CHKERRQ(ierr); 86606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 87b97cf49bSBarry Smith 8894d884c6SBarry Smith PLogObjectDestroy(mat); 8994d884c6SBarry Smith PetscHeaderDestroy(mat); 903a40ed3dSBarry Smith PetscFunctionReturn(0); 91b97cf49bSBarry Smith } 92b97cf49bSBarry Smith 93b97cf49bSBarry Smith 94b97cf49bSBarry Smith #undef __FUNC__ 950752156aSBarry Smith #define __FUNC__ "MatSetOption_MPIAdj" 960752156aSBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 97b97cf49bSBarry Smith { 980752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 99b97cf49bSBarry Smith 100433994e6SBarry Smith PetscFunctionBegin; 101b97cf49bSBarry Smith if (op == MAT_STRUCTURALLY_SYMMETRIC) { 102b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 103b97cf49bSBarry Smith } else { 104981c4779SBarry Smith PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 105b97cf49bSBarry Smith } 1063a40ed3dSBarry Smith PetscFunctionReturn(0); 107b97cf49bSBarry Smith } 108b97cf49bSBarry Smith 109b97cf49bSBarry Smith 110b97cf49bSBarry Smith /* 111b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 112b97cf49bSBarry Smith */ 113b97cf49bSBarry Smith 114b97cf49bSBarry Smith #undef __FUNC__ 1150752156aSBarry Smith #define __FUNC__ "MatMarkDiag_MPIAdj" 1160752156aSBarry Smith int MatMarkDiag_MPIAdj(Mat A) 117b97cf49bSBarry Smith { 1180752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 119b97cf49bSBarry Smith int i,j, *diag, m = a->m; 120b97cf49bSBarry Smith 121433994e6SBarry Smith PetscFunctionBegin; 122b97cf49bSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 123b97cf49bSBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 124b97cf49bSBarry Smith for ( i=0; i<a->m; i++ ) { 125b97cf49bSBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 126b97cf49bSBarry Smith if (a->j[j] == i) { 127b97cf49bSBarry Smith diag[i] = j; 128b97cf49bSBarry Smith break; 129b97cf49bSBarry Smith } 130b97cf49bSBarry Smith } 131b97cf49bSBarry Smith } 132b97cf49bSBarry Smith a->diag = diag; 1333a40ed3dSBarry Smith PetscFunctionReturn(0); 134b97cf49bSBarry Smith } 135b97cf49bSBarry Smith 136b97cf49bSBarry Smith #undef __FUNC__ 1370752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj" 1380752156aSBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n) 139b97cf49bSBarry Smith { 140433994e6SBarry Smith PetscFunctionBegin; 1410752156aSBarry Smith if (m) *m = A->M; 1420752156aSBarry Smith if (n) *n = A->N; 1433a40ed3dSBarry Smith PetscFunctionReturn(0); 144b97cf49bSBarry Smith } 145b97cf49bSBarry Smith 146b97cf49bSBarry Smith #undef __FUNC__ 1470752156aSBarry Smith #define __FUNC__ "MatGetSize_MPIAdj" 1480752156aSBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n) 149b97cf49bSBarry Smith { 1500752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 151433994e6SBarry Smith PetscFunctionBegin; 1520752156aSBarry Smith if (m) *m = a->m; 1530752156aSBarry Smith if (n) *n = A->N; 1543a40ed3dSBarry Smith PetscFunctionReturn(0); 155b97cf49bSBarry Smith } 156b97cf49bSBarry Smith 157b97cf49bSBarry Smith #undef __FUNC__ 1580752156aSBarry Smith #define __FUNC__ "MatGetOwnershipRange_MPIAdj" 1590752156aSBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n) 160b97cf49bSBarry Smith { 1610752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 162433994e6SBarry Smith PetscFunctionBegin; 1630752156aSBarry Smith *m = a->rstart; *n = a->rend; 1643a40ed3dSBarry Smith PetscFunctionReturn(0); 165b97cf49bSBarry Smith } 166b97cf49bSBarry Smith 167b97cf49bSBarry Smith #undef __FUNC__ 1680752156aSBarry Smith #define __FUNC__ "MatGetRow_MPIAdj" 1690752156aSBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 170b97cf49bSBarry Smith { 1710752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *) A->data; 172b97cf49bSBarry Smith int *itmp; 173b97cf49bSBarry Smith 174433994e6SBarry Smith PetscFunctionBegin; 1750752156aSBarry Smith row -= a->rstart; 1760752156aSBarry Smith 177a8c6a408SBarry Smith if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 178b97cf49bSBarry Smith 179b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 180b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 181b97cf49bSBarry Smith if (idx) { 182b97cf49bSBarry Smith itmp = a->j + a->i[row]; 183b97cf49bSBarry Smith if (*nz) { 184b97cf49bSBarry Smith *idx = itmp; 185b97cf49bSBarry Smith } 186b97cf49bSBarry Smith else *idx = 0; 187b97cf49bSBarry Smith } 1883a40ed3dSBarry Smith PetscFunctionReturn(0); 189b97cf49bSBarry Smith } 190b97cf49bSBarry Smith 191b97cf49bSBarry Smith #undef __FUNC__ 1920752156aSBarry Smith #define __FUNC__ "MatRestoreRow_MPIAdj" 1930752156aSBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 194b97cf49bSBarry Smith { 195433994e6SBarry Smith PetscFunctionBegin; 1963a40ed3dSBarry Smith PetscFunctionReturn(0); 197b97cf49bSBarry Smith } 198b97cf49bSBarry Smith 199b97cf49bSBarry Smith #undef __FUNC__ 2000752156aSBarry Smith #define __FUNC__ "MatGetBlockSize_MPIAdj" 2010752156aSBarry Smith int MatGetBlockSize_MPIAdj(Mat A, int *bs) 202b97cf49bSBarry Smith { 203433994e6SBarry Smith PetscFunctionBegin; 204b97cf49bSBarry Smith *bs = 1; 2053a40ed3dSBarry Smith PetscFunctionReturn(0); 206b97cf49bSBarry Smith } 207b97cf49bSBarry Smith 208b97cf49bSBarry Smith 209b97cf49bSBarry Smith #undef __FUNC__ 2100752156aSBarry Smith #define __FUNC__ "MatEqual_MPIAdj" 2110752156aSBarry Smith int MatEqual_MPIAdj(Mat A,Mat B, PetscTruth* flg) 212b97cf49bSBarry Smith { 2130752156aSBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data; 214*0f5bd95cSBarry Smith int ierr; 215*0f5bd95cSBarry Smith PetscTruth flag; 216b97cf49bSBarry Smith 217433994e6SBarry Smith PetscFunctionBegin; 218a8c6a408SBarry Smith if (B->type != MATMPIADJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 219b97cf49bSBarry Smith 220b97cf49bSBarry Smith /* If the matrix dimensions are not equal, or no of nonzeros */ 2210752156aSBarry Smith if ((a->m != b->m ) ||( a->nz != b->nz)) { 222*0f5bd95cSBarry Smith flag = PETSC_FALSE; 223b97cf49bSBarry Smith } 224b97cf49bSBarry Smith 225b97cf49bSBarry Smith /* if the a->i are the same */ 226*0f5bd95cSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int),&flag);CHKERRQ(ierr); 227b97cf49bSBarry Smith 228b97cf49bSBarry Smith /* if a->j are the same */ 229*0f5bd95cSBarry Smith ierr = PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int),&flag);CHKERRQ(ierr); 230b97cf49bSBarry Smith 231ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 2323a40ed3dSBarry Smith PetscFunctionReturn(0); 233b97cf49bSBarry Smith } 234b97cf49bSBarry Smith 235b97cf49bSBarry Smith 236b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2370b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2380b3b1487SBarry Smith MatGetRow_MPIAdj, 2390b3b1487SBarry Smith MatRestoreRow_MPIAdj, 240b97cf49bSBarry Smith 0, 241b97cf49bSBarry Smith 0, 242b97cf49bSBarry Smith 0, 243b97cf49bSBarry 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 0, 2530b3b1487SBarry Smith MatEqual_MPIAdj, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 2580b3b1487SBarry Smith 0, 2590b3b1487SBarry Smith 0, 2600b3b1487SBarry Smith MatSetOption_MPIAdj, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 2660b3b1487SBarry Smith 0, 2670b3b1487SBarry Smith MatGetSize_MPIAdj, 2680b3b1487SBarry Smith MatGetLocalSize_MPIAdj, 2690b3b1487SBarry Smith MatGetOwnershipRange_MPIAdj, 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, 2870b3b1487SBarry Smith 0, 288b97cf49bSBarry Smith 0, 2890752156aSBarry Smith MatGetBlockSize_MPIAdj, 290b97cf49bSBarry Smith 0, 291b97cf49bSBarry Smith 0, 292b97cf49bSBarry Smith 0, 293b97cf49bSBarry Smith 0, 294b97cf49bSBarry Smith 0, 2950752156aSBarry Smith 0, 2960752156aSBarry Smith 0, 2970b3b1487SBarry Smith 0, 2980b3b1487SBarry Smith 0, 2990b3b1487SBarry Smith 0, 3000b3b1487SBarry Smith 0, 3010b3b1487SBarry Smith 0, 3020b3b1487SBarry Smith MatGetMaps_Petsc}; 303b97cf49bSBarry Smith 304b97cf49bSBarry Smith 305b97cf49bSBarry Smith #undef __FUNC__ 3060752156aSBarry Smith #define __FUNC__ "MatCreateMPIAdj" 307b97cf49bSBarry Smith /*@C 3080752156aSBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 309b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 310b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 311b97cf49bSBarry Smith 312ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 313ef5ee4d1SLois Curfman McInnes 314b97cf49bSBarry Smith Input Parameters: 315ef5ee4d1SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 3160752156aSBarry Smith . m - number of local rows 317b97cf49bSBarry Smith . n - number of columns 318b97cf49bSBarry Smith . i - the indices into j for the start of each row 319ef5ee4d1SLois Curfman McInnes - j - the column indices for each row (sorted for each row). 320ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 321b97cf49bSBarry Smith 322b97cf49bSBarry Smith Output Parameter: 323b97cf49bSBarry Smith . A - the matrix 324b97cf49bSBarry Smith 32515091d37SBarry Smith Level: intermediate 32615091d37SBarry Smith 3274bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 3284bc6d8c0SBarry Smith MatSetValues(). 3294bc6d8c0SBarry Smith You must NOT free the ii and jj arrays yourself. PETSc will free them 330b97cf49bSBarry Smith when the matrix is destroyed. 331b97cf49bSBarry Smith 332ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 333b97cf49bSBarry Smith 33491e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering() 335b97cf49bSBarry Smith @*/ 3360752156aSBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j, Mat *A) 337b97cf49bSBarry Smith { 338b97cf49bSBarry Smith Mat B; 3390752156aSBarry Smith Mat_MPIAdj *b; 3400752156aSBarry Smith int ii,ierr, flg,size,rank; 341b97cf49bSBarry Smith 342433994e6SBarry Smith PetscFunctionBegin; 343d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 344d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 345b97cf49bSBarry Smith 346b97cf49bSBarry Smith *A = 0; 3473f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,"Mat",comm,MatDestroy,MatView); 348b97cf49bSBarry Smith PLogObjectCreate(B); 3490752156aSBarry Smith B->data = (void *) (b = PetscNew(Mat_MPIAdj));CHKPTRQ(b); 350549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 351549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 352e1311b90SBarry Smith B->ops->destroy = MatDestroy_MPIAdj; 353e1311b90SBarry Smith B->ops->view = MatView_MPIAdj; 354b97cf49bSBarry Smith B->factor = 0; 355b97cf49bSBarry Smith B->lupivotthreshold = 1.0; 356b97cf49bSBarry Smith B->mapping = 0; 357b97cf49bSBarry Smith B->assembled = PETSC_FALSE; 358b97cf49bSBarry Smith 3590752156aSBarry Smith b->m = m; B->m = m; 360ca161407SBarry Smith ierr = MPI_Allreduce(&m,&B->M,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 3610752156aSBarry Smith B->n = n; B->N = n; 3620752156aSBarry Smith 363c7fcc2eaSBarry Smith /* the information in the maps duplicates the information computed below, eventually 364c7fcc2eaSBarry Smith we should remove the duplicate information that is not contained in the maps */ 365488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,B->M,&B->rmap);CHKERRQ(ierr); 366c7fcc2eaSBarry Smith /* we don't know the "local columns" so just use the row information :-( */ 367488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,B->M,&B->cmap);CHKERRQ(ierr); 368c7fcc2eaSBarry Smith 3690752156aSBarry Smith b->rowners = (int *) PetscMalloc((size+1)*sizeof(int));CHKPTRQ(b->rowners); 3700752156aSBarry Smith PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 371ca161407SBarry Smith ierr = MPI_Allgather(&m,1,MPI_INT,b->rowners+1,1,MPI_INT,comm);CHKERRQ(ierr); 3720752156aSBarry Smith b->rowners[0] = 0; 3730752156aSBarry Smith for ( ii=2; ii<=size; ii++ ) { 3740752156aSBarry Smith b->rowners[ii] += b->rowners[ii-1]; 3750752156aSBarry Smith } 3760752156aSBarry Smith b->rstart = b->rowners[rank]; 3770752156aSBarry Smith b->rend = b->rowners[rank+1]; 378b97cf49bSBarry Smith 379b97cf49bSBarry Smith b->j = j; 380b97cf49bSBarry Smith b->i = i; 381b97cf49bSBarry Smith 382b97cf49bSBarry Smith b->nz = i[m]; 383b97cf49bSBarry Smith b->diag = 0; 384b97cf49bSBarry Smith b->symmetric = PETSC_FALSE; 385b97cf49bSBarry Smith 386b97cf49bSBarry Smith *A = B; 387b97cf49bSBarry Smith 388b97cf49bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 389b97cf49bSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 390b97cf49bSBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 391b97cf49bSBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3923a40ed3dSBarry Smith PetscFunctionReturn(0); 393b97cf49bSBarry Smith } 394b97cf49bSBarry Smith 395b97cf49bSBarry Smith 396b97cf49bSBarry Smith 397433994e6SBarry Smith 398433994e6SBarry Smith 399433994e6SBarry Smith 400433994e6SBarry Smith 401433994e6SBarry Smith 402433994e6SBarry Smith 403