173f4d377SMatthew Knepley /*$Id: mpiadj.c,v 1.66 2001/08/07 03:02:59 balay Exp $*/ 2b97cf49bSBarry Smith 3b97cf49bSBarry Smith /* 4b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 5b97cf49bSBarry Smith */ 63eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h" 7325e03aeSBarry Smith #include "petscsys.h" 8b97cf49bSBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII" 11b0a32e0cSBarry Smith int MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 12b97cf49bSBarry Smith { 133eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 14fb9695e5SSatish Balay int ierr,i,j,m = A->m; 15fb9695e5SSatish Balay char *name; 16ce1411ecSBarry Smith PetscViewerFormat format; 17b97cf49bSBarry Smith 18433994e6SBarry Smith PetscFunctionBegin; 193a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 20b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 21fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 223a40ed3dSBarry Smith PetscFunctionReturn(0); 23fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 24ffac6cdbSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 250752156aSBarry Smith } else { 26b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 27b97cf49bSBarry Smith for (i=0; i<m; i++) { 28b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr); 29b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 30b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr); 310752156aSBarry Smith } 32b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 33b97cf49bSBarry Smith } 34b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 35b97cf49bSBarry Smith } 36b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 373a40ed3dSBarry Smith PetscFunctionReturn(0); 38b97cf49bSBarry Smith } 39b97cf49bSBarry Smith 404a2ae208SSatish Balay #undef __FUNCT__ 414a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj" 42b0a32e0cSBarry Smith int MatView_MPIAdj(Mat A,PetscViewer viewer) 43b97cf49bSBarry Smith { 44b97cf49bSBarry Smith int ierr; 456831982aSBarry Smith PetscTruth isascii; 46b97cf49bSBarry Smith 47433994e6SBarry Smith PetscFunctionBegin; 48b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 490f5bd95cSBarry Smith if (isascii) { 503eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 515cd90555SBarry Smith } else { 5229bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 53b97cf49bSBarry Smith } 543a40ed3dSBarry Smith PetscFunctionReturn(0); 55b97cf49bSBarry Smith } 56b97cf49bSBarry Smith 574a2ae208SSatish Balay #undef __FUNCT__ 584a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 593eda8832SBarry Smith int MatDestroy_MPIAdj(Mat mat) 60b97cf49bSBarry Smith { 613eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 6294d884c6SBarry Smith int ierr; 6394d884c6SBarry Smith 6494d884c6SBarry Smith PetscFunctionBegin; 65aa482453SBarry Smith #if defined(PETSC_USE_LOG) 66b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz); 67b97cf49bSBarry Smith #endif 68606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 690400557aSBarry Smith if (a->freeaij) { 70606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 71606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 721198db28SBarry Smith if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);} 730400557aSBarry Smith } 740400557aSBarry Smith ierr = PetscFree(a->rowners);CHKERRQ(ierr); 751198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 763a40ed3dSBarry Smith PetscFunctionReturn(0); 77b97cf49bSBarry Smith } 78b97cf49bSBarry Smith 794a2ae208SSatish Balay #undef __FUNCT__ 804a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 813eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 82b97cf49bSBarry Smith { 833eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 84b97cf49bSBarry Smith 85433994e6SBarry Smith PetscFunctionBegin; 8612c028f9SKris Buschelman switch (op) { 8777e54ba9SKris Buschelman case MAT_SYMMETRIC: 8812c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 899a4540c5SBarry Smith case MAT_HERMITIAN: 90b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 9112c028f9SKris Buschelman break; 929a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 939a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 949a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 959a4540c5SBarry Smith a->symmetric = PETSC_FALSE; 969a4540c5SBarry Smith break; 979a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 989a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 999a4540c5SBarry Smith break; 10012c028f9SKris Buschelman default: 101b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 10212c028f9SKris Buschelman break; 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 1124a2ae208SSatish Balay #undef __FUNCT__ 1134a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 1143eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A) 115b97cf49bSBarry Smith { 1163eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 11782502324SSatish Balay int i,j,*diag,m = A->m,ierr; 118b97cf49bSBarry Smith 119433994e6SBarry Smith PetscFunctionBegin; 120b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr); 121b0a32e0cSBarry Smith PetscLogObjectMemory(A,(m+1)*sizeof(int)); 122273d9f13SBarry 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 1344a2ae208SSatish Balay #undef __FUNCT__ 1354a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 13687828ca2SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v) 137b97cf49bSBarry Smith { 1383eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 139b97cf49bSBarry Smith int *itmp; 140b97cf49bSBarry Smith 141433994e6SBarry Smith PetscFunctionBegin; 1420752156aSBarry Smith row -= a->rstart; 1430752156aSBarry Smith 144273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 145b97cf49bSBarry Smith 146b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 147b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 148b97cf49bSBarry Smith if (idx) { 149b97cf49bSBarry Smith itmp = a->j + a->i[row]; 150b97cf49bSBarry Smith if (*nz) { 151b97cf49bSBarry Smith *idx = itmp; 152b97cf49bSBarry Smith } 153b97cf49bSBarry Smith else *idx = 0; 154b97cf49bSBarry Smith } 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 156b97cf49bSBarry Smith } 157b97cf49bSBarry Smith 1584a2ae208SSatish Balay #undef __FUNCT__ 1594a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 16087828ca2SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v) 161b97cf49bSBarry Smith { 162433994e6SBarry Smith PetscFunctionBegin; 1633a40ed3dSBarry Smith PetscFunctionReturn(0); 164b97cf49bSBarry Smith } 165b97cf49bSBarry Smith 1664a2ae208SSatish Balay #undef __FUNCT__ 1674a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj" 1683eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs) 169b97cf49bSBarry Smith { 170433994e6SBarry Smith PetscFunctionBegin; 171b97cf49bSBarry Smith *bs = 1; 1723a40ed3dSBarry Smith PetscFunctionReturn(0); 173b97cf49bSBarry Smith } 174b97cf49bSBarry Smith 175b97cf49bSBarry Smith 1764a2ae208SSatish Balay #undef __FUNCT__ 1774a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 1783eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 179b97cf49bSBarry Smith { 1803eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 1810f5bd95cSBarry Smith int ierr; 1820f5bd95cSBarry Smith PetscTruth flag; 183b97cf49bSBarry Smith 184433994e6SBarry Smith PetscFunctionBegin; 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 2004a2ae208SSatish Balay #undef __FUNCT__ 2014a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 2024e7234bfSBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done) 2036cd91f64SBarry Smith { 204d42735a1SBarry Smith int ierr,size,i; 2056cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 2066cd91f64SBarry Smith 2076cd91f64SBarry Smith PetscFunctionBegin; 2086cd91f64SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 2096cd91f64SBarry Smith if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 2106cd91f64SBarry Smith *m = A->m; 2116cd91f64SBarry Smith *ia = a->i; 2126cd91f64SBarry Smith *ja = a->j; 2136cd91f64SBarry Smith *done = PETSC_TRUE; 214d42735a1SBarry Smith if (oshift) { 215d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 216d42735a1SBarry Smith (*ja)[i]++; 217d42735a1SBarry Smith } 218d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 219d42735a1SBarry Smith } 220d42735a1SBarry Smith PetscFunctionReturn(0); 221d42735a1SBarry Smith } 222d42735a1SBarry Smith 2234a2ae208SSatish Balay #undef __FUNCT__ 2244a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 2254e7234bfSBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done) 226d42735a1SBarry Smith { 227d42735a1SBarry Smith int i; 228d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 229d42735a1SBarry Smith 230d42735a1SBarry Smith PetscFunctionBegin; 231c5b3d447SBarry Smith if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()"); 232c5b3d447SBarry Smith if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()"); 233d42735a1SBarry Smith if (oshift) { 234d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 235d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 236d42735a1SBarry Smith (*ja)[i]--; 237d42735a1SBarry Smith } 238d42735a1SBarry Smith } 2396cd91f64SBarry Smith PetscFunctionReturn(0); 2406cd91f64SBarry Smith } 241b97cf49bSBarry Smith 242b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2430b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2443eda8832SBarry Smith MatGetRow_MPIAdj, 2453eda8832SBarry Smith MatRestoreRow_MPIAdj, 246b97cf49bSBarry Smith 0, 24797304618SKris Buschelman /* 4*/ 0, 248b97cf49bSBarry Smith 0, 249b97cf49bSBarry Smith 0, 250b97cf49bSBarry Smith 0, 2510b3b1487SBarry Smith 0, 2520b3b1487SBarry Smith 0, 25397304618SKris Buschelman /*10*/ 0, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 25897304618SKris Buschelman /*15*/ 0, 2593eda8832SBarry Smith MatEqual_MPIAdj, 2600b3b1487SBarry Smith 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 26397304618SKris Buschelman /*20*/ 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 2663eda8832SBarry Smith MatSetOption_MPIAdj, 2670b3b1487SBarry Smith 0, 26897304618SKris Buschelman /*25*/ 0, 2690b3b1487SBarry Smith 0, 2700b3b1487SBarry Smith 0, 2710b3b1487SBarry Smith 0, 2720b3b1487SBarry Smith 0, 27397304618SKris Buschelman /*30*/ 0, 2740b3b1487SBarry Smith 0, 275273d9f13SBarry Smith 0, 276273d9f13SBarry Smith 0, 2770b3b1487SBarry Smith 0, 27897304618SKris Buschelman /*35*/ 0, 2790b3b1487SBarry Smith 0, 2800b3b1487SBarry Smith 0, 2810b3b1487SBarry Smith 0, 2820b3b1487SBarry Smith 0, 28397304618SKris Buschelman /*40*/ 0, 2840b3b1487SBarry Smith 0, 2850b3b1487SBarry Smith 0, 2860b3b1487SBarry Smith 0, 2870b3b1487SBarry Smith 0, 28897304618SKris Buschelman /*45*/ 0, 2890b3b1487SBarry Smith 0, 2900b3b1487SBarry Smith 0, 2910b3b1487SBarry Smith 0, 2920b3b1487SBarry Smith 0, 29397304618SKris Buschelman /*50*/ MatGetBlockSize_MPIAdj, 2946cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 295d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 296b97cf49bSBarry Smith 0, 297b97cf49bSBarry Smith 0, 29897304618SKris Buschelman /*55*/ 0, 299b97cf49bSBarry Smith 0, 3000752156aSBarry Smith 0, 3010752156aSBarry Smith 0, 3020b3b1487SBarry Smith 0, 30397304618SKris Buschelman /*60*/ 0, 304b9b97703SBarry Smith MatDestroy_MPIAdj, 305b9b97703SBarry Smith MatView_MPIAdj, 30697304618SKris Buschelman MatGetPetscMaps_Petsc, 30797304618SKris Buschelman 0, 30897304618SKris Buschelman /*65*/ 0, 30997304618SKris Buschelman 0, 31097304618SKris Buschelman 0, 31197304618SKris Buschelman 0, 31297304618SKris Buschelman 0, 31397304618SKris Buschelman /*70*/ 0, 31497304618SKris Buschelman 0, 31597304618SKris Buschelman 0, 31697304618SKris Buschelman 0, 31797304618SKris Buschelman 0, 31897304618SKris Buschelman /*75*/ 0, 31997304618SKris Buschelman 0, 32097304618SKris Buschelman 0, 32197304618SKris Buschelman 0, 32297304618SKris Buschelman 0, 32397304618SKris Buschelman /*80*/ 0, 32497304618SKris Buschelman 0, 32597304618SKris Buschelman 0, 32697304618SKris Buschelman 0, 32797304618SKris Buschelman /*85*/ 0 32897304618SKris Buschelman }; 329b97cf49bSBarry Smith 330a23d5eceSKris Buschelman EXTERN_C_BEGIN 331a23d5eceSKris Buschelman #undef __FUNCT__ 332a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 333a23d5eceSKris Buschelman int MatMPIAdjSetPreallocation_MPIAdj(Mat B,int *i,int *j,int *values) 334a23d5eceSKris Buschelman { 335a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 336a23d5eceSKris Buschelman int ierr; 337a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g) 338a23d5eceSKris Buschelman int ii; 339a23d5eceSKris Buschelman #endif 340a23d5eceSKris Buschelman 341a23d5eceSKris Buschelman PetscFunctionBegin; 342a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 343a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g) 344a23d5eceSKris Buschelman if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]); 345a23d5eceSKris Buschelman for (ii=1; ii<B->m; ii++) { 346a23d5eceSKris Buschelman if (i[ii] < 0 || i[ii] < i[ii-1]) { 347a23d5eceSKris Buschelman SETERRQ4(1,"i[%d]=%d index is out of range: i[%d]=%d",ii,i[ii],ii-1,i[ii-1]); 348a23d5eceSKris Buschelman } 349a23d5eceSKris Buschelman } 350a23d5eceSKris Buschelman for (ii=0; ii<i[B->m]; ii++) { 351a23d5eceSKris Buschelman if (j[ii] < 0 || j[ii] >= B->N) { 352a23d5eceSKris Buschelman SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]); 353a23d5eceSKris Buschelman } 354a23d5eceSKris Buschelman } 355a23d5eceSKris Buschelman #endif 356a23d5eceSKris Buschelman 357a23d5eceSKris Buschelman b->j = j; 358a23d5eceSKris Buschelman b->i = i; 359a23d5eceSKris Buschelman b->values = values; 360a23d5eceSKris Buschelman 361a23d5eceSKris Buschelman b->nz = i[B->m]; 362a23d5eceSKris Buschelman b->diag = 0; 363a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 364a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 365a23d5eceSKris Buschelman 366a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 367a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 368a23d5eceSKris Buschelman PetscFunctionReturn(0); 369a23d5eceSKris Buschelman } 370a23d5eceSKris Buschelman EXTERN_C_END 371b97cf49bSBarry Smith 3720bad9183SKris Buschelman /*MC 373fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 3740bad9183SKris Buschelman intended for use constructing orderings and partitionings. 3750bad9183SKris Buschelman 3760bad9183SKris Buschelman Level: beginner 3770bad9183SKris Buschelman 3780bad9183SKris Buschelman .seealso: MatCreateMPIAdj 3790bad9183SKris Buschelman M*/ 3800bad9183SKris Buschelman 381273d9f13SBarry Smith EXTERN_C_BEGIN 3824a2ae208SSatish Balay #undef __FUNCT__ 3834a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 384273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B) 385273d9f13SBarry Smith { 386273d9f13SBarry Smith Mat_MPIAdj *b; 387273d9f13SBarry Smith int ii,ierr,size,rank; 388273d9f13SBarry Smith 389273d9f13SBarry Smith PetscFunctionBegin; 390273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 391273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 392273d9f13SBarry Smith 393b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr); 394b0a32e0cSBarry Smith B->data = (void*)b; 395273d9f13SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 396273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 397273d9f13SBarry Smith B->factor = 0; 398273d9f13SBarry Smith B->lupivotthreshold = 1.0; 399273d9f13SBarry Smith B->mapping = 0; 400273d9f13SBarry Smith B->assembled = PETSC_FALSE; 401273d9f13SBarry Smith 402273d9f13SBarry Smith ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr); 403273d9f13SBarry Smith B->n = B->N; 404273d9f13SBarry Smith 405273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 406273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 4078a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 408273d9f13SBarry Smith /* we don't know the "local columns" so just use the row information :-(*/ 4098a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr); 410273d9f13SBarry Smith 411b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr); 412b0a32e0cSBarry Smith PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 413273d9f13SBarry Smith ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr); 414273d9f13SBarry Smith b->rowners[0] = 0; 415273d9f13SBarry Smith for (ii=2; ii<=size; ii++) { 416273d9f13SBarry Smith b->rowners[ii] += b->rowners[ii-1]; 417273d9f13SBarry Smith } 418273d9f13SBarry Smith b->rstart = b->rowners[rank]; 419273d9f13SBarry Smith b->rend = b->rowners[rank+1]; 420a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C", 421a23d5eceSKris Buschelman "MatMPIAdjSetPreallocation_MPIAdj", 422a23d5eceSKris Buschelman MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 423273d9f13SBarry Smith PetscFunctionReturn(0); 424273d9f13SBarry Smith } 425273d9f13SBarry Smith EXTERN_C_END 426273d9f13SBarry Smith 4274a2ae208SSatish Balay #undef __FUNCT__ 4284a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 429a23d5eceSKris Buschelman /*@C 430a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 431a23d5eceSKris Buschelman 432a23d5eceSKris Buschelman Collective on MPI_Comm 433a23d5eceSKris Buschelman 434a23d5eceSKris Buschelman Input Parameters: 435a23d5eceSKris Buschelman + A - the matrix 436a23d5eceSKris Buschelman . i - the indices into j for the start of each row 437a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 438a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 439a23d5eceSKris Buschelman - values - [optional] edge weights 440a23d5eceSKris Buschelman 441a23d5eceSKris Buschelman Level: intermediate 442a23d5eceSKris Buschelman 443a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 444a23d5eceSKris Buschelman @*/ 445273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values) 446273d9f13SBarry Smith { 447a23d5eceSKris Buschelman int ierr,(*f)(Mat,int*,int*,int*); 448273d9f13SBarry Smith 449273d9f13SBarry Smith PetscFunctionBegin; 450a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 451a23d5eceSKris Buschelman if (f) { 452a23d5eceSKris Buschelman ierr = (*f)(B,i,j,values);CHKERRQ(ierr); 453273d9f13SBarry Smith } 454273d9f13SBarry Smith PetscFunctionReturn(0); 455273d9f13SBarry Smith } 456273d9f13SBarry Smith 4574a2ae208SSatish Balay #undef __FUNCT__ 4584a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 459b97cf49bSBarry Smith /*@C 4603eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 461b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 462b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 463b97cf49bSBarry Smith 464ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 465ef5ee4d1SLois Curfman McInnes 466b97cf49bSBarry Smith Input Parameters: 467c2e958feSBarry Smith + comm - MPI communicator 4680752156aSBarry Smith . m - number of local rows 469b97cf49bSBarry Smith . n - number of columns 470b97cf49bSBarry Smith . i - the indices into j for the start of each row 471329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 472ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 473329f5518SBarry Smith - values -[optional] edge weights 474b97cf49bSBarry Smith 475b97cf49bSBarry Smith Output Parameter: 476b97cf49bSBarry Smith . A - the matrix 477b97cf49bSBarry Smith 47815091d37SBarry Smith Level: intermediate 47915091d37SBarry Smith 4804bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 4814bc6d8c0SBarry Smith MatSetValues(). 482c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 483ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 4841198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 485327e1a73SBarry Smith Should not include the matrix diagonals. 486b97cf49bSBarry Smith 487e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 488ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 489ededeb1aSvictorle 490ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 491b97cf49bSBarry Smith 492e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 493b97cf49bSBarry Smith @*/ 4943eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A) 495b97cf49bSBarry Smith { 496273d9f13SBarry Smith int ierr; 497b97cf49bSBarry Smith 498433994e6SBarry Smith PetscFunctionBegin; 499273d9f13SBarry Smith ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr); 500273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 501273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 5023a40ed3dSBarry Smith PetscFunctionReturn(0); 503b97cf49bSBarry Smith } 504b97cf49bSBarry Smith 505273d9f13SBarry Smith EXTERN_C_BEGIN 5064a2ae208SSatish Balay #undef __FUNCT__ 5074a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj" 508676c34cdSKris Buschelman int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *newmat) { 509676c34cdSKris Buschelman Mat B; 5104c49b128SBarry Smith int i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a; 51187828ca2SBarry Smith PetscScalar *ra; 512c2e958feSBarry Smith MPI_Comm comm; 513c2e958feSBarry Smith 514c2e958feSBarry Smith PetscFunctionBegin; 515c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 516c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 517c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 518c2e958feSBarry Smith 519c2e958feSBarry Smith /* count the number of nonzeros per row */ 520c2e958feSBarry Smith for (i=0; i<m; i++) { 521c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 522c2e958feSBarry Smith for (j=0; j<len; j++) { 523c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 524c2e958feSBarry Smith } 525c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 526c2e958feSBarry Smith nzeros += len; 527c2e958feSBarry Smith } 528c2e958feSBarry Smith 529c2e958feSBarry Smith /* malloc space for nonzeros */ 530b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr); 531b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr); 532b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr); 533c2e958feSBarry Smith 534c2e958feSBarry Smith nzeros = 0; 535c2e958feSBarry Smith ia[0] = 0; 536c2e958feSBarry Smith for (i=0; i<m; i++) { 537c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 538c2e958feSBarry Smith cnt = 0; 539c2e958feSBarry Smith for (j=0; j<len; j++) { 540c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 5414c49b128SBarry Smith a[nzeros+cnt] = (int) PetscAbsScalar(ra[j]); 542c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 543c2e958feSBarry Smith } 544c2e958feSBarry Smith } 545c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 546c2e958feSBarry Smith nzeros += cnt; 547c2e958feSBarry Smith ia[i+1] = nzeros; 548c2e958feSBarry Smith } 549c2e958feSBarry Smith 550c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 551*f204ca49SKris Buschelman ierr = MatCreate(comm,m,N,PETSC_DETERMINE,N,&B);CHKERRQ(ierr); 552*f204ca49SKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 553*f204ca49SKris Buschelman ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 554c2e958feSBarry Smith 555676c34cdSKris Buschelman /* Fake support for "inplace" convert. */ 556676c34cdSKris Buschelman if (*newmat == A) { 557676c34cdSKris Buschelman ierr = MatDestroy(A);CHKERRQ(ierr); 558676c34cdSKris Buschelman } 559676c34cdSKris Buschelman *newmat = B; 560c2e958feSBarry Smith PetscFunctionReturn(0); 561c2e958feSBarry Smith } 562273d9f13SBarry Smith EXTERN_C_END 563433994e6SBarry Smith 564433994e6SBarry Smith 565433994e6SBarry Smith 566433994e6SBarry Smith 567433994e6SBarry Smith 568