1*b9b97703SBarry Smith /*$Id: mpiadj.c,v 1.41 2000/05/10 16:40:58 bsmith Exp bsmith $*/ 2b97cf49bSBarry Smith 3b97cf49bSBarry Smith /* 4b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 5b97cf49bSBarry Smith */ 6e090d566SSatish Balay #include "petscsys.h" 73eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h" 8b97cf49bSBarry Smith 9b97cf49bSBarry Smith #undef __FUNC__ 10b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIAdj_ASCII" 11ca44d042SBarry Smith int MatView_MPIAdj_ASCII(Mat A,Viewer viewer) 12b97cf49bSBarry Smith { 133eda8832SBarry 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__ 38b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIAdj" 393eda8832SBarry 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) { 473eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 485cd90555SBarry Smith } else { 493eda8832SBarry 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__ 55b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_MPIAdj" 563eda8832SBarry Smith int MatDestroy_MPIAdj(Mat mat) 57b97cf49bSBarry Smith { 583eda8832SBarry 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);} 80c2e958feSBarry Smith if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);} 81606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 82606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 83606d414cSSatish Balay ierr = PetscFree(a->rowners);CHKERRQ(ierr); 84606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 85b97cf49bSBarry Smith 8694d884c6SBarry Smith PLogObjectDestroy(mat); 8794d884c6SBarry Smith PetscHeaderDestroy(mat); 883a40ed3dSBarry Smith PetscFunctionReturn(0); 89b97cf49bSBarry Smith } 90b97cf49bSBarry Smith 91b97cf49bSBarry Smith 92b97cf49bSBarry Smith #undef __FUNC__ 93b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_MPIAdj" 943eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 95b97cf49bSBarry Smith { 963eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 97b97cf49bSBarry Smith 98433994e6SBarry Smith PetscFunctionBegin; 99b97cf49bSBarry Smith if (op == MAT_STRUCTURALLY_SYMMETRIC) { 100b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 101b97cf49bSBarry Smith } else { 1023eda8832SBarry Smith PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 103b97cf49bSBarry Smith } 1043a40ed3dSBarry Smith PetscFunctionReturn(0); 105b97cf49bSBarry Smith } 106b97cf49bSBarry Smith 107b97cf49bSBarry Smith 108b97cf49bSBarry Smith /* 109b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 110b97cf49bSBarry Smith */ 111b97cf49bSBarry Smith 112b97cf49bSBarry Smith #undef __FUNC__ 113b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMarkDiagonal_MPIAdj" 1143eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A) 115b97cf49bSBarry Smith { 1163eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 117b97cf49bSBarry Smith int i,j,*diag,m = a->m; 118b97cf49bSBarry Smith 119433994e6SBarry Smith PetscFunctionBegin; 120b97cf49bSBarry Smith diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag); 121b97cf49bSBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 122b97cf49bSBarry Smith for (i=0; i<a->m; i++) { 123b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 124b97cf49bSBarry Smith if (a->j[j] == i) { 125b97cf49bSBarry Smith diag[i] = j; 126b97cf49bSBarry Smith break; 127b97cf49bSBarry Smith } 128b97cf49bSBarry Smith } 129b97cf49bSBarry Smith } 130b97cf49bSBarry Smith a->diag = diag; 1313a40ed3dSBarry Smith PetscFunctionReturn(0); 132b97cf49bSBarry Smith } 133b97cf49bSBarry Smith 134b97cf49bSBarry Smith #undef __FUNC__ 135b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSize_MPIAdj" 1363eda8832SBarry Smith int MatGetSize_MPIAdj(Mat A,int *m,int *n) 137b97cf49bSBarry Smith { 138433994e6SBarry Smith PetscFunctionBegin; 1390752156aSBarry Smith if (m) *m = A->M; 1400752156aSBarry Smith if (n) *n = A->N; 1413a40ed3dSBarry Smith PetscFunctionReturn(0); 142b97cf49bSBarry Smith } 143b97cf49bSBarry Smith 144b97cf49bSBarry Smith #undef __FUNC__ 145b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSize_MPIAdj" 1463eda8832SBarry Smith int MatGetLocalSize_MPIAdj(Mat A,int *m,int *n) 147b97cf49bSBarry Smith { 1483eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 149433994e6SBarry Smith PetscFunctionBegin; 1500752156aSBarry Smith if (m) *m = a->m; 1510752156aSBarry Smith if (n) *n = A->N; 1523a40ed3dSBarry Smith PetscFunctionReturn(0); 153b97cf49bSBarry Smith } 154b97cf49bSBarry Smith 155b97cf49bSBarry Smith #undef __FUNC__ 156b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_MPIAdj" 1573eda8832SBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n) 158b97cf49bSBarry Smith { 1593eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 160433994e6SBarry Smith PetscFunctionBegin; 161c2e958feSBarry Smith if (m) *m = a->rstart; 162c2e958feSBarry Smith if (n) *n = a->rend; 1633a40ed3dSBarry Smith PetscFunctionReturn(0); 164b97cf49bSBarry Smith } 165b97cf49bSBarry Smith 166b97cf49bSBarry Smith #undef __FUNC__ 167b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_MPIAdj" 1683eda8832SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 169b97cf49bSBarry Smith { 1703eda8832SBarry 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__ 191b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_MPIAdj" 1923eda8832SBarry 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__ 199b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_MPIAdj" 2003eda8832SBarry 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__ 209b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatEqual_MPIAdj" 2103eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 211b97cf49bSBarry Smith { 2123eda8832SBarry 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; 2173eda8832SBarry 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, 2373eda8832SBarry Smith MatGetRow_MPIAdj, 2383eda8832SBarry 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, 2523eda8832SBarry Smith MatEqual_MPIAdj, 2530b3b1487SBarry Smith 0, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 2580b3b1487SBarry Smith 0, 2593eda8832SBarry Smith MatSetOption_MPIAdj, 2600b3b1487SBarry Smith 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 2663eda8832SBarry Smith MatGetSize_MPIAdj, 2673eda8832SBarry Smith MatGetLocalSize_MPIAdj, 2683eda8832SBarry 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, 2883eda8832SBarry 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, 299*b9b97703SBarry Smith MatDestroy_MPIAdj, 300*b9b97703SBarry Smith MatView_MPIAdj, 3010b3b1487SBarry Smith MatGetMaps_Petsc}; 302b97cf49bSBarry Smith 303b97cf49bSBarry Smith 304b97cf49bSBarry Smith #undef __FUNC__ 305b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIAdj" 306b97cf49bSBarry Smith /*@C 3073eda8832SBarry 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: 314c2e958feSBarry Smith + comm - MPI communicator 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 318329f5518SBarry Smith . 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). 320329f5518SBarry Smith - values -[optional] edge weights 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(). 329c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 330c2e958feSBarry Smith when the matrix is destroyed. And you must allocate them with PetscMalloc() 331b97cf49bSBarry Smith 332ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 333b97cf49bSBarry Smith 33491e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering() 335b97cf49bSBarry Smith @*/ 3363eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A) 337b97cf49bSBarry Smith { 338b97cf49bSBarry Smith Mat B; 3393eda8832SBarry Smith Mat_MPIAdj *b; 340f1af5d2fSBarry Smith int ii,ierr,size,rank; 341f1af5d2fSBarry Smith PetscTruth flg; 342b97cf49bSBarry Smith 343433994e6SBarry Smith PetscFunctionBegin; 344d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 345d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 346b97cf49bSBarry Smith 347b97cf49bSBarry Smith *A = 0; 3483eda8832SBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATMPIADJ,"Mat",comm,MatDestroy,MatView); 349b97cf49bSBarry Smith PLogObjectCreate(B); 3503eda8832SBarry Smith B->data = (void*)(b = PetscNew(Mat_MPIAdj));CHKPTRQ(b); 3513eda8832SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 352549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 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); 3693eda8832SBarry 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 378c2e958feSBarry Smith #if defined(PETSC_BOPT_g) 379c2e958feSBarry Smith if (i[0] != 0) SETERR1(1,1,"First i[] index must be zero, instead it is %d\n",i[0]); 380c2e958feSBarry Smith for (ii=1; ii<m; ii++) { 381c2e958feSBarry Smith if (i[ii] < 0 || i[ii] > i[ii-1]) { 382c2e958feSBarry Smith SETERRQ4(1,1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]); 383c2e958feSBarry Smith } 384c2e958feSBarry Smith } 385c2e958feSBarry Smith for (ii=0; ii<i[m]; ii++) { 386c2e958feSBarry Smith if (i[ii] < 0 || i[ii] >= N) { 387c2e958feSBarry Smith SETERRQ2(1,1,"Column index %d out of range %d\n",ii,i[ii]); 388c2e958feSBarry Smith } 389c2e958feSBarry Smith } 390c2e958feSBarry Smith #endif 391c2e958feSBarry Smith 392b97cf49bSBarry Smith b->j = j; 393b97cf49bSBarry Smith b->i = i; 394329f5518SBarry Smith b->values = values; 395b97cf49bSBarry Smith 396b97cf49bSBarry Smith b->nz = i[m]; 397b97cf49bSBarry Smith b->diag = 0; 398b97cf49bSBarry Smith b->symmetric = PETSC_FALSE; 399b97cf49bSBarry Smith 400b97cf49bSBarry Smith *A = B; 401b97cf49bSBarry Smith 402b97cf49bSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 403b97cf49bSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 404b97cf49bSBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 405b97cf49bSBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4063a40ed3dSBarry Smith PetscFunctionReturn(0); 407b97cf49bSBarry Smith } 408b97cf49bSBarry Smith 409b97cf49bSBarry Smith 410b97cf49bSBarry Smith 411c2e958feSBarry Smith #undef __FUNC__ 412b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatConvert_MPIAdj" 4133eda8832SBarry Smith int MatConvert_MPIAdj(Mat A,MatType type,Mat *B) 414c2e958feSBarry Smith { 4154c49b128SBarry Smith int i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a; 416c2e958feSBarry Smith Scalar *ra; 417c2e958feSBarry Smith MPI_Comm comm; 418c2e958feSBarry Smith 419c2e958feSBarry Smith PetscFunctionBegin; 420c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 421c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 422c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 423c2e958feSBarry Smith 424c2e958feSBarry Smith /* count the number of nonzeros per row */ 425c2e958feSBarry Smith for (i=0; i<m; i++) { 426c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 427c2e958feSBarry Smith for (j=0; j<len; j++) { 428c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 429c2e958feSBarry Smith } 430c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 431c2e958feSBarry Smith nzeros += len; 432c2e958feSBarry Smith } 433c2e958feSBarry Smith 434c2e958feSBarry Smith /* malloc space for nonzeros */ 435c2e958feSBarry Smith a = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(a); 436c2e958feSBarry Smith ia = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(ia); 437c2e958feSBarry Smith ja = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(ja); 438c2e958feSBarry Smith 439c2e958feSBarry Smith nzeros = 0; 440c2e958feSBarry Smith ia[0] = 0; 441c2e958feSBarry Smith for (i=0; i<m; i++) { 442c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 443c2e958feSBarry Smith cnt = 0; 444c2e958feSBarry Smith for (j=0; j<len; j++) { 445c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 4464c49b128SBarry Smith a[nzeros+cnt] = (int) PetscAbsScalar(ra[j]); 447c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 448c2e958feSBarry Smith } 449c2e958feSBarry Smith } 450c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 451c2e958feSBarry Smith nzeros += cnt; 452c2e958feSBarry Smith ia[i+1] = nzeros; 453c2e958feSBarry Smith } 454c2e958feSBarry Smith 455c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 4563eda8832SBarry Smith ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr); 457c2e958feSBarry Smith 458c2e958feSBarry Smith PetscFunctionReturn(0); 459c2e958feSBarry Smith } 460433994e6SBarry Smith 461433994e6SBarry Smith 462433994e6SBarry Smith 463433994e6SBarry Smith 464433994e6SBarry Smith 465433994e6SBarry Smith 466