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*/ 6841d17a1SFande Kong #include <petscsf.h> 7b97cf49bSBarry Smith 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII" 10dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 11b97cf49bSBarry Smith { 123eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 13dfbe8321SBarry Smith PetscErrorCode ierr; 14d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 152dcb1b2aSMatthew Knepley const 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) { 24ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 250752156aSBarry Smith } else { 26d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 277b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 28b97cf49bSBarry Smith for (i=0; i<m; i++) { 29d0f46423SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr); 30b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 3177431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 320752156aSBarry Smith } 33b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 34b97cf49bSBarry Smith } 35d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 36b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 377b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 387b23a99aSBarry Smith } 393a40ed3dSBarry Smith PetscFunctionReturn(0); 40b97cf49bSBarry Smith } 41b97cf49bSBarry Smith 424a2ae208SSatish Balay #undef __FUNCT__ 434a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj" 44dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 45b97cf49bSBarry Smith { 46dfbe8321SBarry Smith PetscErrorCode ierr; 47ace3abfcSBarry Smith PetscBool iascii; 48b97cf49bSBarry Smith 49433994e6SBarry Smith PetscFunctionBegin; 50251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5132077d6dSBarry Smith if (iascii) { 523eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 53b97cf49bSBarry Smith } 543a40ed3dSBarry Smith PetscFunctionReturn(0); 55b97cf49bSBarry Smith } 56b97cf49bSBarry Smith 574a2ae208SSatish Balay #undef __FUNCT__ 584a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 59dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat) 60b97cf49bSBarry Smith { 613eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 62dfbe8321SBarry Smith PetscErrorCode ierr; 6394d884c6SBarry Smith 6494d884c6SBarry Smith PetscFunctionBegin; 65aa482453SBarry Smith #if defined(PETSC_USE_LOG) 66d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz); 67b97cf49bSBarry Smith #endif 6805b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 690400557aSBarry Smith if (a->freeaij) { 7014196391SBarry Smith if (a->freeaijwithfree) { 7114196391SBarry Smith if (a->i) free(a->i); 7214196391SBarry Smith if (a->j) free(a->j); 7314196391SBarry Smith } else { 74606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 75606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 7605b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 7714196391SBarry Smith } 780400557aSBarry Smith } 792b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 80bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 81dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 82bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr); 83bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",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; 124785e854fSJed Brown ierr = PetscMalloc1(m,&a->diag);CHKERRQ(ierr); 1253bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)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; 1422b1d8763SJed Brown PetscErrorCode ierr; 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]; 1502b1d8763SJed Brown if (v) { 1512b1d8763SJed Brown PetscInt j; 1522b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 1532b1d8763SJed Brown ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 1542b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 155785e854fSJed Brown ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr); 156b97cf49bSBarry Smith } 1572b1d8763SJed Brown for (j=0; j<*nz; j++) a->rowvalues[j] = a->values[a->i[row]+j]; 1582b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 159b97cf49bSBarry Smith } 16092b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 1613a40ed3dSBarry Smith PetscFunctionReturn(0); 162b97cf49bSBarry Smith } 163b97cf49bSBarry Smith 1644a2ae208SSatish Balay #undef __FUNCT__ 1654a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 166b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 167b97cf49bSBarry Smith { 168433994e6SBarry Smith PetscFunctionBegin; 1693a40ed3dSBarry Smith PetscFunctionReturn(0); 170b97cf49bSBarry Smith } 171b97cf49bSBarry Smith 1724a2ae208SSatish Balay #undef __FUNCT__ 1734a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 174ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 175b97cf49bSBarry Smith { 1763eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 177dfbe8321SBarry Smith PetscErrorCode ierr; 178ace3abfcSBarry Smith PetscBool flag; 179b97cf49bSBarry Smith 180433994e6SBarry Smith PetscFunctionBegin; 181b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 182d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 1830f5bd95cSBarry Smith flag = PETSC_FALSE; 184b97cf49bSBarry Smith } 185b97cf49bSBarry Smith 186b97cf49bSBarry Smith /* if the a->i are the same */ 187d0f46423SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 188b97cf49bSBarry Smith 189b97cf49bSBarry Smith /* if a->j are the same */ 190b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 191b97cf49bSBarry Smith 192ce94432eSBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 1933a40ed3dSBarry Smith PetscFunctionReturn(0); 194b97cf49bSBarry Smith } 195b97cf49bSBarry Smith 1964a2ae208SSatish Balay #undef __FUNCT__ 1974a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 1981a83f524SJed Brown PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 1996cd91f64SBarry Smith { 200b24ad042SBarry Smith PetscInt i; 2016cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 2021a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 2036cd91f64SBarry Smith 2046cd91f64SBarry Smith PetscFunctionBegin; 205d0f46423SBarry Smith *m = A->rmap->n; 2066cd91f64SBarry Smith *ia = a->i; 2076cd91f64SBarry Smith *ja = a->j; 2086cd91f64SBarry Smith *done = PETSC_TRUE; 209d42735a1SBarry Smith if (oshift) { 210d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 211d42735a1SBarry Smith (*ja)[i]++; 212d42735a1SBarry Smith } 213d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 214d42735a1SBarry Smith } 215d42735a1SBarry Smith PetscFunctionReturn(0); 216d42735a1SBarry Smith } 217d42735a1SBarry Smith 2184a2ae208SSatish Balay #undef __FUNCT__ 2194a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 2201a83f524SJed Brown PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 221d42735a1SBarry Smith { 222b24ad042SBarry Smith PetscInt i; 223d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 2241a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 225d42735a1SBarry Smith 226d42735a1SBarry Smith PetscFunctionBegin; 227e32f2f54SBarry Smith if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 228e32f2f54SBarry Smith if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 229d42735a1SBarry Smith if (oshift) { 230d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 231d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 232d42735a1SBarry Smith (*ja)[i]--; 233d42735a1SBarry Smith } 234d42735a1SBarry Smith } 2356cd91f64SBarry Smith PetscFunctionReturn(0); 2366cd91f64SBarry Smith } 237b97cf49bSBarry Smith 23817667f90SBarry Smith #undef __FUNCT__ 23917667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj" 24019fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 24117667f90SBarry Smith { 24217667f90SBarry Smith Mat B; 24317667f90SBarry Smith PetscErrorCode ierr; 24417667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 24517667f90SBarry Smith const PetscInt *rj; 24617667f90SBarry Smith const PetscScalar *ra; 24717667f90SBarry Smith MPI_Comm comm; 24817667f90SBarry Smith 24917667f90SBarry Smith PetscFunctionBegin; 2500298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 2510298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 2520298fd71SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 25317667f90SBarry Smith 25417667f90SBarry Smith /* count the number of nonzeros per row */ 25517667f90SBarry Smith for (i=0; i<m; i++) { 2560298fd71SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 2575ee9ba1cSJed Brown for (j=0; j<len; j++) { 2585ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 2595ee9ba1cSJed Brown } 26017667f90SBarry Smith nzeros += len; 2610f2063bfSTobin Isaac ierr = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 26217667f90SBarry Smith } 26317667f90SBarry Smith 26417667f90SBarry Smith /* malloc space for nonzeros */ 265854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr); 266854ce69bSBarry Smith ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr); 267854ce69bSBarry Smith ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr); 26817667f90SBarry Smith 26917667f90SBarry Smith nzeros = 0; 27017667f90SBarry Smith ia[0] = 0; 27117667f90SBarry Smith for (i=0; i<m; i++) { 27217667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 27317667f90SBarry Smith cnt = 0; 27417667f90SBarry Smith for (j=0; j<len; j++) { 2755ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 27617667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 27717667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 27817667f90SBarry Smith } 2795ee9ba1cSJed Brown } 28017667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 28117667f90SBarry Smith nzeros += cnt; 28217667f90SBarry Smith ia[i+1] = nzeros; 28317667f90SBarry Smith } 28417667f90SBarry Smith 28517667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 28617667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 28758ec18a5SLisandro Dalcin ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 28817667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 28917667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 29017667f90SBarry Smith 29117667f90SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 29217667f90SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 29317667f90SBarry Smith } else { 29417667f90SBarry Smith *newmat = B; 29517667f90SBarry Smith } 29617667f90SBarry Smith PetscFunctionReturn(0); 29717667f90SBarry Smith } 29817667f90SBarry Smith 299b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 3000b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 3013eda8832SBarry Smith MatGetRow_MPIAdj, 3023eda8832SBarry Smith MatRestoreRow_MPIAdj, 303b97cf49bSBarry Smith 0, 30497304618SKris Buschelman /* 4*/ 0, 305b97cf49bSBarry Smith 0, 306b97cf49bSBarry Smith 0, 307b97cf49bSBarry Smith 0, 3080b3b1487SBarry Smith 0, 3090b3b1487SBarry Smith 0, 31097304618SKris Buschelman /*10*/ 0, 3110b3b1487SBarry Smith 0, 3120b3b1487SBarry Smith 0, 3130b3b1487SBarry Smith 0, 3140b3b1487SBarry Smith 0, 31597304618SKris Buschelman /*15*/ 0, 3163eda8832SBarry Smith MatEqual_MPIAdj, 3170b3b1487SBarry Smith 0, 3180b3b1487SBarry Smith 0, 3190b3b1487SBarry Smith 0, 32097304618SKris Buschelman /*20*/ 0, 3210b3b1487SBarry Smith 0, 3223eda8832SBarry Smith MatSetOption_MPIAdj, 3230b3b1487SBarry Smith 0, 324d519adbfSMatthew Knepley /*24*/ 0, 3250b3b1487SBarry Smith 0, 3260b3b1487SBarry Smith 0, 3270b3b1487SBarry Smith 0, 3280b3b1487SBarry Smith 0, 329d519adbfSMatthew Knepley /*29*/ 0, 3300b3b1487SBarry Smith 0, 331273d9f13SBarry Smith 0, 332273d9f13SBarry Smith 0, 3330b3b1487SBarry Smith 0, 334d519adbfSMatthew Knepley /*34*/ 0, 3350b3b1487SBarry Smith 0, 3360b3b1487SBarry Smith 0, 3370b3b1487SBarry Smith 0, 3380b3b1487SBarry Smith 0, 339d519adbfSMatthew Knepley /*39*/ 0, 3400b3b1487SBarry Smith 0, 3410b3b1487SBarry Smith 0, 3420b3b1487SBarry Smith 0, 3430b3b1487SBarry Smith 0, 344d519adbfSMatthew Knepley /*44*/ 0, 3450b3b1487SBarry Smith 0, 3467d68702bSBarry Smith MatShift_Basic, 3470b3b1487SBarry Smith 0, 3480b3b1487SBarry Smith 0, 349d519adbfSMatthew Knepley /*49*/ 0, 3506cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 351d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 352b97cf49bSBarry Smith 0, 353b97cf49bSBarry Smith 0, 354d519adbfSMatthew Knepley /*54*/ 0, 355b97cf49bSBarry Smith 0, 3560752156aSBarry Smith 0, 3570752156aSBarry Smith 0, 3580b3b1487SBarry Smith 0, 359d519adbfSMatthew Knepley /*59*/ 0, 360b9b97703SBarry Smith MatDestroy_MPIAdj, 361b9b97703SBarry Smith MatView_MPIAdj, 36217667f90SBarry Smith MatConvertFrom_MPIAdj, 36397304618SKris Buschelman 0, 364d519adbfSMatthew Knepley /*64*/ 0, 36597304618SKris Buschelman 0, 36697304618SKris Buschelman 0, 36797304618SKris Buschelman 0, 36897304618SKris Buschelman 0, 369d519adbfSMatthew Knepley /*69*/ 0, 37097304618SKris Buschelman 0, 37197304618SKris Buschelman 0, 37297304618SKris Buschelman 0, 37397304618SKris Buschelman 0, 374d519adbfSMatthew Knepley /*74*/ 0, 37597304618SKris Buschelman 0, 37697304618SKris Buschelman 0, 37797304618SKris Buschelman 0, 37897304618SKris Buschelman 0, 379d519adbfSMatthew Knepley /*79*/ 0, 38097304618SKris Buschelman 0, 38197304618SKris Buschelman 0, 38297304618SKris Buschelman 0, 383865e5f61SKris Buschelman 0, 384d519adbfSMatthew Knepley /*84*/ 0, 385865e5f61SKris Buschelman 0, 386865e5f61SKris Buschelman 0, 387865e5f61SKris Buschelman 0, 388865e5f61SKris Buschelman 0, 389d519adbfSMatthew Knepley /*89*/ 0, 390865e5f61SKris Buschelman 0, 391865e5f61SKris Buschelman 0, 392865e5f61SKris Buschelman 0, 393865e5f61SKris Buschelman 0, 394d519adbfSMatthew Knepley /*94*/ 0, 395865e5f61SKris Buschelman 0, 396865e5f61SKris Buschelman 0, 3973964eb88SJed Brown 0, 3983964eb88SJed Brown 0, 3993964eb88SJed Brown /*99*/ 0, 4003964eb88SJed Brown 0, 4013964eb88SJed Brown 0, 4023964eb88SJed Brown 0, 4033964eb88SJed Brown 0, 4043964eb88SJed Brown /*104*/ 0, 4053964eb88SJed Brown 0, 4063964eb88SJed Brown 0, 4073964eb88SJed Brown 0, 4083964eb88SJed Brown 0, 4093964eb88SJed Brown /*109*/ 0, 4103964eb88SJed Brown 0, 4113964eb88SJed Brown 0, 4123964eb88SJed Brown 0, 4133964eb88SJed Brown 0, 4143964eb88SJed Brown /*114*/ 0, 4153964eb88SJed Brown 0, 4163964eb88SJed Brown 0, 4173964eb88SJed Brown 0, 4183964eb88SJed Brown 0, 4193964eb88SJed Brown /*119*/ 0, 4203964eb88SJed Brown 0, 4213964eb88SJed Brown 0, 4223964eb88SJed Brown 0, 4233964eb88SJed Brown 0, 4243964eb88SJed Brown /*124*/ 0, 4253964eb88SJed Brown 0, 4263964eb88SJed Brown 0, 4273964eb88SJed Brown 0, 4283964eb88SJed Brown 0, 4293964eb88SJed Brown /*129*/ 0, 4303964eb88SJed Brown 0, 4313964eb88SJed Brown 0, 4323964eb88SJed Brown 0, 4333964eb88SJed Brown 0, 4343964eb88SJed Brown /*134*/ 0, 4353964eb88SJed Brown 0, 4363964eb88SJed Brown 0, 4373964eb88SJed Brown 0, 4383964eb88SJed Brown 0, 4393964eb88SJed Brown /*139*/ 0, 440f9426fe0SMark Adams 0, 4413964eb88SJed Brown 0 4423964eb88SJed Brown }; 443b97cf49bSBarry Smith 444a23d5eceSKris Buschelman #undef __FUNCT__ 445a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 446f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 447a23d5eceSKris Buschelman { 448a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 449dfbe8321SBarry Smith PetscErrorCode ierr; 4502515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 451b24ad042SBarry Smith PetscInt ii; 452a23d5eceSKris Buschelman #endif 453a23d5eceSKris Buschelman 454a23d5eceSKris Buschelman PetscFunctionBegin; 45526283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 45626283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 45758ec18a5SLisandro Dalcin 4582515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 459e32f2f54SBarry 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]); 460d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 461e02043d6SBarry 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]); 462a23d5eceSKris Buschelman } 463d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 464e02043d6SBarry 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]); 465a23d5eceSKris Buschelman } 466a23d5eceSKris Buschelman #endif 46758ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 468a23d5eceSKris Buschelman 469a23d5eceSKris Buschelman b->j = j; 470a23d5eceSKris Buschelman b->i = i; 471a23d5eceSKris Buschelman b->values = values; 472a23d5eceSKris Buschelman 473d0f46423SBarry Smith b->nz = i[B->rmap->n]; 474a23d5eceSKris Buschelman b->diag = 0; 475a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 476a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 477a23d5eceSKris Buschelman 478a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 479a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 480a23d5eceSKris Buschelman PetscFunctionReturn(0); 481a23d5eceSKris Buschelman } 482b97cf49bSBarry Smith 483*85374dbcSFande Kong static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat adj,IS rows, IS cols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values); 484*85374dbcSFande Kong 4859a3665c6SJed Brown #undef __FUNCT__ 486*85374dbcSFande Kong #define __FUNCT__ "MatGetSubMatrices_MPIAdj" 487*85374dbcSFande Kong PetscErrorCode MatGetSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 488*85374dbcSFande Kong { 489*85374dbcSFande Kong PetscInt i,irow_localsize,icol_localsize,*sxadj,*sadjncy,*svalues; 490*85374dbcSFande Kong PetscInt *indices,nindices,j,k,loc; 491*85374dbcSFande Kong const PetscInt *irow_indices,*icol_indices; 492*85374dbcSFande Kong PetscErrorCode ierr; 493*85374dbcSFande Kong 494*85374dbcSFande Kong PetscFunctionBegin; 495*85374dbcSFande Kong nindices = 0; 496*85374dbcSFande Kong for(i=0; i<n; i++){ 497*85374dbcSFande Kong ierr = ISGetLocalSize(irow[i],&irow_localsize);CHKERRQ(ierr); 498*85374dbcSFande Kong ierr = ISGetLocalSize(icol[i],&icol_localsize);CHKERRQ(ierr); 499*85374dbcSFande Kong nindices = nindices>(irow_localsize+icol_localsize)? nindices:(irow_localsize+icol_localsize); 500*85374dbcSFande Kong } 501*85374dbcSFande Kong ierr = PetscCalloc1(nindices,&indices);CHKERRQ(ierr); 502*85374dbcSFande Kong for(i=0; i<n; i++){ 503*85374dbcSFande Kong ierr = MatGetSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr); 504*85374dbcSFande Kong ierr = ISGetLocalSize(irow[i],&irow_localsize);CHKERRQ(ierr); 505*85374dbcSFande Kong ierr = ISGetLocalSize(icol[i],&icol_localsize);CHKERRQ(ierr); 506*85374dbcSFande Kong ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr); 507*85374dbcSFande Kong ierr = PetscMemcpy(indices,irow_indices,sizeof(PetscInt)*irow_localsize);CHKERRQ(ierr); 508*85374dbcSFande Kong ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr); 509*85374dbcSFande Kong ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr); 510*85374dbcSFande Kong ierr = PetscMemcpy(indices+irow_localsize,icol_indices,sizeof(PetscInt)*icol_localsize);CHKERRQ(ierr); 511*85374dbcSFande Kong ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr); 512*85374dbcSFande Kong nindices = irow_localsize+icol_localsize; 513*85374dbcSFande Kong ierr = PetscSortRemoveDupsInt(&nindices,indices);CHKERRQ(ierr); 514*85374dbcSFande Kong for(j=0; j<irow_localsize; j++){ 515*85374dbcSFande Kong for(k=sxadj[j]; k<sxadj[j]; k++){ 516*85374dbcSFande Kong ierr = PetscFindInt(sadjncy[k],nindices,indices,&loc);CHKERRQ(ierr); 517*85374dbcSFande Kong #if PETSC_USE_DEBUG 518*85374dbcSFande Kong if(loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %d \n",sadjncy[k]); 519*85374dbcSFande Kong #endif 520*85374dbcSFande Kong sadjncy[k] = loc; 521*85374dbcSFande Kong } 522*85374dbcSFande Kong } 523*85374dbcSFande Kong if(scall==MAT_INITIAL_MATRIX){ 524*85374dbcSFande Kong ierr = MatCreateMPIAdj(PETSC_COMM_SELF,irow_localsize,icol_localsize,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr); 525*85374dbcSFande Kong }else{ 526*85374dbcSFande Kong Mat sadj = *(submat[i]); 527*85374dbcSFande Kong Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 528*85374dbcSFande Kong ierr = PetscMemcpy(sa->i,sxadj,sizeof(PetscInt)*(irow_localsize+1));CHKERRQ(ierr); 529*85374dbcSFande Kong ierr = PetscMemcpy(sa->j,sadjncy,sizeof(PetscInt)*sxadj[irow_localsize]);CHKERRQ(ierr); 530*85374dbcSFande Kong ierr = PetscMemcpy(sa->values,svalues,sizeof(PetscInt)*sxadj[irow_localsize]);CHKERRQ(ierr); 531*85374dbcSFande Kong ierr = PetscFree(sxadj);CHKERRQ(ierr); 532*85374dbcSFande Kong ierr = PetscFree(sadjncy);CHKERRQ(ierr); 533*85374dbcSFande Kong ierr = PetscFree(svalues);CHKERRQ(ierr); 534*85374dbcSFande Kong } 535*85374dbcSFande Kong } 536*85374dbcSFande Kong PetscFunctionReturn(0); 537*85374dbcSFande Kong } 538*85374dbcSFande Kong 539*85374dbcSFande Kong 540*85374dbcSFande Kong #undef __FUNCT__ 541*85374dbcSFande Kong #define __FUNCT__ "MatGetSubMatrix_MPIAdj_data" 542*85374dbcSFande Kong static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat adj,IS rows, IS cols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values) 543841d17a1SFande Kong { 544841d17a1SFande Kong PetscInt rows_localsize,cols_localsize,i,j,nroots,nleaves,owner,localindex,*ncols_send,*ncols_recv; 545*85374dbcSFande Kong PetscInt nlocalrows,*adjncy_recv,Ncols_recv,Ncols_send,*xadj_recv,*values_recv; 546*85374dbcSFande Kong PetscInt *ncols_recv_offsets,loc,rnclos,*sadjncy,*sxadj,*svalues; 547841d17a1SFande Kong const PetscInt *rows_indices,*cols_indices,*xadj, *adjncy; 548841d17a1SFande Kong Mat_MPIAdj *a = (Mat_MPIAdj*)adj->data; 549841d17a1SFande Kong PetscLayout rmap; 550841d17a1SFande Kong MPI_Comm comm; 551841d17a1SFande Kong PetscSF sf; 552841d17a1SFande Kong PetscSFNode *iremote; 553841d17a1SFande Kong PetscBool done; 554841d17a1SFande Kong PetscErrorCode ierr; 555841d17a1SFande Kong 556841d17a1SFande Kong PetscFunctionBegin; 557841d17a1SFande Kong /* communicator */ 558841d17a1SFande Kong ierr = PetscObjectGetComm((PetscObject)adj,&comm);CHKERRQ(ierr); 559841d17a1SFande Kong /* Layouts */ 560841d17a1SFande Kong ierr = MatGetLayouts(adj,&rmap,PETSC_NULL);CHKERRQ(ierr); 561841d17a1SFande Kong /* get rows information */ 562841d17a1SFande Kong ierr = ISGetLocalSize(rows,&rows_localsize);CHKERRQ(ierr); 563841d17a1SFande Kong ierr = ISGetIndices(rows,&rows_indices);CHKERRQ(ierr); 564841d17a1SFande Kong ierr = PetscCalloc1(rows_localsize,&iremote);CHKERRQ(ierr); 565841d17a1SFande Kong 566841d17a1SFande Kong nleaves = rows_localsize; 567841d17a1SFande Kong for(i=0; i<rows_localsize; i++){ 568841d17a1SFande Kong ierr = PetscLayoutFindOwnerIndex(rmap,rows_indices[i],&owner,&localindex);CHKERRQ(ierr); 569841d17a1SFande Kong iremote[i].rank = owner; 570841d17a1SFande Kong iremote[i].index = localindex; 571841d17a1SFande Kong } 572841d17a1SFande Kong ierr = MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlocalrows,&xadj,&adjncy,&done);CHKERRQ(ierr); 573*85374dbcSFande Kong ierr = PetscCalloc4(nlocalrows,&ncols_send,rows_localsize,&xadj_recv,rows_localsize+1,&ncols_recv_offsets,rows_localsize,&ncols_recv);CHKERRQ(ierr); 574841d17a1SFande Kong nroots = nlocalrows; 575841d17a1SFande Kong for(i=0; i<nlocalrows; i++){ 576841d17a1SFande Kong ncols_send[i] = xadj[i+1]-xadj[i]; 577841d17a1SFande Kong } 578841d17a1SFande Kong ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 579*85374dbcSFande Kong ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 580841d17a1SFande Kong ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 581841d17a1SFande Kong ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 582841d17a1SFande Kong ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 583841d17a1SFande Kong ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 584841d17a1SFande Kong ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 585841d17a1SFande Kong ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 586841d17a1SFande Kong ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 587841d17a1SFande Kong Ncols_recv =0; 588841d17a1SFande Kong for(i=0; i<rows_localsize; i++){ 589841d17a1SFande Kong Ncols_recv += ncols_recv[i]; 590*85374dbcSFande Kong ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i]; 591841d17a1SFande Kong } 592841d17a1SFande Kong Ncols_send = 0; 593841d17a1SFande Kong for(i=0; i<nlocalrows; i++){ 594841d17a1SFande Kong Ncols_send += ncols_send[i]; 595841d17a1SFande Kong } 596841d17a1SFande Kong ierr = PetscCalloc1(Ncols_recv,&iremote);CHKERRQ(ierr); 597841d17a1SFande Kong ierr = PetscCalloc2(Ncols_recv,&adjncy_recv,Ncols_recv,&values_recv);CHKERRQ(ierr); 598841d17a1SFande Kong nleaves = Ncols_recv; 599841d17a1SFande Kong for(i=0; i<rows_localsize; i++){ 600841d17a1SFande Kong ierr = PetscLayoutFindOwner(rmap,rows_indices[i],&owner);CHKERRQ(ierr); 601841d17a1SFande Kong iremote[i].rank = owner; 602841d17a1SFande Kong for(j=0; j<ncols_recv[i]; j++){ 603841d17a1SFande Kong iremote[i].index =xadj_recv[i]+j; 604841d17a1SFande Kong } 605841d17a1SFande Kong } 606*85374dbcSFande Kong ierr = ISRestoreIndices(rows,&rows_indices);CHKERRQ(ierr); 607841d17a1SFande Kong nroots = Ncols_send; 608841d17a1SFande Kong ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 609*85374dbcSFande Kong ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 610841d17a1SFande Kong ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 611841d17a1SFande Kong ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 612841d17a1SFande Kong ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 613841d17a1SFande Kong ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 614841d17a1SFande Kong ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 615841d17a1SFande Kong ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 616841d17a1SFande Kong ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 617841d17a1SFande Kong ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlocalrows,&xadj,&adjncy,&done);CHKERRQ(ierr); 618*85374dbcSFande Kong ierr = ISGetLocalSize(cols,&cols_localsize);CHKERRQ(ierr); 619*85374dbcSFande Kong ierr = ISGetIndices(cols,&cols_indices);CHKERRQ(ierr); 620*85374dbcSFande Kong rnclos = 0; 621*85374dbcSFande Kong for(i=0; i<rows_localsize; i++){ 622*85374dbcSFande Kong for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 623*85374dbcSFande Kong ierr = PetscFindInt(adjncy_recv[j], cols_localsize, cols_indices, &loc);CHKERRQ(ierr); 624*85374dbcSFande Kong if(loc<0){ 625*85374dbcSFande Kong adjncy_recv[j] = -1; 626*85374dbcSFande Kong values_recv[j] = -1; 627*85374dbcSFande Kong ncols_recv[i]--; 628*85374dbcSFande Kong }else{ 629*85374dbcSFande Kong rnclos++; 630*85374dbcSFande Kong } 631*85374dbcSFande Kong } 632*85374dbcSFande Kong } 633*85374dbcSFande Kong ierr = ISRestoreIndices(cols,&cols_indices);CHKERRQ(ierr); 634*85374dbcSFande Kong ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr); 635*85374dbcSFande Kong ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr); 636*85374dbcSFande Kong ierr = PetscCalloc1(rows_localsize+1,&sxadj);CHKERRQ(ierr); 637*85374dbcSFande Kong rnclos = 0; 638*85374dbcSFande Kong for(i=0; i<rows_localsize; i++){ 639*85374dbcSFande Kong for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 640*85374dbcSFande Kong if(adjncy_recv[j]<0) continue; 641*85374dbcSFande Kong sadjncy[rnclos] = adjncy_recv[j]; 642*85374dbcSFande Kong svalues[rnclos] = values_recv[j]; 643*85374dbcSFande Kong rnclos++; 644*85374dbcSFande Kong } 645*85374dbcSFande Kong } 646*85374dbcSFande Kong for(i=0; i<rows_localsize; i++){ 647*85374dbcSFande Kong sxadj[i+1] = sxadj[i]+ncols_recv[i]; 648*85374dbcSFande Kong } 649*85374dbcSFande Kong if(sadj_xadj) { *sadj_xadj = sxadj;}else { ierr = PetscFree(sxadj);CHKERRQ(ierr);} 650*85374dbcSFande Kong if(sadj_adjncy){ *sadj_adjncy = sadjncy;}else{ ierr = PetscFree(sadjncy);CHKERRQ(ierr);} 651*85374dbcSFande Kong if(sadj_values){ *sadj_values = svalues;}else{ ierr = PetscFree(svalues);CHKERRQ(ierr);} 652*85374dbcSFande Kong ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr); 653*85374dbcSFande Kong ierr = PetscFree2(adjncy_recv,values_recv);CHKERRQ(ierr); 654841d17a1SFande Kong PetscFunctionReturn(0); 655841d17a1SFande Kong } 656841d17a1SFande Kong 657841d17a1SFande Kong 658841d17a1SFande Kong #undef __FUNCT__ 6599a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj" 660f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 6619a3665c6SJed Brown { 6629a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 6639a3665c6SJed Brown PetscErrorCode ierr; 6649a3665c6SJed Brown const PetscInt *ranges; 6659a3665c6SJed Brown MPI_Comm acomm,bcomm; 6669a3665c6SJed Brown MPI_Group agroup,bgroup; 6679a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 6689a3665c6SJed Brown 6699a3665c6SJed Brown PetscFunctionBegin; 6700298fd71SBarry Smith *B = NULL; 671ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr); 6729a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr); 6739a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr); 6749a3665c6SJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 6759a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6769a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 6779a3665c6SJed Brown } 6789a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 6799a3665c6SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 6809a3665c6SJed Brown *B = A; 6819a3665c6SJed Brown PetscFunctionReturn(0); 6829a3665c6SJed Brown } 6839a3665c6SJed Brown 684785e854fSJed Brown ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr); 6859a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 6869a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 6879a3665c6SJed Brown } 6889a3665c6SJed Brown ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr); 6899a3665c6SJed Brown ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr); 6909a3665c6SJed Brown ierr = PetscFree(ranks);CHKERRQ(ierr); 6919a3665c6SJed Brown ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr); 6929a3665c6SJed Brown ierr = MPI_Group_free(&agroup);CHKERRQ(ierr); 6939a3665c6SJed Brown ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr); 6949a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 6959a3665c6SJed Brown PetscInt m,N; 6969a3665c6SJed Brown Mat_MPIAdj *b; 6970298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 6980298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 6999a3665c6SJed Brown ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 7009a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 7019a3665c6SJed Brown b->freeaij = PETSC_FALSE; 7029a3665c6SJed Brown ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr); 7039a3665c6SJed Brown } 7049a3665c6SJed Brown PetscFunctionReturn(0); 7059a3665c6SJed Brown } 7069a3665c6SJed Brown 7079a3665c6SJed Brown #undef __FUNCT__ 7089a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat" 7099a3665c6SJed Brown /*@ 7109a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 7119a3665c6SJed Brown 7129a3665c6SJed Brown Collective 7139a3665c6SJed Brown 7149a3665c6SJed Brown Input Arguments: 7159a3665c6SJed Brown . A - original MPIAdj matrix 7169a3665c6SJed Brown 7179a3665c6SJed Brown Output Arguments: 7180298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 7199a3665c6SJed Brown 7209a3665c6SJed Brown Level: developer 7219a3665c6SJed Brown 7229a3665c6SJed Brown Note: 7239a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 7249a3665c6SJed Brown 7259a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 7269a3665c6SJed Brown 7279a3665c6SJed Brown .seealso: MatCreateMPIAdj() 7289a3665c6SJed Brown @*/ 7299a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 7309a3665c6SJed Brown { 7319a3665c6SJed Brown PetscErrorCode ierr; 7329a3665c6SJed Brown 7339a3665c6SJed Brown PetscFunctionBegin; 7349a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 7359a3665c6SJed Brown ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 7369a3665c6SJed Brown PetscFunctionReturn(0); 7379a3665c6SJed Brown } 7389a3665c6SJed Brown 7390bad9183SKris Buschelman /*MC 740fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 7410bad9183SKris Buschelman intended for use constructing orderings and partitionings. 7420bad9183SKris Buschelman 7430bad9183SKris Buschelman Level: beginner 7440bad9183SKris Buschelman 7450bad9183SKris Buschelman .seealso: MatCreateMPIAdj 7460bad9183SKris Buschelman M*/ 7470bad9183SKris Buschelman 7484a2ae208SSatish Balay #undef __FUNCT__ 7494a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 7508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 751273d9f13SBarry Smith { 752273d9f13SBarry Smith Mat_MPIAdj *b; 7536849ba73SBarry Smith PetscErrorCode ierr; 754273d9f13SBarry Smith 755273d9f13SBarry Smith PetscFunctionBegin; 756b00a9115SJed Brown ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 757b0a32e0cSBarry Smith B->data = (void*)b; 758273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 759273d9f13SBarry Smith B->assembled = PETSC_FALSE; 760273d9f13SBarry Smith 761bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 762bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 76317667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 764273d9f13SBarry Smith PetscFunctionReturn(0); 765273d9f13SBarry Smith } 766273d9f13SBarry Smith 7674a2ae208SSatish Balay #undef __FUNCT__ 7684a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 769a23d5eceSKris Buschelman /*@C 770a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 771a23d5eceSKris Buschelman 7723f9fe445SBarry Smith Logically Collective on MPI_Comm 773a23d5eceSKris Buschelman 774a23d5eceSKris Buschelman Input Parameters: 775a23d5eceSKris Buschelman + A - the matrix 776a23d5eceSKris Buschelman . i - the indices into j for the start of each row 777a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 778a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 779a23d5eceSKris Buschelman - values - [optional] edge weights 780a23d5eceSKris Buschelman 781a23d5eceSKris Buschelman Level: intermediate 782a23d5eceSKris Buschelman 783a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 784a23d5eceSKris Buschelman @*/ 7857087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 786273d9f13SBarry Smith { 7874ac538c5SBarry Smith PetscErrorCode ierr; 788273d9f13SBarry Smith 789273d9f13SBarry Smith PetscFunctionBegin; 7904ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 791273d9f13SBarry Smith PetscFunctionReturn(0); 792273d9f13SBarry Smith } 793273d9f13SBarry Smith 7944a2ae208SSatish Balay #undef __FUNCT__ 7954a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 796b97cf49bSBarry Smith /*@C 7973eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 798b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 799b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 800b97cf49bSBarry Smith 801ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 802ef5ee4d1SLois Curfman McInnes 803b97cf49bSBarry Smith Input Parameters: 804c2e958feSBarry Smith + comm - MPI communicator 8050752156aSBarry Smith . m - number of local rows 80658ec18a5SLisandro Dalcin . N - number of global columns 807b97cf49bSBarry Smith . i - the indices into j for the start of each row 808329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 809ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 810329f5518SBarry Smith - values -[optional] edge weights 811b97cf49bSBarry Smith 812b97cf49bSBarry Smith Output Parameter: 813b97cf49bSBarry Smith . A - the matrix 814b97cf49bSBarry Smith 81515091d37SBarry Smith Level: intermediate 81615091d37SBarry Smith 8174bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 8184bc6d8c0SBarry Smith MatSetValues(). 819c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 820ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 8211198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 822327e1a73SBarry Smith Should not include the matrix diagonals. 823b97cf49bSBarry Smith 824e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 825ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 826ededeb1aSvictorle 827ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 828b97cf49bSBarry Smith 829e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 830b97cf49bSBarry Smith @*/ 8317087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 832b97cf49bSBarry Smith { 833dfbe8321SBarry Smith PetscErrorCode ierr; 834b97cf49bSBarry Smith 835433994e6SBarry Smith PetscFunctionBegin; 836f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 83758ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 838273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 839273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 8403a40ed3dSBarry Smith PetscFunctionReturn(0); 841b97cf49bSBarry Smith } 842b97cf49bSBarry Smith 843c2e958feSBarry Smith 844433994e6SBarry Smith 845433994e6SBarry Smith 846433994e6SBarry Smith 847433994e6SBarry Smith 848433994e6SBarry Smith 849