1*6cd91f64SBarry Smith /*$Id: mpiadj.c,v 1.50 2000/10/30 03:28:19 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); 22ffac6cdbSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 23ffac6cdbSBarry 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 200*6cd91f64SBarry Smith #undef __FUNC__ 201*6cd91f64SBarry Smith #define __FUNC__ /*<a name="MatGetRowIJ_MPIAdj"></a>*/"MatGetRowIJ_MPIAdj" 202*6cd91f64SBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done) 203*6cd91f64SBarry Smith { 204*6cd91f64SBarry Smith int ierr,size; 205*6cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 206*6cd91f64SBarry Smith 207*6cd91f64SBarry Smith PetscFunctionBegin; 208*6cd91f64SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 209*6cd91f64SBarry Smith if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 210*6cd91f64SBarry Smith if (oshift) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 211*6cd91f64SBarry Smith *m = A->m; 212*6cd91f64SBarry Smith *ia = a->i; 213*6cd91f64SBarry Smith *ja = a->j; 214*6cd91f64SBarry Smith *done = PETSC_TRUE; 215*6cd91f64SBarry Smith PetscFunctionReturn(0); 216*6cd91f64SBarry Smith } 217b97cf49bSBarry Smith 218b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2190b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2203eda8832SBarry Smith MatGetRow_MPIAdj, 2213eda8832SBarry Smith MatRestoreRow_MPIAdj, 222b97cf49bSBarry Smith 0, 223b97cf49bSBarry Smith 0, 224b97cf49bSBarry Smith 0, 225b97cf49bSBarry Smith 0, 2260b3b1487SBarry Smith 0, 2270b3b1487SBarry Smith 0, 2280b3b1487SBarry Smith 0, 2290b3b1487SBarry Smith 0, 2300b3b1487SBarry Smith 0, 2310b3b1487SBarry Smith 0, 2320b3b1487SBarry Smith 0, 2330b3b1487SBarry Smith 0, 2340b3b1487SBarry Smith 0, 2353eda8832SBarry Smith MatEqual_MPIAdj, 2360b3b1487SBarry Smith 0, 2370b3b1487SBarry Smith 0, 2380b3b1487SBarry Smith 0, 2390b3b1487SBarry Smith 0, 2400b3b1487SBarry Smith 0, 2410b3b1487SBarry Smith 0, 2423eda8832SBarry Smith MatSetOption_MPIAdj, 2430b3b1487SBarry Smith 0, 2440b3b1487SBarry Smith 0, 2450b3b1487SBarry Smith 0, 2460b3b1487SBarry Smith 0, 2470b3b1487SBarry Smith 0, 2480b3b1487SBarry Smith 0, 249273d9f13SBarry Smith 0, 250273d9f13SBarry Smith 0, 2513eda8832SBarry Smith MatGetOwnershipRange_MPIAdj, 2520b3b1487SBarry Smith 0, 2530b3b1487SBarry Smith 0, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 2580b3b1487SBarry Smith 0, 2590b3b1487SBarry Smith 0, 2600b3b1487SBarry Smith 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 2660b3b1487SBarry Smith 0, 2670b3b1487SBarry Smith 0, 2680b3b1487SBarry Smith 0, 2690b3b1487SBarry Smith 0, 270b97cf49bSBarry Smith 0, 2713eda8832SBarry Smith MatGetBlockSize_MPIAdj, 272*6cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 273b97cf49bSBarry Smith 0, 274b97cf49bSBarry Smith 0, 275b97cf49bSBarry Smith 0, 276b97cf49bSBarry Smith 0, 2770752156aSBarry Smith 0, 2780752156aSBarry Smith 0, 2790b3b1487SBarry Smith 0, 2800b3b1487SBarry Smith 0, 2810b3b1487SBarry Smith 0, 282b9b97703SBarry Smith MatDestroy_MPIAdj, 283b9b97703SBarry Smith MatView_MPIAdj, 2840b3b1487SBarry Smith MatGetMaps_Petsc}; 285b97cf49bSBarry Smith 286b97cf49bSBarry Smith 287273d9f13SBarry Smith EXTERN_C_BEGIN 288273d9f13SBarry Smith #undef __FUNC__ 289273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_MPIAdj" 290273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B) 291273d9f13SBarry Smith { 292273d9f13SBarry Smith Mat_MPIAdj *b; 293273d9f13SBarry Smith int ii,ierr,size,rank; 294273d9f13SBarry Smith 295273d9f13SBarry Smith PetscFunctionBegin; 296273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 297273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 298273d9f13SBarry Smith 299273d9f13SBarry Smith B->data = (void*)(b = PetscNew(Mat_MPIAdj));CHKPTRQ(b); 300273d9f13SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 301273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 302273d9f13SBarry Smith B->factor = 0; 303273d9f13SBarry Smith B->lupivotthreshold = 1.0; 304273d9f13SBarry Smith B->mapping = 0; 305273d9f13SBarry Smith B->assembled = PETSC_FALSE; 306273d9f13SBarry Smith 307273d9f13SBarry Smith ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr); 308273d9f13SBarry Smith B->n = B->N; 309273d9f13SBarry Smith 310273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 311273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 312273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 313273d9f13SBarry Smith /* we don't know the "local columns" so just use the row information :-(*/ 314273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr); 315273d9f13SBarry Smith 316273d9f13SBarry Smith b->rowners = (int*)PetscMalloc((size+1)*sizeof(int));CHKPTRQ(b->rowners); 317273d9f13SBarry Smith PLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 318273d9f13SBarry Smith ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr); 319273d9f13SBarry Smith b->rowners[0] = 0; 320273d9f13SBarry Smith for (ii=2; ii<=size; ii++) { 321273d9f13SBarry Smith b->rowners[ii] += b->rowners[ii-1]; 322273d9f13SBarry Smith } 323273d9f13SBarry Smith b->rstart = b->rowners[rank]; 324273d9f13SBarry Smith b->rend = b->rowners[rank+1]; 325273d9f13SBarry Smith 326273d9f13SBarry Smith PetscFunctionReturn(0); 327273d9f13SBarry Smith } 328273d9f13SBarry Smith EXTERN_C_END 329273d9f13SBarry Smith 330273d9f13SBarry Smith #undef __FUNC__ 331273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMPIAdjSetPreallocation" 332273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values) 333273d9f13SBarry Smith { 334273d9f13SBarry Smith Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 335273d9f13SBarry Smith int ierr; 336273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g) 337273d9f13SBarry Smith int ii; 338273d9f13SBarry Smith #endif 339273d9f13SBarry Smith 340273d9f13SBarry Smith PetscFunctionBegin; 341273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 342273d9f13SBarry Smith #if defined(PETSC_USE_BOPT_g) 343273d9f13SBarry Smith if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]); 344273d9f13SBarry Smith for (ii=1; ii<B->m; ii++) { 345273d9f13SBarry Smith if (i[ii] < 0 || i[ii] < i[ii-1]) { 346273d9f13SBarry Smith SETERRQ4(1,"i[%d] index is out of range: i[%d]",ii,i[ii],ii-1,i[ii-1]); 347273d9f13SBarry Smith } 348273d9f13SBarry Smith } 349273d9f13SBarry Smith for (ii=0; ii<i[B->m]; ii++) { 350273d9f13SBarry Smith if (j[ii] < 0 || j[ii] >= B->N) { 351273d9f13SBarry Smith SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]); 352273d9f13SBarry Smith } 353273d9f13SBarry Smith } 354273d9f13SBarry Smith #endif 355273d9f13SBarry Smith 356273d9f13SBarry Smith b->j = j; 357273d9f13SBarry Smith b->i = i; 358273d9f13SBarry Smith b->values = values; 359273d9f13SBarry Smith 360273d9f13SBarry Smith b->nz = i[B->m]; 361273d9f13SBarry Smith b->diag = 0; 362273d9f13SBarry Smith b->symmetric = PETSC_FALSE; 363273d9f13SBarry Smith b->freeaij = PETSC_TRUE; 364273d9f13SBarry Smith 365273d9f13SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 366273d9f13SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 367273d9f13SBarry Smith PetscFunctionReturn(0); 368273d9f13SBarry Smith } 369273d9f13SBarry Smith 370b97cf49bSBarry Smith #undef __FUNC__ 371b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateMPIAdj" 372b97cf49bSBarry Smith /*@C 3733eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 374b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 375b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 376b97cf49bSBarry Smith 377ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 378ef5ee4d1SLois Curfman McInnes 379b97cf49bSBarry Smith Input Parameters: 380c2e958feSBarry Smith + comm - MPI communicator 3810752156aSBarry Smith . m - number of local rows 382b97cf49bSBarry Smith . n - number of columns 383b97cf49bSBarry Smith . i - the indices into j for the start of each row 384329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 385ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 386329f5518SBarry Smith - values -[optional] edge weights 387b97cf49bSBarry Smith 388b97cf49bSBarry Smith Output Parameter: 389b97cf49bSBarry Smith . A - the matrix 390b97cf49bSBarry Smith 39115091d37SBarry Smith Level: intermediate 39215091d37SBarry Smith 3934bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 3944bc6d8c0SBarry Smith MatSetValues(). 395c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 3961198db28SBarry Smith when the matrix is destroyed. And you must allocate them with PetscMalloc(). If you 3971198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 398327e1a73SBarry Smith Should not include the matrix diagonals. 399b97cf49bSBarry Smith 400ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 401b97cf49bSBarry Smith 40291e9ee9fSBarry Smith .seealso: MatCreate(), MatCreateSeqAdj(), MatGetOrdering() 403b97cf49bSBarry Smith @*/ 4043eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A) 405b97cf49bSBarry Smith { 406273d9f13SBarry Smith int ierr; 407b97cf49bSBarry Smith 408433994e6SBarry Smith PetscFunctionBegin; 409273d9f13SBarry Smith ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr); 410273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 411273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 4123a40ed3dSBarry Smith PetscFunctionReturn(0); 413b97cf49bSBarry Smith } 414b97cf49bSBarry Smith 415273d9f13SBarry Smith EXTERN_C_BEGIN 416c2e958feSBarry Smith #undef __FUNC__ 417273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatConvertTo_MPIAdj" 418273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B) 419c2e958feSBarry Smith { 4204c49b128SBarry Smith int i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a; 421c2e958feSBarry Smith Scalar *ra; 422c2e958feSBarry Smith MPI_Comm comm; 423c2e958feSBarry Smith 424c2e958feSBarry Smith PetscFunctionBegin; 425c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 426c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 427c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 428c2e958feSBarry Smith 429c2e958feSBarry Smith /* count the number of nonzeros per row */ 430c2e958feSBarry Smith for (i=0; i<m; i++) { 431c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 432c2e958feSBarry Smith for (j=0; j<len; j++) { 433c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 434c2e958feSBarry Smith } 435c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 436c2e958feSBarry Smith nzeros += len; 437c2e958feSBarry Smith } 438c2e958feSBarry Smith 439c2e958feSBarry Smith /* malloc space for nonzeros */ 440c2e958feSBarry Smith a = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(a); 441c2e958feSBarry Smith ia = (int*)PetscMalloc((N+1)*sizeof(int));CHKPTRQ(ia); 442c2e958feSBarry Smith ja = (int*)PetscMalloc((nzeros+1)*sizeof(int));CHKPTRQ(ja); 443c2e958feSBarry Smith 444c2e958feSBarry Smith nzeros = 0; 445c2e958feSBarry Smith ia[0] = 0; 446c2e958feSBarry Smith for (i=0; i<m; i++) { 447c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 448c2e958feSBarry Smith cnt = 0; 449c2e958feSBarry Smith for (j=0; j<len; j++) { 450c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 4514c49b128SBarry Smith a[nzeros+cnt] = (int) PetscAbsScalar(ra[j]); 452c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 453c2e958feSBarry Smith } 454c2e958feSBarry Smith } 455c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 456c2e958feSBarry Smith nzeros += cnt; 457c2e958feSBarry Smith ia[i+1] = nzeros; 458c2e958feSBarry Smith } 459c2e958feSBarry Smith 460c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 4613eda8832SBarry Smith ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr); 462c2e958feSBarry Smith 463c2e958feSBarry Smith PetscFunctionReturn(0); 464c2e958feSBarry Smith } 465273d9f13SBarry Smith EXTERN_C_END 466433994e6SBarry Smith 467433994e6SBarry Smith 468433994e6SBarry Smith 469433994e6SBarry Smith 470433994e6SBarry Smith 471