1*273d9f13SBarry Smith /*$Id: mpiadj.c,v 1.48 2000/09/28 21:11:35 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; 14*273d9f13SBarry 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 { 4929bbc08cSBarry Smith SETERRQ1(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; 62aa482453SBarry Smith #if defined(PETSC_USE_LOG) 6394d884c6SBarry Smith PLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz); 64b97cf49bSBarry Smith #endif 65606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 660400557aSBarry Smith if (a->freeaij) { 67606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 68606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 691198db28SBarry Smith if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);} 700400557aSBarry Smith } 710400557aSBarry Smith ierr = PetscFree(a->rowners);CHKERRQ(ierr); 721198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 733a40ed3dSBarry Smith PetscFunctionReturn(0); 74b97cf49bSBarry Smith } 75b97cf49bSBarry Smith 76b97cf49bSBarry Smith #undef __FUNC__ 77b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_MPIAdj" 783eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 79b97cf49bSBarry Smith { 803eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 81b97cf49bSBarry Smith 82433994e6SBarry Smith PetscFunctionBegin; 83b97cf49bSBarry Smith if (op == MAT_STRUCTURALLY_SYMMETRIC) { 84b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 85b97cf49bSBarry Smith } else { 863eda8832SBarry Smith PLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 87b97cf49bSBarry Smith } 883a40ed3dSBarry Smith PetscFunctionReturn(0); 89b97cf49bSBarry Smith } 90b97cf49bSBarry Smith 91b97cf49bSBarry Smith 92b97cf49bSBarry Smith /* 93b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 94b97cf49bSBarry Smith */ 95b97cf49bSBarry Smith 96b97cf49bSBarry Smith #undef __FUNC__ 97b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMarkDiagonal_MPIAdj" 983eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A) 99b97cf49bSBarry Smith { 1003eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 101*273d9f13SBarry Smith int i,j,*diag,m = A->m; 102b97cf49bSBarry Smith 103433994e6SBarry Smith PetscFunctionBegin; 104b97cf49bSBarry Smith diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag); 105b97cf49bSBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 106*273d9f13SBarry Smith for (i=0; i<A->m; i++) { 107b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 108b97cf49bSBarry Smith if (a->j[j] == i) { 109b97cf49bSBarry Smith diag[i] = j; 110b97cf49bSBarry Smith break; 111b97cf49bSBarry Smith } 112b97cf49bSBarry Smith } 113b97cf49bSBarry Smith } 114b97cf49bSBarry Smith a->diag = diag; 1153a40ed3dSBarry Smith PetscFunctionReturn(0); 116b97cf49bSBarry Smith } 117b97cf49bSBarry Smith 118b97cf49bSBarry Smith #undef __FUNC__ 119b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_MPIAdj" 1203eda8832SBarry Smith int MatGetOwnershipRange_MPIAdj(Mat A,int *m,int *n) 121b97cf49bSBarry Smith { 1223eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 123433994e6SBarry Smith PetscFunctionBegin; 124c2e958feSBarry Smith if (m) *m = a->rstart; 125c2e958feSBarry Smith if (n) *n = a->rend; 1263a40ed3dSBarry Smith PetscFunctionReturn(0); 127b97cf49bSBarry Smith } 128b97cf49bSBarry Smith 129b97cf49bSBarry Smith #undef __FUNC__ 130b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_MPIAdj" 1313eda8832SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 132b97cf49bSBarry Smith { 1333eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 134b97cf49bSBarry Smith int *itmp; 135b97cf49bSBarry Smith 136433994e6SBarry Smith PetscFunctionBegin; 1370752156aSBarry Smith row -= a->rstart; 1380752156aSBarry Smith 139*273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 140b97cf49bSBarry Smith 141b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 142b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 143b97cf49bSBarry Smith if (idx) { 144b97cf49bSBarry Smith itmp = a->j + a->i[row]; 145b97cf49bSBarry Smith if (*nz) { 146b97cf49bSBarry Smith *idx = itmp; 147b97cf49bSBarry Smith } 148b97cf49bSBarry Smith else *idx = 0; 149b97cf49bSBarry Smith } 1503a40ed3dSBarry Smith PetscFunctionReturn(0); 151b97cf49bSBarry Smith } 152b97cf49bSBarry Smith 153b97cf49bSBarry Smith #undef __FUNC__ 154b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_MPIAdj" 1553eda8832SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,Scalar **v) 156b97cf49bSBarry Smith { 157433994e6SBarry Smith PetscFunctionBegin; 1583a40ed3dSBarry Smith PetscFunctionReturn(0); 159b97cf49bSBarry Smith } 160b97cf49bSBarry Smith 161b97cf49bSBarry Smith #undef __FUNC__ 162b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_MPIAdj" 1633eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs) 164b97cf49bSBarry Smith { 165433994e6SBarry Smith PetscFunctionBegin; 166b97cf49bSBarry Smith *bs = 1; 1673a40ed3dSBarry Smith PetscFunctionReturn(0); 168b97cf49bSBarry Smith } 169b97cf49bSBarry Smith 170b97cf49bSBarry Smith 171b97cf49bSBarry Smith #undef __FUNC__ 172b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatEqual_MPIAdj" 1733eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 174b97cf49bSBarry Smith { 1753eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 1760f5bd95cSBarry Smith int ierr; 1770f5bd95cSBarry Smith PetscTruth flag; 178b97cf49bSBarry Smith 179433994e6SBarry Smith PetscFunctionBegin; 180*273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr); 181*273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type"); 182b97cf49bSBarry Smith 183b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 184*273d9f13SBarry Smith if ((A->m != B->m) ||(a->nz != b->nz)) { 1850f5bd95cSBarry Smith flag = PETSC_FALSE; 186b97cf49bSBarry Smith } 187b97cf49bSBarry Smith 188b97cf49bSBarry Smith /* if the a->i are the same */ 189*273d9f13SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr); 190b97cf49bSBarry Smith 191b97cf49bSBarry Smith /* if a->j are the same */ 1920f5bd95cSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr); 193b97cf49bSBarry Smith 194ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 1953a40ed3dSBarry Smith PetscFunctionReturn(0); 196b97cf49bSBarry Smith } 197b97cf49bSBarry Smith 198b97cf49bSBarry Smith 199b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2000b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2013eda8832SBarry Smith MatGetRow_MPIAdj, 2023eda8832SBarry Smith MatRestoreRow_MPIAdj, 203b97cf49bSBarry Smith 0, 204b97cf49bSBarry Smith 0, 205b97cf49bSBarry Smith 0, 206b97cf49bSBarry Smith 0, 2070b3b1487SBarry Smith 0, 2080b3b1487SBarry Smith 0, 2090b3b1487SBarry Smith 0, 2100b3b1487SBarry Smith 0, 2110b3b1487SBarry Smith 0, 2120b3b1487SBarry Smith 0, 2130b3b1487SBarry Smith 0, 2140b3b1487SBarry Smith 0, 2150b3b1487SBarry Smith 0, 2163eda8832SBarry Smith MatEqual_MPIAdj, 2170b3b1487SBarry Smith 0, 2180b3b1487SBarry Smith 0, 2190b3b1487SBarry Smith 0, 2200b3b1487SBarry Smith 0, 2210b3b1487SBarry Smith 0, 2220b3b1487SBarry Smith 0, 2233eda8832SBarry Smith MatSetOption_MPIAdj, 2240b3b1487SBarry Smith 0, 2250b3b1487SBarry Smith 0, 2260b3b1487SBarry Smith 0, 2270b3b1487SBarry Smith 0, 2280b3b1487SBarry Smith 0, 2290b3b1487SBarry Smith 0, 230*273d9f13SBarry Smith 0, 231*273d9f13SBarry Smith 0, 2323eda8832SBarry Smith MatGetOwnershipRange_MPIAdj, 2330b3b1487SBarry Smith 0, 2340b3b1487SBarry Smith 0, 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, 251b97cf49bSBarry Smith 0, 2523eda8832SBarry Smith MatGetBlockSize_MPIAdj, 253b97cf49bSBarry Smith 0, 254b97cf49bSBarry Smith 0, 255b97cf49bSBarry Smith 0, 256b97cf49bSBarry Smith 0, 257b97cf49bSBarry Smith 0, 2580752156aSBarry Smith 0, 2590752156aSBarry Smith 0, 2600b3b1487SBarry Smith 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 263b9b97703SBarry Smith MatDestroy_MPIAdj, 264b9b97703SBarry Smith MatView_MPIAdj, 2650b3b1487SBarry Smith MatGetMaps_Petsc}; 266b97cf49bSBarry Smith 267b97cf49bSBarry Smith 268*273d9f13SBarry Smith EXTERN_C_BEGIN 269*273d9f13SBarry Smith #undef __FUNC__ 270*273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_MPIAdj" 271*273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B) 272*273d9f13SBarry Smith { 273*273d9f13SBarry Smith Mat_MPIAdj *b; 274*273d9f13SBarry Smith int ii,ierr,size,rank; 275*273d9f13SBarry Smith 276*273d9f13SBarry Smith PetscFunctionBegin; 277*273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 278*273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 279*273d9f13SBarry Smith 280*273d9f13SBarry Smith B->data = (void*)(b = PetscNew(Mat_MPIAdj));CHKPTRQ(b); 281*273d9f13SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 282*273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 283*273d9f13SBarry Smith B->factor = 0; 284*273d9f13SBarry Smith B->lupivotthreshold = 1.0; 285*273d9f13SBarry Smith B->mapping = 0; 286*273d9f13SBarry Smith B->assembled = PETSC_FALSE; 287*273d9f13SBarry Smith 288*273d9f13SBarry Smith ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr); 289*273d9f13SBarry Smith B->n = B->N; 290*273d9f13SBarry Smith 291*273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 292*273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 293*273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 294*273d9f13SBarry Smith /* we don't know the "local columns" so just use the row information :-(*/ 295*273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr); 296*273d9f13SBarry Smith 297*273d9f13SBarry Smith b->rowners = (int*)PetscMalloc((size+1)*sizeof(int));CHKPTRQ(b->rowners); 298*273d9f13SBarry Smith PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 299*273d9f13SBarry Smith ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr); 300*273d9f13SBarry Smith b->rowners[0] = 0; 301*273d9f13SBarry Smith for (ii=2; ii<=size; ii++) { 302*273d9f13SBarry Smith b->rowners[ii] += b->rowners[ii-1]; 303*273d9f13SBarry Smith } 304*273d9f13SBarry Smith b->rstart = b->rowners[rank]; 305*273d9f13SBarry Smith b->rend = b->rowners[rank+1]; 306*273d9f13SBarry Smith 307*273d9f13SBarry Smith PetscFunctionReturn(0); 308*273d9f13SBarry Smith } 309*273d9f13SBarry Smith EXTERN_C_END 310*273d9f13SBarry Smith 311*273d9f13SBarry Smith #undef __FUNC__ 312*273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMPIAdjSetPreallocation" 313*273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values) 314*273d9f13SBarry Smith { 315*273d9f13SBarry Smith Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 316*273d9f13SBarry Smith int ierr; 317*273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g) 318*273d9f13SBarry Smith int ii; 319*273d9f13SBarry Smith #endif 320*273d9f13SBarry Smith 321*273d9f13SBarry Smith PetscFunctionBegin; 322*273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 323*273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g) 324*273d9f13SBarry Smith if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]); 325*273d9f13SBarry Smith for (ii=1; ii<B->m; ii++) { 326*273d9f13SBarry Smith if (i[ii] < 0 || i[ii] < i[ii-1]) { 327*273d9f13SBarry Smith SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]); 328*273d9f13SBarry Smith } 329*273d9f13SBarry Smith } 330*273d9f13SBarry Smith for (ii=0; ii<i[B->m]; ii++) { 331*273d9f13SBarry Smith if (j[ii] < 0 || j[ii] >= B->N) { 332*273d9f13SBarry Smith SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]); 333*273d9f13SBarry Smith } 334*273d9f13SBarry Smith } 335*273d9f13SBarry Smith #endif 336*273d9f13SBarry Smith 337*273d9f13SBarry Smith b->j = j; 338*273d9f13SBarry Smith b->i = i; 339*273d9f13SBarry Smith b->values = values; 340*273d9f13SBarry Smith 341*273d9f13SBarry Smith b->nz = i[B->m]; 342*273d9f13SBarry Smith b->diag = 0; 343*273d9f13SBarry Smith b->symmetric = PETSC_FALSE; 344*273d9f13SBarry Smith b->freeaij = PETSC_TRUE; 345*273d9f13SBarry Smith 346*273d9f13SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 347*273d9f13SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 348*273d9f13SBarry Smith PetscFunctionReturn(0); 349*273d9f13SBarry Smith } 350*273d9f13SBarry Smith 351b97cf49bSBarry Smith #undef __FUNC__ 352b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIAdj" 353b97cf49bSBarry Smith /*@C 3543eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 355b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 356b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 357b97cf49bSBarry Smith 358ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 359ef5ee4d1SLois Curfman McInnes 360b97cf49bSBarry Smith Input Parameters: 361c2e958feSBarry Smith + comm - MPI communicator 3620752156aSBarry Smith . m - number of local rows 363b97cf49bSBarry Smith . n - number of columns 364b97cf49bSBarry Smith . i - the indices into j for the start of each row 365329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 366ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 367329f5518SBarry Smith - values -[optional] edge weights 368b97cf49bSBarry Smith 369b97cf49bSBarry Smith Output Parameter: 370b97cf49bSBarry Smith . A - the matrix 371b97cf49bSBarry Smith 37215091d37SBarry Smith Level: intermediate 37315091d37SBarry Smith 3744bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 3754bc6d8c0SBarry Smith MatSetValues(). 376c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 3771198db28SBarry Smith when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you 3781198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 379327e1a73SBarry Smith Should not include the matrix diagonals. 380b97cf49bSBarry Smith 381ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 382b97cf49bSBarry Smith 38391e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering() 384b97cf49bSBarry Smith @*/ 3853eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A) 386b97cf49bSBarry Smith { 387*273d9f13SBarry Smith int ierr; 388b97cf49bSBarry Smith 389433994e6SBarry Smith PetscFunctionBegin; 390*273d9f13SBarry Smith ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr); 391*273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 392*273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 3933a40ed3dSBarry Smith PetscFunctionReturn(0); 394b97cf49bSBarry Smith } 395b97cf49bSBarry Smith 396*273d9f13SBarry Smith EXTERN_C_BEGIN 397c2e958feSBarry Smith #undef __FUNC__ 398*273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatConvertTo_MPIAdj" 399*273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B) 400c2e958feSBarry Smith { 4014c49b128SBarry Smith int i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a; 402c2e958feSBarry Smith Scalar *ra; 403c2e958feSBarry Smith MPI_Comm comm; 404c2e958feSBarry Smith 405c2e958feSBarry Smith PetscFunctionBegin; 406c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 407c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 408c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 409c2e958feSBarry Smith 410c2e958feSBarry Smith /* count the number of nonzeros per row */ 411c2e958feSBarry Smith for (i=0; i<m; i++) { 412c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 413c2e958feSBarry Smith for (j=0; j<len; j++) { 414c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 415c2e958feSBarry Smith } 416c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 417c2e958feSBarry Smith nzeros += len; 418c2e958feSBarry Smith } 419c2e958feSBarry Smith 420c2e958feSBarry Smith /* malloc space for nonzeros */ 421c2e958feSBarry Smith a = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(a); 422c2e958feSBarry Smith ia = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(ia); 423c2e958feSBarry Smith ja = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(ja); 424c2e958feSBarry Smith 425c2e958feSBarry Smith nzeros = 0; 426c2e958feSBarry Smith ia[0] = 0; 427c2e958feSBarry Smith for (i=0; i<m; i++) { 428c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 429c2e958feSBarry Smith cnt = 0; 430c2e958feSBarry Smith for (j=0; j<len; j++) { 431c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 4324c49b128SBarry Smith a[nzeros+cnt] = (int) PetscAbsScalar(ra[j]); 433c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 434c2e958feSBarry Smith } 435c2e958feSBarry Smith } 436c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 437c2e958feSBarry Smith nzeros += cnt; 438c2e958feSBarry Smith ia[i+1] = nzeros; 439c2e958feSBarry Smith } 440c2e958feSBarry Smith 441c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 4423eda8832SBarry Smith ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr); 443c2e958feSBarry Smith 444c2e958feSBarry Smith PetscFunctionReturn(0); 445c2e958feSBarry Smith } 446*273d9f13SBarry Smith EXTERN_C_END 447433994e6SBarry Smith 448433994e6SBarry Smith 449433994e6SBarry Smith 450433994e6SBarry Smith 451433994e6SBarry Smith 452