1be1d678aSKris Buschelman 2b97cf49bSBarry Smith /* 3b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 4b97cf49bSBarry Smith */ 5*c6db04a5SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h> 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); 26b97cf49bSBarry Smith for (i=0; i<m; i++) { 27d0f46423SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr); 28b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 2977431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 300752156aSBarry Smith } 31b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 32b97cf49bSBarry Smith } 33d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 34b97cf49bSBarry Smith } 35b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 363a40ed3dSBarry Smith PetscFunctionReturn(0); 37b97cf49bSBarry Smith } 38b97cf49bSBarry Smith 394a2ae208SSatish Balay #undef __FUNCT__ 404a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj" 41dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 42b97cf49bSBarry Smith { 43dfbe8321SBarry Smith PetscErrorCode ierr; 44ace3abfcSBarry Smith PetscBool iascii; 45b97cf49bSBarry Smith 46433994e6SBarry Smith PetscFunctionBegin; 472692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4832077d6dSBarry Smith if (iascii) { 493eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 505cd90555SBarry Smith } else { 51e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 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 } 781198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 79dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 80901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 813a40ed3dSBarry Smith PetscFunctionReturn(0); 82b97cf49bSBarry Smith } 83b97cf49bSBarry Smith 844a2ae208SSatish Balay #undef __FUNCT__ 854a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 86ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 87b97cf49bSBarry Smith { 883eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 8963ba0a88SBarry Smith PetscErrorCode ierr; 90b97cf49bSBarry Smith 91433994e6SBarry Smith PetscFunctionBegin; 9212c028f9SKris Buschelman switch (op) { 9377e54ba9SKris Buschelman case MAT_SYMMETRIC: 9412c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 959a4540c5SBarry Smith case MAT_HERMITIAN: 964e0d8c25SBarry Smith a->symmetric = flg; 979a4540c5SBarry Smith break; 989a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 999a4540c5SBarry Smith break; 10012c028f9SKris Buschelman default: 101290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 10212c028f9SKris Buschelman break; 103b97cf49bSBarry Smith } 1043a40ed3dSBarry Smith PetscFunctionReturn(0); 105b97cf49bSBarry Smith } 106b97cf49bSBarry Smith 107b97cf49bSBarry Smith 108b97cf49bSBarry Smith /* 109b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 110b97cf49bSBarry Smith */ 111b97cf49bSBarry Smith 1124a2ae208SSatish Balay #undef __FUNCT__ 1134a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 114dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A) 115b97cf49bSBarry Smith { 1163eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 1176849ba73SBarry Smith PetscErrorCode ierr; 118d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n; 119b97cf49bSBarry Smith 120433994e6SBarry Smith PetscFunctionBegin; 12109f38230SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 12209f38230SBarry Smith ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr); 123d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 124b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 125b97cf49bSBarry Smith if (a->j[j] == i) { 12609f38230SBarry Smith a->diag[i] = j; 127b97cf49bSBarry Smith break; 128b97cf49bSBarry Smith } 129b97cf49bSBarry Smith } 130b97cf49bSBarry Smith } 1313a40ed3dSBarry Smith PetscFunctionReturn(0); 132b97cf49bSBarry Smith } 133b97cf49bSBarry Smith 1344a2ae208SSatish Balay #undef __FUNCT__ 1354a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 136b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 137b97cf49bSBarry Smith { 1383eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 139b24ad042SBarry Smith PetscInt *itmp; 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]; 147b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 148b97cf49bSBarry Smith if (idx) { 149b97cf49bSBarry Smith itmp = a->j + a->i[row]; 150b97cf49bSBarry Smith if (*nz) { 151b97cf49bSBarry Smith *idx = itmp; 152b97cf49bSBarry Smith } 153b97cf49bSBarry Smith else *idx = 0; 154b97cf49bSBarry Smith } 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 156b97cf49bSBarry Smith } 157b97cf49bSBarry Smith 1584a2ae208SSatish Balay #undef __FUNCT__ 1594a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 160b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 161b97cf49bSBarry Smith { 162433994e6SBarry Smith PetscFunctionBegin; 1633a40ed3dSBarry Smith PetscFunctionReturn(0); 164b97cf49bSBarry Smith } 165b97cf49bSBarry Smith 1664a2ae208SSatish Balay #undef __FUNCT__ 1674a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 168ace3abfcSBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 169b97cf49bSBarry Smith { 1703eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 171dfbe8321SBarry Smith PetscErrorCode ierr; 172ace3abfcSBarry Smith PetscBool flag; 173b97cf49bSBarry Smith 174433994e6SBarry Smith PetscFunctionBegin; 175b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 176d0f46423SBarry Smith if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 1770f5bd95cSBarry Smith flag = PETSC_FALSE; 178b97cf49bSBarry Smith } 179b97cf49bSBarry Smith 180b97cf49bSBarry Smith /* if the a->i are the same */ 181d0f46423SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 182b97cf49bSBarry Smith 183b97cf49bSBarry Smith /* if a->j are the same */ 184b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 185b97cf49bSBarry Smith 1867adad957SLisandro Dalcin ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);CHKERRQ(ierr); 1873a40ed3dSBarry Smith PetscFunctionReturn(0); 188b97cf49bSBarry Smith } 189b97cf49bSBarry Smith 1904a2ae208SSatish Balay #undef __FUNCT__ 1914a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 192ace3abfcSBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool *done) 1936cd91f64SBarry Smith { 194b24ad042SBarry Smith PetscInt i; 1956cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 1966cd91f64SBarry Smith 1976cd91f64SBarry Smith PetscFunctionBegin; 198d0f46423SBarry Smith *m = A->rmap->n; 1996cd91f64SBarry Smith *ia = a->i; 2006cd91f64SBarry Smith *ja = a->j; 2016cd91f64SBarry Smith *done = PETSC_TRUE; 202d42735a1SBarry Smith if (oshift) { 203d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 204d42735a1SBarry Smith (*ja)[i]++; 205d42735a1SBarry Smith } 206d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 207d42735a1SBarry Smith } 208d42735a1SBarry Smith PetscFunctionReturn(0); 209d42735a1SBarry Smith } 210d42735a1SBarry Smith 2114a2ae208SSatish Balay #undef __FUNCT__ 2124a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 213ace3abfcSBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscBool *done) 214d42735a1SBarry Smith { 215b24ad042SBarry Smith PetscInt i; 216d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 217d42735a1SBarry Smith 218d42735a1SBarry Smith PetscFunctionBegin; 219e32f2f54SBarry Smith if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 220e32f2f54SBarry Smith if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 221d42735a1SBarry Smith if (oshift) { 222d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 223d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 224d42735a1SBarry Smith (*ja)[i]--; 225d42735a1SBarry Smith } 226d42735a1SBarry Smith } 2276cd91f64SBarry Smith PetscFunctionReturn(0); 2286cd91f64SBarry Smith } 229b97cf49bSBarry Smith 23017667f90SBarry Smith #undef __FUNCT__ 23117667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj" 2327087cfbeSBarry Smith PetscErrorCode MatConvertFrom_MPIAdj(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 23317667f90SBarry Smith { 23417667f90SBarry Smith Mat B; 23517667f90SBarry Smith PetscErrorCode ierr; 23617667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 23717667f90SBarry Smith const PetscInt *rj; 23817667f90SBarry Smith const PetscScalar *ra; 23917667f90SBarry Smith MPI_Comm comm; 24017667f90SBarry Smith 24117667f90SBarry Smith PetscFunctionBegin; 24217667f90SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 24317667f90SBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 24417667f90SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 24517667f90SBarry Smith 24617667f90SBarry Smith /* count the number of nonzeros per row */ 24717667f90SBarry Smith for (i=0; i<m; i++) { 24817667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 24917667f90SBarry Smith for (j=0; j<len; j++) { 25017667f90SBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 25117667f90SBarry Smith } 25217667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 25317667f90SBarry Smith nzeros += len; 25417667f90SBarry Smith } 25517667f90SBarry Smith 25617667f90SBarry Smith /* malloc space for nonzeros */ 25717667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr); 25817667f90SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr); 25917667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr); 26017667f90SBarry Smith 26117667f90SBarry Smith nzeros = 0; 26217667f90SBarry Smith ia[0] = 0; 26317667f90SBarry Smith for (i=0; i<m; i++) { 26417667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 26517667f90SBarry Smith cnt = 0; 26617667f90SBarry Smith for (j=0; j<len; j++) { 26717667f90SBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 26817667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 26917667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 27017667f90SBarry Smith } 27117667f90SBarry Smith } 27217667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 27317667f90SBarry Smith nzeros += cnt; 27417667f90SBarry Smith ia[i+1] = nzeros; 27517667f90SBarry Smith } 27617667f90SBarry Smith 27717667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 27817667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 27958ec18a5SLisandro Dalcin ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 28017667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 28117667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 28217667f90SBarry Smith 28317667f90SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 28417667f90SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 28517667f90SBarry Smith } else { 28617667f90SBarry Smith *newmat = B; 28717667f90SBarry Smith } 28817667f90SBarry Smith PetscFunctionReturn(0); 28917667f90SBarry Smith } 29017667f90SBarry Smith 291b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2920b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2933eda8832SBarry Smith MatGetRow_MPIAdj, 2943eda8832SBarry Smith MatRestoreRow_MPIAdj, 295b97cf49bSBarry Smith 0, 29697304618SKris Buschelman /* 4*/ 0, 297b97cf49bSBarry Smith 0, 298b97cf49bSBarry Smith 0, 299b97cf49bSBarry Smith 0, 3000b3b1487SBarry Smith 0, 3010b3b1487SBarry Smith 0, 30297304618SKris Buschelman /*10*/ 0, 3030b3b1487SBarry Smith 0, 3040b3b1487SBarry Smith 0, 3050b3b1487SBarry Smith 0, 3060b3b1487SBarry Smith 0, 30797304618SKris Buschelman /*15*/ 0, 3083eda8832SBarry Smith MatEqual_MPIAdj, 3090b3b1487SBarry Smith 0, 3100b3b1487SBarry Smith 0, 3110b3b1487SBarry Smith 0, 31297304618SKris Buschelman /*20*/ 0, 3130b3b1487SBarry Smith 0, 3143eda8832SBarry Smith MatSetOption_MPIAdj, 3150b3b1487SBarry Smith 0, 316d519adbfSMatthew Knepley /*24*/ 0, 3170b3b1487SBarry Smith 0, 3180b3b1487SBarry Smith 0, 3190b3b1487SBarry Smith 0, 3200b3b1487SBarry Smith 0, 321d519adbfSMatthew Knepley /*29*/ 0, 3220b3b1487SBarry Smith 0, 323273d9f13SBarry Smith 0, 324273d9f13SBarry Smith 0, 3250b3b1487SBarry Smith 0, 326d519adbfSMatthew Knepley /*34*/ 0, 3270b3b1487SBarry Smith 0, 3280b3b1487SBarry Smith 0, 3290b3b1487SBarry Smith 0, 3300b3b1487SBarry Smith 0, 331d519adbfSMatthew Knepley /*39*/ 0, 3320b3b1487SBarry Smith 0, 3330b3b1487SBarry Smith 0, 3340b3b1487SBarry Smith 0, 3350b3b1487SBarry Smith 0, 336d519adbfSMatthew Knepley /*44*/ 0, 3370b3b1487SBarry Smith 0, 3380b3b1487SBarry Smith 0, 3390b3b1487SBarry Smith 0, 3400b3b1487SBarry Smith 0, 341d519adbfSMatthew Knepley /*49*/ 0, 3426cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 343d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 344b97cf49bSBarry Smith 0, 345b97cf49bSBarry Smith 0, 346d519adbfSMatthew Knepley /*54*/ 0, 347b97cf49bSBarry Smith 0, 3480752156aSBarry Smith 0, 3490752156aSBarry Smith 0, 3500b3b1487SBarry Smith 0, 351d519adbfSMatthew Knepley /*59*/ 0, 352b9b97703SBarry Smith MatDestroy_MPIAdj, 353b9b97703SBarry Smith MatView_MPIAdj, 35417667f90SBarry Smith MatConvertFrom_MPIAdj, 35597304618SKris Buschelman 0, 356d519adbfSMatthew Knepley /*64*/ 0, 35797304618SKris Buschelman 0, 35897304618SKris Buschelman 0, 35997304618SKris Buschelman 0, 36097304618SKris Buschelman 0, 361d519adbfSMatthew Knepley /*69*/ 0, 36297304618SKris Buschelman 0, 36397304618SKris Buschelman 0, 36497304618SKris Buschelman 0, 36597304618SKris Buschelman 0, 366d519adbfSMatthew Knepley /*74*/ 0, 36797304618SKris Buschelman 0, 36897304618SKris Buschelman 0, 36997304618SKris Buschelman 0, 37097304618SKris Buschelman 0, 371d519adbfSMatthew Knepley /*79*/ 0, 37297304618SKris Buschelman 0, 37397304618SKris Buschelman 0, 37497304618SKris Buschelman 0, 375865e5f61SKris Buschelman 0, 376d519adbfSMatthew Knepley /*84*/ 0, 377865e5f61SKris Buschelman 0, 378865e5f61SKris Buschelman 0, 379865e5f61SKris Buschelman 0, 380865e5f61SKris Buschelman 0, 381d519adbfSMatthew Knepley /*89*/ 0, 382865e5f61SKris Buschelman 0, 383865e5f61SKris Buschelman 0, 384865e5f61SKris Buschelman 0, 385865e5f61SKris Buschelman 0, 386d519adbfSMatthew Knepley /*94*/ 0, 387865e5f61SKris Buschelman 0, 388865e5f61SKris Buschelman 0, 389865e5f61SKris Buschelman 0}; 390b97cf49bSBarry Smith 391a23d5eceSKris Buschelman EXTERN_C_BEGIN 392a23d5eceSKris Buschelman #undef __FUNCT__ 393a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 3947087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 395a23d5eceSKris Buschelman { 396a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 397dfbe8321SBarry Smith PetscErrorCode ierr; 3982515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 399b24ad042SBarry Smith PetscInt ii; 400a23d5eceSKris Buschelman #endif 401a23d5eceSKris Buschelman 402a23d5eceSKris Buschelman PetscFunctionBegin; 40326283091SBarry Smith ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr); 40426283091SBarry Smith ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr); 40526283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 40626283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 40758ec18a5SLisandro Dalcin 4082515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 409e32f2f54SBarry 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]); 410d0f46423SBarry Smith for (ii=1; ii<B->rmap->n; ii++) { 411a23d5eceSKris Buschelman if (i[ii] < 0 || i[ii] < i[ii-1]) { 412e32f2f54SBarry Smith 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 } 414a23d5eceSKris Buschelman } 415d0f46423SBarry Smith for (ii=0; ii<i[B->rmap->n]; ii++) { 416d0f46423SBarry Smith if (j[ii] < 0 || j[ii] >= B->cmap->N) { 417e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]); 418a23d5eceSKris Buschelman } 419a23d5eceSKris Buschelman } 420a23d5eceSKris Buschelman #endif 42158ec18a5SLisandro Dalcin B->preallocated = PETSC_TRUE; 422a23d5eceSKris Buschelman 423a23d5eceSKris Buschelman b->j = j; 424a23d5eceSKris Buschelman b->i = i; 425a23d5eceSKris Buschelman b->values = values; 426a23d5eceSKris Buschelman 427d0f46423SBarry Smith b->nz = i[B->rmap->n]; 428a23d5eceSKris Buschelman b->diag = 0; 429a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 430a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 431a23d5eceSKris Buschelman 432a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 433a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 434a23d5eceSKris Buschelman PetscFunctionReturn(0); 435a23d5eceSKris Buschelman } 436a23d5eceSKris Buschelman EXTERN_C_END 437b97cf49bSBarry Smith 4380bad9183SKris Buschelman /*MC 439fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 4400bad9183SKris Buschelman intended for use constructing orderings and partitionings. 4410bad9183SKris Buschelman 4420bad9183SKris Buschelman Level: beginner 4430bad9183SKris Buschelman 4440bad9183SKris Buschelman .seealso: MatCreateMPIAdj 4450bad9183SKris Buschelman M*/ 4460bad9183SKris Buschelman 447273d9f13SBarry Smith EXTERN_C_BEGIN 4484a2ae208SSatish Balay #undef __FUNCT__ 4494a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 4507087cfbeSBarry Smith PetscErrorCode MatCreate_MPIAdj(Mat B) 451273d9f13SBarry Smith { 452273d9f13SBarry Smith Mat_MPIAdj *b; 4536849ba73SBarry Smith PetscErrorCode ierr; 454273d9f13SBarry Smith 455273d9f13SBarry Smith PetscFunctionBegin; 45638f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_MPIAdj,&b);CHKERRQ(ierr); 457b0a32e0cSBarry Smith B->data = (void*)b; 458273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 459273d9f13SBarry Smith B->assembled = PETSC_FALSE; 460273d9f13SBarry Smith 461a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C", 462a23d5eceSKris Buschelman "MatMPIAdjSetPreallocation_MPIAdj", 463a23d5eceSKris Buschelman MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 46417667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 465273d9f13SBarry Smith PetscFunctionReturn(0); 466273d9f13SBarry Smith } 467273d9f13SBarry Smith EXTERN_C_END 468273d9f13SBarry Smith 4694a2ae208SSatish Balay #undef __FUNCT__ 4704a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 471a23d5eceSKris Buschelman /*@C 472a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 473a23d5eceSKris Buschelman 4743f9fe445SBarry Smith Logically Collective on MPI_Comm 475a23d5eceSKris Buschelman 476a23d5eceSKris Buschelman Input Parameters: 477a23d5eceSKris Buschelman + A - the matrix 478a23d5eceSKris Buschelman . i - the indices into j for the start of each row 479a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 480a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 481a23d5eceSKris Buschelman - values - [optional] edge weights 482a23d5eceSKris Buschelman 483a23d5eceSKris Buschelman Level: intermediate 484a23d5eceSKris Buschelman 485a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 486a23d5eceSKris Buschelman @*/ 4877087cfbeSBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 488273d9f13SBarry Smith { 4894ac538c5SBarry Smith PetscErrorCode ierr; 490273d9f13SBarry Smith 491273d9f13SBarry Smith PetscFunctionBegin; 4924ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 493273d9f13SBarry Smith PetscFunctionReturn(0); 494273d9f13SBarry Smith } 495273d9f13SBarry Smith 4964a2ae208SSatish Balay #undef __FUNCT__ 4974a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 498b97cf49bSBarry Smith /*@C 4993eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 500b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 501b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 502b97cf49bSBarry Smith 503ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 504ef5ee4d1SLois Curfman McInnes 505b97cf49bSBarry Smith Input Parameters: 506c2e958feSBarry Smith + comm - MPI communicator 5070752156aSBarry Smith . m - number of local rows 50858ec18a5SLisandro Dalcin . N - number of global columns 509b97cf49bSBarry Smith . i - the indices into j for the start of each row 510329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 511ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 512329f5518SBarry Smith - values -[optional] edge weights 513b97cf49bSBarry Smith 514b97cf49bSBarry Smith Output Parameter: 515b97cf49bSBarry Smith . A - the matrix 516b97cf49bSBarry Smith 51715091d37SBarry Smith Level: intermediate 51815091d37SBarry Smith 5194bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 5204bc6d8c0SBarry Smith MatSetValues(). 521c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 522ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 5231198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 524327e1a73SBarry Smith Should not include the matrix diagonals. 525b97cf49bSBarry Smith 526e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 527ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 528ededeb1aSvictorle 529ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 530b97cf49bSBarry Smith 531e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 532b97cf49bSBarry Smith @*/ 5337087cfbeSBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 534b97cf49bSBarry Smith { 535dfbe8321SBarry Smith PetscErrorCode ierr; 536b97cf49bSBarry Smith 537433994e6SBarry Smith PetscFunctionBegin; 538f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 53958ec18a5SLisandro Dalcin ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 540273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 541273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 5423a40ed3dSBarry Smith PetscFunctionReturn(0); 543b97cf49bSBarry Smith } 544b97cf49bSBarry Smith 545c2e958feSBarry Smith 546433994e6SBarry Smith 547433994e6SBarry Smith 548433994e6SBarry Smith 549433994e6SBarry Smith 550433994e6SBarry Smith 551