1*ffac6cdbSBarry Smith /*$Id: mpiadj.c,v 1.49 2000/10/24 20:25: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; 14273d9f13SBarry 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); 22*ffac6cdbSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 23*ffac6cdbSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 240752156aSBarry Smith } else { 256831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 26b97cf49bSBarry Smith for (i=0; i<m; i++) { 276831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr); 28b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 296831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr); 300752156aSBarry Smith } 316831982aSBarry Smith ierr = ViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 32b97cf49bSBarry Smith } 336831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 34b97cf49bSBarry Smith } 356831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 363a40ed3dSBarry Smith PetscFunctionReturn(0); 37b97cf49bSBarry Smith } 38b97cf49bSBarry Smith 39b97cf49bSBarry Smith #undef __FUNC__ 40b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_MPIAdj" 413eda8832SBarry Smith int MatView_MPIAdj(Mat A,Viewer viewer) 42b97cf49bSBarry Smith { 43b97cf49bSBarry Smith int ierr; 446831982aSBarry Smith PetscTruth isascii; 45b97cf49bSBarry Smith 46433994e6SBarry Smith PetscFunctionBegin; 476831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 480f5bd95cSBarry Smith if (isascii) { 493eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 505cd90555SBarry Smith } else { 5129bbc08cSBarry Smith SETERRQ1(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__ 57b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_MPIAdj" 583eda8832SBarry Smith int MatDestroy_MPIAdj(Mat mat) 59b97cf49bSBarry Smith { 603eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 6194d884c6SBarry Smith int ierr; 6294d884c6SBarry Smith 6394d884c6SBarry Smith PetscFunctionBegin; 64aa482453SBarry Smith #if defined(PETSC_USE_LOG) 6594d884c6SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz); 66b97cf49bSBarry Smith #endif 67606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 680400557aSBarry Smith if (a->freeaij) { 69606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 70606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 711198db28SBarry Smith if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);} 720400557aSBarry Smith } 730400557aSBarry Smith ierr = PetscFree(a->rowners);CHKERRQ(ierr); 741198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 753a40ed3dSBarry Smith PetscFunctionReturn(0); 76b97cf49bSBarry Smith } 77b97cf49bSBarry Smith 78b97cf49bSBarry Smith #undef __FUNC__ 79b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_MPIAdj" 803eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 81b97cf49bSBarry Smith { 823eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 83b97cf49bSBarry Smith 84433994e6SBarry Smith PetscFunctionBegin; 85b97cf49bSBarry Smith if (op == MAT_STRUCTURALLY_SYMMETRIC) { 86b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 87b97cf49bSBarry Smith } else { 883eda8832SBarry Smith PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 89b97cf49bSBarry Smith } 903a40ed3dSBarry Smith PetscFunctionReturn(0); 91b97cf49bSBarry Smith } 92b97cf49bSBarry Smith 93b97cf49bSBarry Smith 94b97cf49bSBarry Smith /* 95b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 96b97cf49bSBarry Smith */ 97b97cf49bSBarry Smith 98b97cf49bSBarry Smith #undef __FUNC__ 99b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMarkDiagonal_MPIAdj" 1003eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A) 101b97cf49bSBarry Smith { 1023eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 103273d9f13SBarry Smith int i,j,*diag,m = A->m; 104b97cf49bSBarry Smith 105433994e6SBarry Smith PetscFunctionBegin; 106b97cf49bSBarry Smith diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag); 107b97cf49bSBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 108273d9f13SBarry Smith for (i=0; i<A->m; i++) { 109b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 110b97cf49bSBarry Smith if (a->j[j] == i) { 111b97cf49bSBarry Smith diag[i] = j; 112b97cf49bSBarry Smith break; 113b97cf49bSBarry Smith } 114b97cf49bSBarry Smith } 115b97cf49bSBarry Smith } 116b97cf49bSBarry Smith a->diag = diag; 1173a40ed3dSBarry Smith PetscFunctionReturn(0); 118b97cf49bSBarry Smith } 119b97cf49bSBarry Smith 120b97cf49bSBarry Smith #undef __FUNC__ 121b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_MPIAdj" 1223eda8832SBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n) 123b97cf49bSBarry Smith { 1243eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 125433994e6SBarry Smith PetscFunctionBegin; 126c2e958feSBarry Smith if (m) *m = a->rstart; 127c2e958feSBarry Smith if (n) *n = a->rend; 1283a40ed3dSBarry Smith PetscFunctionReturn(0); 129b97cf49bSBarry Smith } 130b97cf49bSBarry Smith 131b97cf49bSBarry Smith #undef __FUNC__ 132b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_MPIAdj" 1333eda8832SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 134b97cf49bSBarry Smith { 1353eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 136b97cf49bSBarry Smith int *itmp; 137b97cf49bSBarry Smith 138433994e6SBarry Smith PetscFunctionBegin; 1390752156aSBarry Smith row -= a->rstart; 1400752156aSBarry Smith 141273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 142b97cf49bSBarry Smith 143b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 144b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 145b97cf49bSBarry Smith if (idx) { 146b97cf49bSBarry Smith itmp = a->j + a->i[row]; 147b97cf49bSBarry Smith if (*nz) { 148b97cf49bSBarry Smith *idx = itmp; 149b97cf49bSBarry Smith } 150b97cf49bSBarry Smith else *idx = 0; 151b97cf49bSBarry Smith } 1523a40ed3dSBarry Smith PetscFunctionReturn(0); 153b97cf49bSBarry Smith } 154b97cf49bSBarry Smith 155b97cf49bSBarry Smith #undef __FUNC__ 156b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_MPIAdj" 1573eda8832SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 158b97cf49bSBarry Smith { 159433994e6SBarry Smith PetscFunctionBegin; 1603a40ed3dSBarry Smith PetscFunctionReturn(0); 161b97cf49bSBarry Smith } 162b97cf49bSBarry Smith 163b97cf49bSBarry Smith #undef __FUNC__ 164b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_MPIAdj" 1653eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs) 166b97cf49bSBarry Smith { 167433994e6SBarry Smith PetscFunctionBegin; 168b97cf49bSBarry Smith *bs = 1; 1693a40ed3dSBarry Smith PetscFunctionReturn(0); 170b97cf49bSBarry Smith } 171b97cf49bSBarry Smith 172b97cf49bSBarry Smith 173b97cf49bSBarry Smith #undef __FUNC__ 174b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatEqual_MPIAdj" 1753eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 176b97cf49bSBarry Smith { 1773eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 1780f5bd95cSBarry Smith int ierr; 1790f5bd95cSBarry Smith PetscTruth flag; 180b97cf49bSBarry Smith 181433994e6SBarry Smith PetscFunctionBegin; 182273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr); 183273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type"); 184b97cf49bSBarry Smith 185b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 186273d9f13SBarry Smith if ((A->m != B->m) ||(a->nz != b->nz)) { 1870f5bd95cSBarry Smith flag = PETSC_FALSE; 188b97cf49bSBarry Smith } 189b97cf49bSBarry Smith 190b97cf49bSBarry Smith /* if the a->i are the same */ 191273d9f13SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr); 192b97cf49bSBarry Smith 193b97cf49bSBarry Smith /* if a->j are the same */ 1940f5bd95cSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr); 195b97cf49bSBarry Smith 196ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 1973a40ed3dSBarry Smith PetscFunctionReturn(0); 198b97cf49bSBarry Smith } 199b97cf49bSBarry Smith 200b97cf49bSBarry Smith 201b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2020b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2033eda8832SBarry Smith MatGetRow_MPIAdj, 2043eda8832SBarry Smith MatRestoreRow_MPIAdj, 205b97cf49bSBarry Smith 0, 206b97cf49bSBarry Smith 0, 207b97cf49bSBarry Smith 0, 208b97cf49bSBarry Smith 0, 2090b3b1487SBarry Smith 0, 2100b3b1487SBarry Smith 0, 2110b3b1487SBarry Smith 0, 2120b3b1487SBarry Smith 0, 2130b3b1487SBarry Smith 0, 2140b3b1487SBarry Smith 0, 2150b3b1487SBarry Smith 0, 2160b3b1487SBarry Smith 0, 2170b3b1487SBarry Smith 0, 2183eda8832SBarry Smith MatEqual_MPIAdj, 2190b3b1487SBarry Smith 0, 2200b3b1487SBarry Smith 0, 2210b3b1487SBarry Smith 0, 2220b3b1487SBarry Smith 0, 2230b3b1487SBarry Smith 0, 2240b3b1487SBarry Smith 0, 2253eda8832SBarry Smith MatSetOption_MPIAdj, 2260b3b1487SBarry Smith 0, 2270b3b1487SBarry Smith 0, 2280b3b1487SBarry Smith 0, 2290b3b1487SBarry Smith 0, 2300b3b1487SBarry Smith 0, 2310b3b1487SBarry Smith 0, 232273d9f13SBarry Smith 0, 233273d9f13SBarry Smith 0, 2343eda8832SBarry Smith MatGetOwnershipRange_MPIAdj, 2350b3b1487SBarry Smith 0, 2360b3b1487SBarry Smith 0, 2370b3b1487SBarry Smith 0, 2380b3b1487SBarry Smith 0, 2390b3b1487SBarry Smith 0, 2400b3b1487SBarry 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 0, 2510b3b1487SBarry Smith 0, 2520b3b1487SBarry Smith 0, 253b97cf49bSBarry Smith 0, 2543eda8832SBarry Smith MatGetBlockSize_MPIAdj, 255b97cf49bSBarry Smith 0, 256b97cf49bSBarry Smith 0, 257b97cf49bSBarry Smith 0, 258b97cf49bSBarry Smith 0, 259b97cf49bSBarry Smith 0, 2600752156aSBarry Smith 0, 2610752156aSBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 265b9b97703SBarry Smith MatDestroy_MPIAdj, 266b9b97703SBarry Smith MatView_MPIAdj, 2670b3b1487SBarry Smith MatGetMaps_Petsc}; 268b97cf49bSBarry Smith 269b97cf49bSBarry Smith 270273d9f13SBarry Smith EXTERN_C_BEGIN 271273d9f13SBarry Smith #undef __FUNC__ 272273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_MPIAdj" 273273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B) 274273d9f13SBarry Smith { 275273d9f13SBarry Smith Mat_MPIAdj *b; 276273d9f13SBarry Smith int ii,ierr,size,rank; 277273d9f13SBarry Smith 278273d9f13SBarry Smith PetscFunctionBegin; 279273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 280273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 281273d9f13SBarry Smith 282273d9f13SBarry Smith B->data = (void*)(b = PetscNew(Mat_MPIAdj));CHKPTRQ(b); 283273d9f13SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 284273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 285273d9f13SBarry Smith B->factor = 0; 286273d9f13SBarry Smith B->lupivotthreshold = 1.0; 287273d9f13SBarry Smith B->mapping = 0; 288273d9f13SBarry Smith B->assembled = PETSC_FALSE; 289273d9f13SBarry Smith 290273d9f13SBarry Smith ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr); 291273d9f13SBarry Smith B->n = B->N; 292273d9f13SBarry Smith 293273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 294273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 295273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 296273d9f13SBarry Smith /* we don't know the "local columns" so just use the row information :-(*/ 297273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr); 298273d9f13SBarry Smith 299273d9f13SBarry Smith b->rowners = (int*)PetscMalloc((size+1)*sizeof(int));CHKPTRQ(b->rowners); 300273d9f13SBarry Smith PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 301273d9f13SBarry Smith ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr); 302273d9f13SBarry Smith b->rowners[0] = 0; 303273d9f13SBarry Smith for (ii=2; ii<=size; ii++) { 304273d9f13SBarry Smith b->rowners[ii] += b->rowners[ii-1]; 305273d9f13SBarry Smith } 306273d9f13SBarry Smith b->rstart = b->rowners[rank]; 307273d9f13SBarry Smith b->rend = b->rowners[rank+1]; 308273d9f13SBarry Smith 309273d9f13SBarry Smith PetscFunctionReturn(0); 310273d9f13SBarry Smith } 311273d9f13SBarry Smith EXTERN_C_END 312273d9f13SBarry Smith 313273d9f13SBarry Smith #undef __FUNC__ 314273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMPIAdjSetPreallocation" 315273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values) 316273d9f13SBarry Smith { 317273d9f13SBarry Smith Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 318273d9f13SBarry Smith int ierr; 319273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g) 320273d9f13SBarry Smith int ii; 321273d9f13SBarry Smith #endif 322273d9f13SBarry Smith 323273d9f13SBarry Smith PetscFunctionBegin; 324273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 325273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g) 326273d9f13SBarry Smith if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]); 327273d9f13SBarry Smith for (ii=1; ii<B->m; ii++) { 328273d9f13SBarry Smith if (i[ii] < 0 || i[ii] < i[ii-1]) { 329273d9f13SBarry Smith SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]); 330273d9f13SBarry Smith } 331273d9f13SBarry Smith } 332273d9f13SBarry Smith for (ii=0; ii<i[B->m]; ii++) { 333273d9f13SBarry Smith if (j[ii] < 0 || j[ii] >= B->N) { 334273d9f13SBarry Smith SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]); 335273d9f13SBarry Smith } 336273d9f13SBarry Smith } 337273d9f13SBarry Smith #endif 338273d9f13SBarry Smith 339273d9f13SBarry Smith b->j = j; 340273d9f13SBarry Smith b->i = i; 341273d9f13SBarry Smith b->values = values; 342273d9f13SBarry Smith 343273d9f13SBarry Smith b->nz = i[B->m]; 344273d9f13SBarry Smith b->diag = 0; 345273d9f13SBarry Smith b->symmetric = PETSC_FALSE; 346273d9f13SBarry Smith b->freeaij = PETSC_TRUE; 347273d9f13SBarry Smith 348273d9f13SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 349273d9f13SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 350273d9f13SBarry Smith PetscFunctionReturn(0); 351273d9f13SBarry Smith } 352273d9f13SBarry Smith 353b97cf49bSBarry Smith #undef __FUNC__ 354b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIAdj" 355b97cf49bSBarry Smith /*@C 3563eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 357b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 358b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 359b97cf49bSBarry Smith 360ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 361ef5ee4d1SLois Curfman McInnes 362b97cf49bSBarry Smith Input Parameters: 363c2e958feSBarry Smith + comm - MPI communicator 3640752156aSBarry Smith . m - number of local rows 365b97cf49bSBarry Smith . n - number of columns 366b97cf49bSBarry Smith . i - the indices into j for the start of each row 367329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 368ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 369329f5518SBarry Smith - values -[optional] edge weights 370b97cf49bSBarry Smith 371b97cf49bSBarry Smith Output Parameter: 372b97cf49bSBarry Smith . A - the matrix 373b97cf49bSBarry Smith 37415091d37SBarry Smith Level: intermediate 37515091d37SBarry Smith 3764bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 3774bc6d8c0SBarry Smith MatSetValues(). 378c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 3791198db28SBarry Smith when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you 3801198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 381327e1a73SBarry Smith Should not include the matrix diagonals. 382b97cf49bSBarry Smith 383ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 384b97cf49bSBarry Smith 38591e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering() 386b97cf49bSBarry Smith @*/ 3873eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A) 388b97cf49bSBarry Smith { 389273d9f13SBarry Smith int ierr; 390b97cf49bSBarry Smith 391433994e6SBarry Smith PetscFunctionBegin; 392273d9f13SBarry Smith ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr); 393273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 394273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 3953a40ed3dSBarry Smith PetscFunctionReturn(0); 396b97cf49bSBarry Smith } 397b97cf49bSBarry Smith 398273d9f13SBarry Smith EXTERN_C_BEGIN 399c2e958feSBarry Smith #undef __FUNC__ 400273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatConvertTo_MPIAdj" 401273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B) 402c2e958feSBarry Smith { 4034c49b128SBarry Smith int i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a; 404c2e958feSBarry Smith Scalar *ra; 405c2e958feSBarry Smith MPI_Comm comm; 406c2e958feSBarry Smith 407c2e958feSBarry Smith PetscFunctionBegin; 408c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 409c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 410c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 411c2e958feSBarry Smith 412c2e958feSBarry Smith /* count the number of nonzeros per row */ 413c2e958feSBarry Smith for (i=0; i<m; i++) { 414c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 415c2e958feSBarry Smith for (j=0; j<len; j++) { 416c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 417c2e958feSBarry Smith } 418c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 419c2e958feSBarry Smith nzeros += len; 420c2e958feSBarry Smith } 421c2e958feSBarry Smith 422c2e958feSBarry Smith /* malloc space for nonzeros */ 423c2e958feSBarry Smith a = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(a); 424c2e958feSBarry Smith ia = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(ia); 425c2e958feSBarry Smith ja = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(ja); 426c2e958feSBarry Smith 427c2e958feSBarry Smith nzeros = 0; 428c2e958feSBarry Smith ia[0] = 0; 429c2e958feSBarry Smith for (i=0; i<m; i++) { 430c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 431c2e958feSBarry Smith cnt = 0; 432c2e958feSBarry Smith for (j=0; j<len; j++) { 433c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 4344c49b128SBarry Smith a[nzeros+cnt] = (int) PetscAbsScalar(ra[j]); 435c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 436c2e958feSBarry Smith } 437c2e958feSBarry Smith } 438c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 439c2e958feSBarry Smith nzeros += cnt; 440c2e958feSBarry Smith ia[i+1] = nzeros; 441c2e958feSBarry Smith } 442c2e958feSBarry Smith 443c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 4443eda8832SBarry Smith ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr); 445c2e958feSBarry Smith 446c2e958feSBarry Smith PetscFunctionReturn(0); 447c2e958feSBarry Smith } 448273d9f13SBarry Smith EXTERN_C_END 449433994e6SBarry Smith 450433994e6SBarry Smith 451433994e6SBarry Smith 452433994e6SBarry Smith 453433994e6SBarry Smith 454