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) { 23ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)A),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; 49251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5032077d6dSBarry Smith if (iascii) { 513eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 52b97cf49bSBarry Smith } 533a40ed3dSBarry Smith PetscFunctionReturn(0); 54b97cf49bSBarry Smith } 55b97cf49bSBarry Smith 564a2ae208SSatish Balay #undef __FUNCT__ 574a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 58dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat) 59b97cf49bSBarry Smith { 603eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 61dfbe8321SBarry Smith PetscErrorCode ierr; 6294d884c6SBarry Smith 6394d884c6SBarry Smith PetscFunctionBegin; 64aa482453SBarry Smith #if defined(PETSC_USE_LOG) 65d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz); 66b97cf49bSBarry Smith #endif 6705b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 680400557aSBarry Smith if (a->freeaij) { 6914196391SBarry Smith if (a->freeaijwithfree) { 7014196391SBarry Smith if (a->i) free(a->i); 7114196391SBarry Smith if (a->j) free(a->j); 7214196391SBarry Smith } else { 73606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 74606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 7505b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 7614196391SBarry Smith } 770400557aSBarry Smith } 78bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 79dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 80bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr); 81bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr); 823a40ed3dSBarry Smith PetscFunctionReturn(0); 83b97cf49bSBarry Smith } 84b97cf49bSBarry Smith 854a2ae208SSatish Balay #undef __FUNCT__ 864a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 87ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 88b97cf49bSBarry Smith { 893eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 9063ba0a88SBarry Smith PetscErrorCode ierr; 91b97cf49bSBarry Smith 92433994e6SBarry Smith PetscFunctionBegin; 9312c028f9SKris Buschelman switch (op) { 9477e54ba9SKris Buschelman case MAT_SYMMETRIC: 9512c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 969a4540c5SBarry Smith case MAT_HERMITIAN: 974e0d8c25SBarry Smith a->symmetric = flg; 989a4540c5SBarry Smith break; 999a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1009a4540c5SBarry Smith break; 10112c028f9SKris Buschelman default: 102290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 10312c028f9SKris Buschelman break; 104b97cf49bSBarry Smith } 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 106b97cf49bSBarry Smith } 107b97cf49bSBarry Smith 108b97cf49bSBarry Smith 109b97cf49bSBarry Smith /* 110b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 111b97cf49bSBarry Smith */ 112b97cf49bSBarry Smith 1134a2ae208SSatish Balay #undef __FUNCT__ 1144a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 115dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A) 116b97cf49bSBarry Smith { 1173eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 1186849ba73SBarry Smith PetscErrorCode ierr; 119d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 120b97cf49bSBarry Smith 121433994e6SBarry Smith PetscFunctionBegin; 12209f38230SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 12309f38230SBarry Smith ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr); 124d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 125b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 126b97cf49bSBarry Smith if (a->j[j] == i) { 12709f38230SBarry Smith a->diag[i] = j; 128b97cf49bSBarry Smith break; 129b97cf49bSBarry Smith } 130b97cf49bSBarry Smith } 131b97cf49bSBarry Smith } 1323a40ed3dSBarry Smith PetscFunctionReturn(0); 133b97cf49bSBarry Smith } 134b97cf49bSBarry Smith 1354a2ae208SSatish Balay #undef __FUNCT__ 1364a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 137b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 138b97cf49bSBarry Smith { 1393eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 140b97cf49bSBarry Smith 141433994e6SBarry Smith PetscFunctionBegin; 142d0f46423SBarry Smith row -= A->rmap->rstart; 1430752156aSBarry Smith 144e32f2f54SBarry Smith if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 145b97cf49bSBarry Smith 146b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 1470298fd71SBarry Smith if (v) *v = NULL; 148*92b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 1493a40ed3dSBarry Smith PetscFunctionReturn(0); 150b97cf49bSBarry Smith } 151b97cf49bSBarry Smith 1524a2ae208SSatish Balay #undef __FUNCT__ 1534a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 154b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 155b97cf49bSBarry Smith { 156433994e6SBarry Smith PetscFunctionBegin; 1573a40ed3dSBarry Smith PetscFunctionReturn(0); 158b97cf49bSBarry Smith } 159b97cf49bSBarry Smith 1604a2ae208SSatish Balay #undef __FUNCT__ 1614a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 162ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 163b97cf49bSBarry Smith { 1643eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 165dfbe8321SBarry Smith PetscErrorCode ierr; 166ace3abfcSBarry Smith PetscBool flag; 167b97cf49bSBarry Smith 168433994e6SBarry Smith PetscFunctionBegin; 169b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 170d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 1710f5bd95cSBarry Smith flag = PETSC_FALSE; 172b97cf49bSBarry Smith } 173b97cf49bSBarry Smith 174b97cf49bSBarry Smith /* if the a->i are the same */ 175d0f46423SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 176b97cf49bSBarry Smith 177b97cf49bSBarry Smith /* if a->j are the same */ 178b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 179b97cf49bSBarry Smith 180ce94432eSBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 1813a40ed3dSBarry Smith PetscFunctionReturn(0); 182b97cf49bSBarry Smith } 183b97cf49bSBarry Smith 1844a2ae208SSatish Balay #undef __FUNCT__ 1854a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 1861a83f524SJed Brown PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 1876cd91f64SBarry Smith { 188b24ad042SBarry Smith PetscInt i; 1896cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 1901a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 1916cd91f64SBarry Smith 1926cd91f64SBarry Smith PetscFunctionBegin; 193d0f46423SBarry Smith *m = A->rmap->n; 1946cd91f64SBarry Smith *ia = a->i; 1956cd91f64SBarry Smith *ja = a->j; 1966cd91f64SBarry Smith *done = PETSC_TRUE; 197d42735a1SBarry Smith if (oshift) { 198d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 199d42735a1SBarry Smith (*ja)[i]++; 200d42735a1SBarry Smith } 201d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 202d42735a1SBarry Smith } 203d42735a1SBarry Smith PetscFunctionReturn(0); 204d42735a1SBarry Smith } 205d42735a1SBarry Smith 2064a2ae208SSatish Balay #undef __FUNCT__ 2074a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 2081a83f524SJed Brown PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 209d42735a1SBarry Smith { 210b24ad042SBarry Smith PetscInt i; 211d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 2121a83f524SJed Brown PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 213d42735a1SBarry Smith 214d42735a1SBarry Smith PetscFunctionBegin; 215e32f2f54SBarry Smith if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 216e32f2f54SBarry Smith if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 217d42735a1SBarry Smith if (oshift) { 218d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 219d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 220d42735a1SBarry Smith (*ja)[i]--; 221d42735a1SBarry Smith } 222d42735a1SBarry Smith } 2236cd91f64SBarry Smith PetscFunctionReturn(0); 2246cd91f64SBarry Smith } 225b97cf49bSBarry Smith 22617667f90SBarry Smith #undef __FUNCT__ 22717667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj" 22819fd82e9SBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 22917667f90SBarry Smith { 23017667f90SBarry Smith Mat B; 23117667f90SBarry Smith PetscErrorCode ierr; 23217667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 23317667f90SBarry Smith const PetscInt *rj; 23417667f90SBarry Smith const PetscScalar *ra; 23517667f90SBarry Smith MPI_Comm comm; 23617667f90SBarry Smith 23717667f90SBarry Smith PetscFunctionBegin; 2380298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 2390298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 2400298fd71SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 24117667f90SBarry Smith 24217667f90SBarry Smith /* count the number of nonzeros per row */ 24317667f90SBarry Smith for (i=0; i<m; i++) { 2440298fd71SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 2455ee9ba1cSJed Brown for (j=0; j<len; j++) { 2465ee9ba1cSJed Brown if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 2475ee9ba1cSJed Brown } 24817667f90SBarry Smith nzeros += len; 2490f2063bfSTobin Isaac ierr = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 25017667f90SBarry Smith } 25117667f90SBarry Smith 25217667f90SBarry Smith /* malloc space for nonzeros */ 25317667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr); 25417667f90SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr); 25517667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr); 25617667f90SBarry Smith 25717667f90SBarry Smith nzeros = 0; 25817667f90SBarry Smith ia[0] = 0; 25917667f90SBarry Smith for (i=0; i<m; i++) { 26017667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 26117667f90SBarry Smith cnt = 0; 26217667f90SBarry Smith for (j=0; j<len; j++) { 2635ee9ba1cSJed Brown if (rj[j] != i+rstart) { /* if not diagonal */ 26417667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 26517667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 26617667f90SBarry Smith } 2675ee9ba1cSJed Brown } 26817667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 26917667f90SBarry Smith nzeros += cnt; 27017667f90SBarry Smith ia[i+1] = nzeros; 27117667f90SBarry Smith } 27217667f90SBarry Smith 27317667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 27417667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 27558ec18a5SLisandro Dalcin ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 27617667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 27717667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 27817667f90SBarry Smith 27917667f90SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 28017667f90SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 28117667f90SBarry Smith } else { 28217667f90SBarry Smith *newmat = B; 28317667f90SBarry Smith } 28417667f90SBarry Smith PetscFunctionReturn(0); 28517667f90SBarry Smith } 28617667f90SBarry Smith 287b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2880b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2893eda8832SBarry Smith MatGetRow_MPIAdj, 2903eda8832SBarry Smith MatRestoreRow_MPIAdj, 291b97cf49bSBarry Smith 0, 29297304618SKris Buschelman /* 4*/ 0, 293b97cf49bSBarry Smith 0, 294b97cf49bSBarry Smith 0, 295b97cf49bSBarry Smith 0, 2960b3b1487SBarry Smith 0, 2970b3b1487SBarry Smith 0, 29897304618SKris Buschelman /*10*/ 0, 2990b3b1487SBarry Smith 0, 3000b3b1487SBarry Smith 0, 3010b3b1487SBarry Smith 0, 3020b3b1487SBarry Smith 0, 30397304618SKris Buschelman /*15*/ 0, 3043eda8832SBarry Smith MatEqual_MPIAdj, 3050b3b1487SBarry Smith 0, 3060b3b1487SBarry Smith 0, 3070b3b1487SBarry Smith 0, 30897304618SKris Buschelman /*20*/ 0, 3090b3b1487SBarry Smith 0, 3103eda8832SBarry Smith MatSetOption_MPIAdj, 3110b3b1487SBarry Smith 0, 312d519adbfSMatthew Knepley /*24*/ 0, 3130b3b1487SBarry Smith 0, 3140b3b1487SBarry Smith 0, 3150b3b1487SBarry Smith 0, 3160b3b1487SBarry Smith 0, 317d519adbfSMatthew Knepley /*29*/ 0, 3180b3b1487SBarry Smith 0, 319273d9f13SBarry Smith 0, 320273d9f13SBarry Smith 0, 3210b3b1487SBarry Smith 0, 322d519adbfSMatthew Knepley /*34*/ 0, 3230b3b1487SBarry Smith 0, 3240b3b1487SBarry Smith 0, 3250b3b1487SBarry Smith 0, 3260b3b1487SBarry Smith 0, 327d519adbfSMatthew Knepley /*39*/ 0, 3280b3b1487SBarry Smith 0, 3290b3b1487SBarry Smith 0, 3300b3b1487SBarry Smith 0, 3310b3b1487SBarry Smith 0, 332d519adbfSMatthew Knepley /*44*/ 0, 3330b3b1487SBarry Smith 0, 3340b3b1487SBarry Smith 0, 3350b3b1487SBarry Smith 0, 3360b3b1487SBarry Smith 0, 337d519adbfSMatthew Knepley /*49*/ 0, 3386cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 339d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 340b97cf49bSBarry Smith 0, 341b97cf49bSBarry Smith 0, 342d519adbfSMatthew Knepley /*54*/ 0, 343b97cf49bSBarry Smith 0, 3440752156aSBarry Smith 0, 3450752156aSBarry Smith 0, 3460b3b1487SBarry Smith 0, 347d519adbfSMatthew Knepley /*59*/ 0, 348b9b97703SBarry Smith MatDestroy_MPIAdj, 349b9b97703SBarry Smith MatView_MPIAdj, 35017667f90SBarry Smith MatConvertFrom_MPIAdj, 35197304618SKris Buschelman 0, 352d519adbfSMatthew Knepley /*64*/ 0, 35397304618SKris Buschelman 0, 35497304618SKris Buschelman 0, 35597304618SKris Buschelman 0, 35697304618SKris Buschelman 0, 357d519adbfSMatthew Knepley /*69*/ 0, 35897304618SKris Buschelman 0, 35997304618SKris Buschelman 0, 36097304618SKris Buschelman 0, 36197304618SKris Buschelman 0, 362d519adbfSMatthew Knepley /*74*/ 0, 36397304618SKris Buschelman 0, 36497304618SKris Buschelman 0, 36597304618SKris Buschelman 0, 36697304618SKris Buschelman 0, 367d519adbfSMatthew Knepley /*79*/ 0, 36897304618SKris Buschelman 0, 36997304618SKris Buschelman 0, 37097304618SKris Buschelman 0, 371865e5f61SKris Buschelman 0, 372d519adbfSMatthew Knepley /*84*/ 0, 373865e5f61SKris Buschelman 0, 374865e5f61SKris Buschelman 0, 375865e5f61SKris Buschelman 0, 376865e5f61SKris Buschelman 0, 377d519adbfSMatthew Knepley /*89*/ 0, 378865e5f61SKris Buschelman 0, 379865e5f61SKris Buschelman 0, 380865e5f61SKris Buschelman 0, 381865e5f61SKris Buschelman 0, 382d519adbfSMatthew Knepley /*94*/ 0, 383865e5f61SKris Buschelman 0, 384865e5f61SKris Buschelman 0, 3853964eb88SJed Brown 0, 3863964eb88SJed Brown 0, 3873964eb88SJed Brown /*99*/ 0, 3883964eb88SJed Brown 0, 3893964eb88SJed Brown 0, 3903964eb88SJed Brown 0, 3913964eb88SJed Brown 0, 3923964eb88SJed Brown /*104*/ 0, 3933964eb88SJed Brown 0, 3943964eb88SJed Brown 0, 3953964eb88SJed Brown 0, 3963964eb88SJed Brown 0, 3973964eb88SJed Brown /*109*/ 0, 3983964eb88SJed Brown 0, 3993964eb88SJed Brown 0, 4003964eb88SJed Brown 0, 4013964eb88SJed Brown 0, 4023964eb88SJed Brown /*114*/ 0, 4033964eb88SJed Brown 0, 4043964eb88SJed Brown 0, 4053964eb88SJed Brown 0, 4063964eb88SJed Brown 0, 4073964eb88SJed Brown /*119*/ 0, 4083964eb88SJed Brown 0, 4093964eb88SJed Brown 0, 4103964eb88SJed Brown 0, 4113964eb88SJed Brown 0, 4123964eb88SJed Brown /*124*/ 0, 4133964eb88SJed Brown 0, 4143964eb88SJed Brown 0, 4153964eb88SJed Brown 0, 4163964eb88SJed Brown 0, 4173964eb88SJed Brown /*129*/ 0, 4183964eb88SJed Brown 0, 4193964eb88SJed Brown 0, 4203964eb88SJed Brown 0, 4213964eb88SJed Brown 0, 4223964eb88SJed Brown /*134*/ 0, 4233964eb88SJed Brown 0, 4243964eb88SJed Brown 0, 4253964eb88SJed Brown 0, 4263964eb88SJed Brown 0, 4273964eb88SJed Brown /*139*/ 0, 4283964eb88SJed Brown 0 4293964eb88SJed Brown }; 430b97cf49bSBarry Smith 431a23d5eceSKris Buschelman #undef __FUNCT__ 432a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 433f7a08781SBarry Smith static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 434a23d5eceSKris Buschelman { 435a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 436dfbe8321SBarry Smith PetscErrorCode ierr; 4372515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 438b24ad042SBarry Smith PetscInt ii; 439a23d5eceSKris Buschelman #endif 440a23d5eceSKris Buschelman 441a23d5eceSKris Buschelman PetscFunctionBegin; 44226283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 44326283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 44458ec18a5SLisandro Dalcin 4452515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 446e32f2f54SBarry 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]); 447d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 448e02043d6SBarry 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]); 449a23d5eceSKris Buschelman } 450d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 451e02043d6SBarry 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]); 452a23d5eceSKris Buschelman } 453a23d5eceSKris Buschelman #endif 45458ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 455a23d5eceSKris Buschelman 456a23d5eceSKris Buschelman b->j = j; 457a23d5eceSKris Buschelman b->i = i; 458a23d5eceSKris Buschelman b->values = values; 459a23d5eceSKris Buschelman 460d0f46423SBarry Smith b->nz = i[B->rmap->n]; 461a23d5eceSKris Buschelman b->diag = 0; 462a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 463a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 464a23d5eceSKris Buschelman 465a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 466a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 467a23d5eceSKris Buschelman PetscFunctionReturn(0); 468a23d5eceSKris Buschelman } 469b97cf49bSBarry Smith 4709a3665c6SJed Brown #undef __FUNCT__ 4719a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj" 472f7a08781SBarry Smith static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 4739a3665c6SJed Brown { 4749a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 4759a3665c6SJed Brown PetscErrorCode ierr; 4769a3665c6SJed Brown const PetscInt *ranges; 4779a3665c6SJed Brown MPI_Comm acomm,bcomm; 4789a3665c6SJed Brown MPI_Group agroup,bgroup; 4799a3665c6SJed Brown PetscMPIInt i,rank,size,nranks,*ranks; 4809a3665c6SJed Brown 4819a3665c6SJed Brown PetscFunctionBegin; 4820298fd71SBarry Smith *B = NULL; 483ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr); 4849a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr); 4859a3665c6SJed Brown ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr); 4869a3665c6SJed Brown ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 4879a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 4889a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) nranks++; 4899a3665c6SJed Brown } 4909a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 4919a3665c6SJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 4929a3665c6SJed Brown *B = A; 4939a3665c6SJed Brown PetscFunctionReturn(0); 4949a3665c6SJed Brown } 4959a3665c6SJed Brown 4969a3665c6SJed Brown ierr = PetscMalloc(nranks*sizeof(PetscMPIInt),&ranks);CHKERRQ(ierr); 4979a3665c6SJed Brown for (i=0,nranks=0; i<size; i++) { 4989a3665c6SJed Brown if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 4999a3665c6SJed Brown } 5009a3665c6SJed Brown ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr); 5019a3665c6SJed Brown ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr); 5029a3665c6SJed Brown ierr = PetscFree(ranks);CHKERRQ(ierr); 5039a3665c6SJed Brown ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr); 5049a3665c6SJed Brown ierr = MPI_Group_free(&agroup);CHKERRQ(ierr); 5059a3665c6SJed Brown ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr); 5069a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 5079a3665c6SJed Brown PetscInt m,N; 5089a3665c6SJed Brown Mat_MPIAdj *b; 5090298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 5100298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 5119a3665c6SJed Brown ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 5129a3665c6SJed Brown b = (Mat_MPIAdj*)(*B)->data; 5139a3665c6SJed Brown b->freeaij = PETSC_FALSE; 5149a3665c6SJed Brown ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr); 5159a3665c6SJed Brown } 5169a3665c6SJed Brown PetscFunctionReturn(0); 5179a3665c6SJed Brown } 5189a3665c6SJed Brown 5199a3665c6SJed Brown #undef __FUNCT__ 5209a3665c6SJed Brown #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat" 5219a3665c6SJed Brown /*@ 5229a3665c6SJed Brown MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 5239a3665c6SJed Brown 5249a3665c6SJed Brown Collective 5259a3665c6SJed Brown 5269a3665c6SJed Brown Input Arguments: 5279a3665c6SJed Brown . A - original MPIAdj matrix 5289a3665c6SJed Brown 5299a3665c6SJed Brown Output Arguments: 5300298fd71SBarry Smith . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 5319a3665c6SJed Brown 5329a3665c6SJed Brown Level: developer 5339a3665c6SJed Brown 5349a3665c6SJed Brown Note: 5359a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 5369a3665c6SJed Brown 5379a3665c6SJed Brown The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 5389a3665c6SJed Brown 5399a3665c6SJed Brown .seealso: MatCreateMPIAdj() 5409a3665c6SJed Brown @*/ 5419a3665c6SJed Brown PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 5429a3665c6SJed Brown { 5439a3665c6SJed Brown PetscErrorCode ierr; 5449a3665c6SJed Brown 5459a3665c6SJed Brown PetscFunctionBegin; 5469a3665c6SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 5479a3665c6SJed Brown ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 5489a3665c6SJed Brown PetscFunctionReturn(0); 5499a3665c6SJed Brown } 5509a3665c6SJed Brown 5510bad9183SKris Buschelman /*MC 552fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 5530bad9183SKris Buschelman intended for use constructing orderings and partitionings. 5540bad9183SKris Buschelman 5550bad9183SKris Buschelman Level: beginner 5560bad9183SKris Buschelman 5570bad9183SKris Buschelman .seealso: MatCreateMPIAdj 5580bad9183SKris Buschelman M*/ 5590bad9183SKris Buschelman 5604a2ae208SSatish Balay #undef __FUNCT__ 5614a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 5628cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 563273d9f13SBarry Smith { 564273d9f13SBarry Smith Mat_MPIAdj *b; 5656849ba73SBarry Smith PetscErrorCode ierr; 566273d9f13SBarry Smith 567273d9f13SBarry Smith PetscFunctionBegin; 56838f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr); 569b0a32e0cSBarry Smith B->data = (void*)b; 570273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 571273d9f13SBarry Smith B->assembled = PETSC_FALSE; 572273d9f13SBarry Smith 573bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 574bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 57517667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 576273d9f13SBarry Smith PetscFunctionReturn(0); 577273d9f13SBarry Smith } 578273d9f13SBarry Smith 5794a2ae208SSatish Balay #undef __FUNCT__ 5804a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 581a23d5eceSKris Buschelman /*@C 582a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 583a23d5eceSKris Buschelman 5843f9fe445SBarry Smith Logically Collective on MPI_Comm 585a23d5eceSKris Buschelman 586a23d5eceSKris Buschelman Input Parameters: 587a23d5eceSKris Buschelman + A - the matrix 588a23d5eceSKris Buschelman . i - the indices into j for the start of each row 589a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 590a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 591a23d5eceSKris Buschelman - values - [optional] edge weights 592a23d5eceSKris Buschelman 593a23d5eceSKris Buschelman Level: intermediate 594a23d5eceSKris Buschelman 595a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 596a23d5eceSKris Buschelman @*/ 5977087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 598273d9f13SBarry Smith { 5994ac538c5SBarry Smith PetscErrorCode ierr; 600273d9f13SBarry Smith 601273d9f13SBarry Smith PetscFunctionBegin; 6024ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 603273d9f13SBarry Smith PetscFunctionReturn(0); 604273d9f13SBarry Smith } 605273d9f13SBarry Smith 6064a2ae208SSatish Balay #undef __FUNCT__ 6074a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 608b97cf49bSBarry Smith /*@C 6093eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 610b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 611b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 612b97cf49bSBarry Smith 613ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 614ef5ee4d1SLois Curfman McInnes 615b97cf49bSBarry Smith Input Parameters: 616c2e958feSBarry Smith + comm - MPI communicator 6170752156aSBarry Smith . m - number of local rows 61858ec18a5SLisandro Dalcin . N - number of global columns 619b97cf49bSBarry Smith . i - the indices into j for the start of each row 620329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 621ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 622329f5518SBarry Smith - values -[optional] edge weights 623b97cf49bSBarry Smith 624b97cf49bSBarry Smith Output Parameter: 625b97cf49bSBarry Smith . A - the matrix 626b97cf49bSBarry Smith 62715091d37SBarry Smith Level: intermediate 62815091d37SBarry Smith 6294bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 6304bc6d8c0SBarry Smith MatSetValues(). 631c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 632ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 6331198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 634327e1a73SBarry Smith Should not include the matrix diagonals. 635b97cf49bSBarry Smith 636e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 637ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 638ededeb1aSvictorle 639ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 640b97cf49bSBarry Smith 641e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 642b97cf49bSBarry Smith @*/ 6437087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 644b97cf49bSBarry Smith { 645dfbe8321SBarry Smith PetscErrorCode ierr; 646b97cf49bSBarry Smith 647433994e6SBarry Smith PetscFunctionBegin; 648f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 64958ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 650273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 651273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 6523a40ed3dSBarry Smith PetscFunctionReturn(0); 653b97cf49bSBarry Smith } 654b97cf49bSBarry Smith 655c2e958feSBarry Smith 656433994e6SBarry Smith 657433994e6SBarry Smith 658433994e6SBarry Smith 659433994e6SBarry Smith 660433994e6SBarry Smith 661