1be1d678aSKris Buschelman 2b97cf49bSBarry Smith /* 3b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 4b97cf49bSBarry Smith */ 5e37d9405SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/ 6b97cf49bSBarry Smith 74a2ae208SSatish Balay #undef __FUNCT__ 84a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII" 9dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 10b97cf49bSBarry Smith { 113eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 12dfbe8321SBarry Smith PetscErrorCode ierr; 13d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 142dcb1b2aSMatthew Knepley const char *name; 15ce1411ecSBarry Smith PetscViewerFormat format; 16b97cf49bSBarry Smith 17433994e6SBarry Smith PetscFunctionBegin; 183a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 19b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 20fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 213a40ed3dSBarry Smith PetscFunctionReturn(0); 22fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 23e3c5b3baSBarry Smith SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MATLAB format not supported"); 240752156aSBarry Smith } else { 25d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 267b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 27b97cf49bSBarry Smith for (i=0; i<m; i++) { 28d0f46423SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr); 29b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 3077431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 310752156aSBarry Smith } 32b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 33b97cf49bSBarry Smith } 34d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 35b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 367b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 377b23a99aSBarry Smith } 383a40ed3dSBarry Smith PetscFunctionReturn(0); 39b97cf49bSBarry Smith } 40b97cf49bSBarry Smith 414a2ae208SSatish Balay #undef __FUNCT__ 424a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj" 43dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 44b97cf49bSBarry Smith { 45dfbe8321SBarry Smith PetscErrorCode ierr; 46ace3abfcSBarry Smith PetscBool iascii; 47b97cf49bSBarry Smith 48433994e6SBarry Smith PetscFunctionBegin; 49*251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5032077d6dSBarry Smith if (iascii) { 513eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 525cd90555SBarry Smith } else { 53e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 54b97cf49bSBarry Smith } 553a40ed3dSBarry Smith PetscFunctionReturn(0); 56b97cf49bSBarry Smith } 57b97cf49bSBarry Smith 584a2ae208SSatish Balay #undef __FUNCT__ 594a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 60dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat) 61b97cf49bSBarry Smith { 623eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 63dfbe8321SBarry Smith PetscErrorCode ierr; 6494d884c6SBarry Smith 6594d884c6SBarry Smith PetscFunctionBegin; 66aa482453SBarry Smith #if defined(PETSC_USE_LOG) 67d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz); 68b97cf49bSBarry Smith #endif 6905b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 700400557aSBarry Smith if (a->freeaij) { 7114196391SBarry Smith if (a->freeaijwithfree) { 7214196391SBarry Smith if (a->i) free(a->i); 7314196391SBarry Smith if (a->j) free(a->j); 7414196391SBarry Smith } else { 75606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 76606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 7705b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 7814196391SBarry Smith } 790400557aSBarry Smith } 80bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 81dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 82901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 839a3665c6SJed Brown ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C","",PETSC_NULL);CHKERRQ(ierr); 843a40ed3dSBarry Smith PetscFunctionReturn(0); 85b97cf49bSBarry Smith } 86b97cf49bSBarry Smith 874a2ae208SSatish Balay #undef __FUNCT__ 884a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 89ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 90b97cf49bSBarry Smith { 913eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 9263ba0a88SBarry Smith PetscErrorCode ierr; 93b97cf49bSBarry Smith 94433994e6SBarry Smith PetscFunctionBegin; 9512c028f9SKris Buschelman switch (op) { 9677e54ba9SKris Buschelman case MAT_SYMMETRIC: 9712c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 989a4540c5SBarry Smith case MAT_HERMITIAN: 994e0d8c25SBarry Smith a->symmetric = flg; 1009a4540c5SBarry Smith break; 1019a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1029a4540c5SBarry Smith break; 10312c028f9SKris Buschelman default: 104290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 10512c028f9SKris Buschelman break; 106b97cf49bSBarry Smith } 1073a40ed3dSBarry Smith PetscFunctionReturn(0); 108b97cf49bSBarry Smith } 109b97cf49bSBarry Smith 110b97cf49bSBarry Smith 111b97cf49bSBarry Smith /* 112b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 113b97cf49bSBarry Smith */ 114b97cf49bSBarry Smith 1154a2ae208SSatish Balay #undef __FUNCT__ 1164a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 117dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A) 118b97cf49bSBarry Smith { 1193eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 1206849ba73SBarry Smith PetscErrorCode ierr; 121d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 122b97cf49bSBarry Smith 123433994e6SBarry Smith PetscFunctionBegin; 12409f38230SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 12509f38230SBarry Smith ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr); 126d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 127b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 128b97cf49bSBarry Smith if (a->j[j] == i) { 12909f38230SBarry Smith a->diag[i] = j; 130b97cf49bSBarry Smith break; 131b97cf49bSBarry Smith } 132b97cf49bSBarry Smith } 133b97cf49bSBarry Smith } 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 135b97cf49bSBarry Smith } 136b97cf49bSBarry Smith 1374a2ae208SSatish Balay #undef __FUNCT__ 1384a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 139b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 140b97cf49bSBarry Smith { 1413eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 142b24ad042SBarry Smith PetscInt *itmp; 143b97cf49bSBarry Smith 144433994e6SBarry Smith PetscFunctionBegin; 145d0f46423SBarry Smith row -= A->rmap->rstart; 1460752156aSBarry Smith 147e32f2f54SBarry Smith if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 148b97cf49bSBarry Smith 149b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 150b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 151b97cf49bSBarry Smith if (idx) { 152b97cf49bSBarry Smith itmp = a->j + a->i[row]; 153b97cf49bSBarry Smith if (*nz) { 154b97cf49bSBarry Smith *idx = itmp; 155b97cf49bSBarry Smith } 156b97cf49bSBarry Smith else *idx = 0; 157b97cf49bSBarry Smith } 1583a40ed3dSBarry Smith PetscFunctionReturn(0); 159b97cf49bSBarry Smith } 160b97cf49bSBarry Smith 1614a2ae208SSatish Balay #undef __FUNCT__ 1624a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 163b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 164b97cf49bSBarry Smith { 165433994e6SBarry Smith PetscFunctionBegin; 1663a40ed3dSBarry Smith PetscFunctionReturn(0); 167b97cf49bSBarry Smith } 168b97cf49bSBarry Smith 1694a2ae208SSatish Balay #undef __FUNCT__ 1704a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 171ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 172b97cf49bSBarry Smith { 1733eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 174dfbe8321SBarry Smith PetscErrorCode ierr; 175ace3abfcSBarry Smith PetscBool flag; 176b97cf49bSBarry Smith 177433994e6SBarry Smith PetscFunctionBegin; 178b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 179d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 1800f5bd95cSBarry Smith flag = PETSC_FALSE; 181b97cf49bSBarry Smith } 182b97cf49bSBarry Smith 183b97cf49bSBarry Smith /* if the a->i are the same */ 184d0f46423SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 185b97cf49bSBarry Smith 186b97cf49bSBarry Smith /* if a->j are the same */ 187b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 188b97cf49bSBarry Smith 1897adad957SLisandro Dalcin ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr); 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 191b97cf49bSBarry Smith } 192b97cf49bSBarry Smith 1934a2ae208SSatish Balay #undef __FUNCT__ 1944a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 195ace3abfcSBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool *done) 1966cd91f64SBarry Smith { 197b24ad042SBarry Smith PetscInt i; 1986cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 1996cd91f64SBarry Smith 2006cd91f64SBarry Smith PetscFunctionBegin; 201d0f46423SBarry Smith *m = A->rmap->n; 2026cd91f64SBarry Smith *ia = a->i; 2036cd91f64SBarry Smith *ja = a->j; 2046cd91f64SBarry Smith *done = PETSC_TRUE; 205d42735a1SBarry Smith if (oshift) { 206d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 207d42735a1SBarry Smith (*ja)[i]++; 208d42735a1SBarry Smith } 209d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 210d42735a1SBarry Smith } 211d42735a1SBarry Smith PetscFunctionReturn(0); 212d42735a1SBarry Smith } 213d42735a1SBarry Smith 2144a2ae208SSatish Balay #undef __FUNCT__ 2154a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 216ace3abfcSBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool *done) 217d42735a1SBarry Smith { 218b24ad042SBarry Smith PetscInt i; 219d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 220d42735a1SBarry Smith 221d42735a1SBarry Smith PetscFunctionBegin; 222e32f2f54SBarry Smith if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 223e32f2f54SBarry Smith if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 224d42735a1SBarry Smith if (oshift) { 225d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 226d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 227d42735a1SBarry Smith (*ja)[i]--; 228d42735a1SBarry Smith } 229d42735a1SBarry Smith } 2306cd91f64SBarry Smith PetscFunctionReturn(0); 2316cd91f64SBarry Smith } 232b97cf49bSBarry Smith 23317667f90SBarry Smith #undef __FUNCT__ 23417667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj" 2357087cfbeSBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 23617667f90SBarry Smith { 23717667f90SBarry Smith Mat B; 23817667f90SBarry Smith PetscErrorCode ierr; 23917667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 24017667f90SBarry Smith const PetscInt *rj; 24117667f90SBarry Smith const PetscScalar *ra; 24217667f90SBarry Smith MPI_Comm comm; 24317667f90SBarry Smith 24417667f90SBarry Smith PetscFunctionBegin; 24517667f90SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 24617667f90SBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 24717667f90SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 24817667f90SBarry Smith 24917667f90SBarry Smith /* count the number of nonzeros per row */ 25017667f90SBarry Smith for (i=0; i<m; i++) { 25117667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 25217667f90SBarry Smith for (j=0; j<len; j++) { 25317667f90SBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 25417667f90SBarry Smith } 25517667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 25617667f90SBarry Smith nzeros += len; 25717667f90SBarry Smith } 25817667f90SBarry Smith 25917667f90SBarry Smith /* malloc space for nonzeros */ 26017667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr); 26117667f90SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr); 26217667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr); 26317667f90SBarry Smith 26417667f90SBarry Smith nzeros = 0; 26517667f90SBarry Smith ia[0] = 0; 26617667f90SBarry Smith for (i=0; i<m; i++) { 26717667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 26817667f90SBarry Smith cnt = 0; 26917667f90SBarry Smith for (j=0; j<len; j++) { 27017667f90SBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 27117667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 27217667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 27317667f90SBarry Smith } 27417667f90SBarry Smith } 27517667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 27617667f90SBarry Smith nzeros += cnt; 27717667f90SBarry Smith ia[i+1] = nzeros; 27817667f90SBarry Smith } 27917667f90SBarry Smith 28017667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 28117667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 28258ec18a5SLisandro Dalcin ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 28317667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 28417667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 28517667f90SBarry Smith 28617667f90SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 28717667f90SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 28817667f90SBarry Smith } else { 28917667f90SBarry Smith *newmat = B; 29017667f90SBarry Smith } 29117667f90SBarry Smith PetscFunctionReturn(0); 29217667f90SBarry Smith } 29317667f90SBarry Smith 294b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2950b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2963eda8832SBarry Smith MatGetRow_MPIAdj, 2973eda8832SBarry Smith MatRestoreRow_MPIAdj, 298b97cf49bSBarry Smith 0, 29997304618SKris Buschelman /* 4*/ 0, 300b97cf49bSBarry Smith 0, 301b97cf49bSBarry Smith 0, 302b97cf49bSBarry Smith 0, 3030b3b1487SBarry Smith 0, 3040b3b1487SBarry Smith 0, 30597304618SKris Buschelman /*10*/ 0, 3060b3b1487SBarry Smith 0, 3070b3b1487SBarry Smith 0, 3080b3b1487SBarry Smith 0, 3090b3b1487SBarry Smith 0, 31097304618SKris Buschelman /*15*/ 0, 3113eda8832SBarry Smith MatEqual_MPIAdj, 3120b3b1487SBarry Smith 0, 3130b3b1487SBarry Smith 0, 3140b3b1487SBarry Smith 0, 31597304618SKris Buschelman /*20*/ 0, 3160b3b1487SBarry Smith 0, 3173eda8832SBarry Smith MatSetOption_MPIAdj, 3180b3b1487SBarry Smith 0, 319d519adbfSMatthew Knepley /*24*/ 0, 3200b3b1487SBarry Smith 0, 3210b3b1487SBarry Smith 0, 3220b3b1487SBarry Smith 0, 3230b3b1487SBarry Smith 0, 324d519adbfSMatthew Knepley /*29*/ 0, 3250b3b1487SBarry Smith 0, 326273d9f13SBarry Smith 0, 327273d9f13SBarry Smith 0, 3280b3b1487SBarry Smith 0, 329d519adbfSMatthew Knepley /*34*/ 0, 3300b3b1487SBarry Smith 0, 3310b3b1487SBarry Smith 0, 3320b3b1487SBarry Smith 0, 3330b3b1487SBarry Smith 0, 334d519adbfSMatthew Knepley /*39*/ 0, 3350b3b1487SBarry Smith 0, 3360b3b1487SBarry Smith 0, 3370b3b1487SBarry Smith 0, 3380b3b1487SBarry Smith 0, 339d519adbfSMatthew Knepley /*44*/ 0, 3400b3b1487SBarry Smith 0, 3410b3b1487SBarry Smith 0, 3420b3b1487SBarry Smith 0, 3430b3b1487SBarry Smith 0, 344d519adbfSMatthew Knepley /*49*/ 0, 3456cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 346d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 347b97cf49bSBarry Smith 0, 348b97cf49bSBarry Smith 0, 349d519adbfSMatthew Knepley /*54*/ 0, 350b97cf49bSBarry Smith 0, 3510752156aSBarry Smith 0, 3520752156aSBarry Smith 0, 3530b3b1487SBarry Smith 0, 354d519adbfSMatthew Knepley /*59*/ 0, 355b9b97703SBarry Smith MatDestroy_MPIAdj, 356b9b97703SBarry Smith MatView_MPIAdj, 35717667f90SBarry Smith MatConvertFrom_MPIAdj, 35897304618SKris Buschelman 0, 359d519adbfSMatthew Knepley /*64*/ 0, 36097304618SKris Buschelman 0, 36197304618SKris Buschelman 0, 36297304618SKris Buschelman 0, 36397304618SKris Buschelman 0, 364d519adbfSMatthew Knepley /*69*/ 0, 36597304618SKris Buschelman 0, 36697304618SKris Buschelman 0, 36797304618SKris Buschelman 0, 36897304618SKris Buschelman 0, 369d519adbfSMatthew Knepley /*74*/ 0, 37097304618SKris Buschelman 0, 37197304618SKris Buschelman 0, 37297304618SKris Buschelman 0, 37397304618SKris Buschelman 0, 374d519adbfSMatthew Knepley /*79*/ 0, 37597304618SKris Buschelman 0, 37697304618SKris Buschelman 0, 37797304618SKris Buschelman 0, 378865e5f61SKris Buschelman 0, 379d519adbfSMatthew Knepley /*84*/ 0, 380865e5f61SKris Buschelman 0, 381865e5f61SKris Buschelman 0, 382865e5f61SKris Buschelman 0, 383865e5f61SKris Buschelman 0, 384d519adbfSMatthew Knepley /*89*/ 0, 385865e5f61SKris Buschelman 0, 386865e5f61SKris Buschelman 0, 387865e5f61SKris Buschelman 0, 388865e5f61SKris Buschelman 0, 389d519adbfSMatthew Knepley /*94*/ 0, 390865e5f61SKris Buschelman 0, 391865e5f61SKris Buschelman 0, 392865e5f61SKris Buschelman 0}; 393b97cf49bSBarry Smith 394a23d5eceSKris Buschelman EXTERN_C_BEGIN 395a23d5eceSKris Buschelman #undef __FUNCT__ 396a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 3977087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 398a23d5eceSKris Buschelman { 399a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 400dfbe8321SBarry Smith PetscErrorCode ierr; 4012515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 402b24ad042SBarry Smith PetscInt ii; 403a23d5eceSKris Buschelman #endif 404a23d5eceSKris Buschelman 405a23d5eceSKris Buschelman PetscFunctionBegin; 40626283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 40726283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 40858ec18a5SLisandro Dalcin 4092515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 410e32f2f54SBarry Smith if (i[0] != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]); 411d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 412e02043d6SBarry Smith if (i[ii] < 0 || i[ii] < i[ii-1]) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]); 413a23d5eceSKris Buschelman } 414d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 415e02043d6SBarry Smith if (j[ii] < 0 || j[ii] >= B->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]); 416a23d5eceSKris Buschelman } 417a23d5eceSKris Buschelman #endif 41858ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 419a23d5eceSKris Buschelman 420a23d5eceSKris Buschelman b->j = j; 421a23d5eceSKris Buschelman b->i = i; 422a23d5eceSKris Buschelman b->values = values; 423a23d5eceSKris Buschelman 424d0f46423SBarry Smith b->nz = i[B->rmap->n]; 425a23d5eceSKris Buschelman b->diag = 0; 426a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 427a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 428a23d5eceSKris Buschelman 429a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 430a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 431a23d5eceSKris Buschelman PetscFunctionReturn(0); 432a23d5eceSKris Buschelman } 433a23d5eceSKris Buschelman EXTERN_C_END 434b97cf49bSBarry Smith 4359a3665c6SJed Brown #undef __FUNCT__ 4369a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj" 4379a3665c6SJed Brown PETSC_EXTERN_C PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 4389a3665c6SJed Brown { 4399a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 4409a3665c6SJed Brown PetscErrorCode ierr; 4419a3665c6SJed Brown const PetscInt *ranges; 4429a3665c6SJed Brown MPI_Comm acomm,bcomm; 4439a3665c6SJed Brown MPI_Group agroup,bgroup; 4449a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 4459a3665c6SJed Brown 4469a3665c6SJed Brown PetscFunctionBegin; 4479a3665c6SJed Brown *B = PETSC_NULL; 4489a3665c6SJed Brown acomm = ((PetscObject)A)->comm; 4499a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr); 4509a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr); 4519a3665c6SJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 4529a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 4539a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 4549a3665c6SJed Brown } 4559a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 4569a3665c6SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 4579a3665c6SJed Brown *B = A; 4589a3665c6SJed Brown PetscFunctionReturn(0); 4599a3665c6SJed Brown } 4609a3665c6SJed Brown 4619a3665c6SJed Brown ierr = PetscMalloc(nranks*sizeof(PetscMPIInt),&ranks);CHKERRQ(ierr); 4629a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 4639a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 4649a3665c6SJed Brown } 4659a3665c6SJed Brown ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr); 4669a3665c6SJed Brown ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr); 4679a3665c6SJed Brown ierr = PetscFree(ranks);CHKERRQ(ierr); 4689a3665c6SJed Brown ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr); 4699a3665c6SJed Brown ierr = MPI_Group_free(&agroup);CHKERRQ(ierr); 4709a3665c6SJed Brown ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr); 4719a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 4729a3665c6SJed Brown PetscInt m,N; 4739a3665c6SJed Brown Mat_MPIAdj *b; 4749a3665c6SJed Brown ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 4759a3665c6SJed Brown ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 4769a3665c6SJed Brown ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 4779a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 4789a3665c6SJed Brown b->freeaij = PETSC_FALSE; 4799a3665c6SJed Brown ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr); 4809a3665c6SJed Brown } 4819a3665c6SJed Brown PetscFunctionReturn(0); 4829a3665c6SJed Brown } 4839a3665c6SJed Brown 4849a3665c6SJed Brown #undef __FUNCT__ 4859a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat" 4869a3665c6SJed Brown /*@ 4879a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 4889a3665c6SJed Brown 4899a3665c6SJed Brown Collective 4909a3665c6SJed Brown 4919a3665c6SJed Brown Input Arguments: 4929a3665c6SJed Brown . A - original MPIAdj matrix 4939a3665c6SJed Brown 4949a3665c6SJed Brown Output Arguments: 4959a3665c6SJed Brown . B - matrix on subcommunicator, PETSC_NULL on ranks that owned zero rows of A 4969a3665c6SJed Brown 4979a3665c6SJed Brown Level: developer 4989a3665c6SJed Brown 4999a3665c6SJed Brown Note: 5009a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 5019a3665c6SJed Brown 5029a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 5039a3665c6SJed Brown 5049a3665c6SJed Brown .seealso: MatCreateMPIAdj() 5059a3665c6SJed Brown @*/ 5069a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 5079a3665c6SJed Brown { 5089a3665c6SJed Brown PetscErrorCode ierr; 5099a3665c6SJed Brown 5109a3665c6SJed Brown PetscFunctionBegin; 5119a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 5129a3665c6SJed Brown ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 5139a3665c6SJed Brown PetscFunctionReturn(0); 5149a3665c6SJed Brown } 5159a3665c6SJed Brown 5160bad9183SKris Buschelman /*MC 517fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 5180bad9183SKris Buschelman intended for use constructing orderings and partitionings. 5190bad9183SKris Buschelman 5200bad9183SKris Buschelman Level: beginner 5210bad9183SKris Buschelman 5220bad9183SKris Buschelman .seealso: MatCreateMPIAdj 5230bad9183SKris Buschelman M*/ 5240bad9183SKris Buschelman 525273d9f13SBarry Smith EXTERN_C_BEGIN 5264a2ae208SSatish Balay #undef __FUNCT__ 5274a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 5287087cfbeSBarry Smith PetscErrorCode MatCreate_MPIAdj(Mat B) 529273d9f13SBarry Smith { 530273d9f13SBarry Smith Mat_MPIAdj *b; 5316849ba73SBarry Smith PetscErrorCode ierr; 532273d9f13SBarry Smith 533273d9f13SBarry Smith PetscFunctionBegin; 53438f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr); 535b0a32e0cSBarry Smith B->data = (void*)b; 536273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 537273d9f13SBarry Smith B->assembled = PETSC_FALSE; 538273d9f13SBarry Smith 539a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C", 540a23d5eceSKris Buschelman "MatMPIAdjSetPreallocation_MPIAdj", 541a23d5eceSKris Buschelman MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 5429a3665c6SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C", 5439a3665c6SJed Brown "MatMPIAdjCreateNonemptySubcommMat_MPIAdj", 5449a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 54517667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 546273d9f13SBarry Smith PetscFunctionReturn(0); 547273d9f13SBarry Smith } 548273d9f13SBarry Smith EXTERN_C_END 549273d9f13SBarry Smith 5504a2ae208SSatish Balay #undef __FUNCT__ 5514a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 552a23d5eceSKris Buschelman /*@C 553a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 554a23d5eceSKris Buschelman 5553f9fe445SBarry Smith Logically Collective on MPI_Comm 556a23d5eceSKris Buschelman 557a23d5eceSKris Buschelman Input Parameters: 558a23d5eceSKris Buschelman + A - the matrix 559a23d5eceSKris Buschelman . i - the indices into j for the start of each row 560a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 561a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 562a23d5eceSKris Buschelman - values - [optional] edge weights 563a23d5eceSKris Buschelman 564a23d5eceSKris Buschelman Level: intermediate 565a23d5eceSKris Buschelman 566a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 567a23d5eceSKris Buschelman @*/ 5687087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 569273d9f13SBarry Smith { 5704ac538c5SBarry Smith PetscErrorCode ierr; 571273d9f13SBarry Smith 572273d9f13SBarry Smith PetscFunctionBegin; 5734ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 574273d9f13SBarry Smith PetscFunctionReturn(0); 575273d9f13SBarry Smith } 576273d9f13SBarry Smith 5774a2ae208SSatish Balay #undef __FUNCT__ 5784a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 579b97cf49bSBarry Smith /*@C 5803eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 581b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 582b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 583b97cf49bSBarry Smith 584ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 585ef5ee4d1SLois Curfman McInnes 586b97cf49bSBarry Smith Input Parameters: 587c2e958feSBarry Smith + comm - MPI communicator 5880752156aSBarry Smith . m - number of local rows 58958ec18a5SLisandro Dalcin . N - number of global columns 590b97cf49bSBarry Smith . i - the indices into j for the start of each row 591329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 592ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 593329f5518SBarry Smith - values -[optional] edge weights 594b97cf49bSBarry Smith 595b97cf49bSBarry Smith Output Parameter: 596b97cf49bSBarry Smith . A - the matrix 597b97cf49bSBarry Smith 59815091d37SBarry Smith Level: intermediate 59915091d37SBarry Smith 6004bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 6014bc6d8c0SBarry Smith MatSetValues(). 602c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 603ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 6041198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 605327e1a73SBarry Smith Should not include the matrix diagonals. 606b97cf49bSBarry Smith 607e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 608ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 609ededeb1aSvictorle 610ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 611b97cf49bSBarry Smith 612e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 613b97cf49bSBarry Smith @*/ 6147087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 615b97cf49bSBarry Smith { 616dfbe8321SBarry Smith PetscErrorCode ierr; 617b97cf49bSBarry Smith 618433994e6SBarry Smith PetscFunctionBegin; 619f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 62058ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 621273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 622273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 6233a40ed3dSBarry Smith PetscFunctionReturn(0); 624b97cf49bSBarry Smith } 625b97cf49bSBarry Smith 626c2e958feSBarry Smith 627433994e6SBarry Smith 628433994e6SBarry Smith 629433994e6SBarry Smith 630433994e6SBarry Smith 631433994e6SBarry Smith 632