1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3b97cf49bSBarry Smith /* 4b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 5b97cf49bSBarry Smith */ 63eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h" 7325e03aeSBarry Smith #include "petscsys.h" 8b97cf49bSBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII" 11dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 12b97cf49bSBarry Smith { 133eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 14dfbe8321SBarry Smith PetscErrorCode ierr; 15899cda47SBarry Smith PetscInt i,j,m = A->rmap.n; 162dcb1b2aSMatthew Knepley const char *name; 17ce1411ecSBarry Smith PetscViewerFormat format; 18b97cf49bSBarry Smith 19433994e6SBarry Smith PetscFunctionBegin; 203a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 21b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 22fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 233a40ed3dSBarry Smith PetscFunctionReturn(0); 24fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 25ffac6cdbSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 260752156aSBarry Smith } else { 27b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 28b97cf49bSBarry Smith for (i=0; i<m; i++) { 29899cda47SBarry 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 } 35b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 36b97cf49bSBarry Smith } 37b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 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; 4632077d6dSBarry Smith PetscTruth iascii; 47b97cf49bSBarry Smith 48433994e6SBarry Smith PetscFunctionBegin; 4932077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 5032077d6dSBarry Smith if (iascii) { 513eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 525cd90555SBarry Smith } else { 5379a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 54b97cf49bSBarry Smith } 553a40ed3dSBarry Smith PetscFunctionReturn(0); 56b97cf49bSBarry Smith } 57b97cf49bSBarry Smith 584a2ae208SSatish Balay #undef __FUNCT__ 594a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 60dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat) 61b97cf49bSBarry Smith { 623eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 63dfbe8321SBarry Smith PetscErrorCode ierr; 6494d884c6SBarry Smith 6594d884c6SBarry Smith PetscFunctionBegin; 66aa482453SBarry Smith #if defined(PETSC_USE_LOG) 67899cda47SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap.n,mat->cmap.n,a->nz); 68b97cf49bSBarry Smith #endif 6905b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 700400557aSBarry Smith if (a->freeaij) { 71606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 72606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 7305b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 740400557aSBarry Smith } 751198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 76901853e0SKris Buschelman 77901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 783a40ed3dSBarry Smith PetscFunctionReturn(0); 79b97cf49bSBarry Smith } 80b97cf49bSBarry Smith 814a2ae208SSatish Balay #undef __FUNCT__ 824a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 83dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op) 84b97cf49bSBarry Smith { 853eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 8663ba0a88SBarry Smith PetscErrorCode ierr; 87b97cf49bSBarry Smith 88433994e6SBarry Smith PetscFunctionBegin; 8912c028f9SKris Buschelman switch (op) { 9077e54ba9SKris Buschelman case MAT_SYMMETRIC: 9112c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 929a4540c5SBarry Smith case MAT_HERMITIAN: 93b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 9412c028f9SKris Buschelman break; 959a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 969a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 979a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 989a4540c5SBarry Smith a->symmetric = PETSC_FALSE; 999a4540c5SBarry Smith break; 1009a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1019a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 1029a4540c5SBarry Smith break; 10312c028f9SKris Buschelman default: 104*290bbb0aSBarry 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; 12109f38230SBarry Smith PetscInt i,j,m = A->rmap.n; 122b97cf49bSBarry Smith 123433994e6SBarry Smith PetscFunctionBegin; 12409f38230SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 12509f38230SBarry Smith ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr); 126899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 127b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 128b97cf49bSBarry Smith if (a->j[j] == i) { 12909f38230SBarry Smith a->diag[i] = j; 130b97cf49bSBarry Smith break; 131b97cf49bSBarry Smith } 132b97cf49bSBarry Smith } 133b97cf49bSBarry Smith } 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 135b97cf49bSBarry Smith } 136b97cf49bSBarry Smith 1374a2ae208SSatish Balay #undef __FUNCT__ 1384a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 139b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 140b97cf49bSBarry Smith { 1413eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 142b24ad042SBarry Smith PetscInt *itmp; 143b97cf49bSBarry Smith 144433994e6SBarry Smith PetscFunctionBegin; 145899cda47SBarry Smith row -= A->rmap.rstart; 1460752156aSBarry Smith 147899cda47SBarry Smith if (row < 0 || row >= A->rmap.n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 148b97cf49bSBarry Smith 149b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 150b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 151b97cf49bSBarry Smith if (idx) { 152b97cf49bSBarry Smith itmp = a->j + a->i[row]; 153b97cf49bSBarry Smith if (*nz) { 154b97cf49bSBarry Smith *idx = itmp; 155b97cf49bSBarry Smith } 156b97cf49bSBarry Smith else *idx = 0; 157b97cf49bSBarry Smith } 1583a40ed3dSBarry Smith PetscFunctionReturn(0); 159b97cf49bSBarry Smith } 160b97cf49bSBarry Smith 1614a2ae208SSatish Balay #undef __FUNCT__ 1624a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 163b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 164b97cf49bSBarry Smith { 165433994e6SBarry Smith PetscFunctionBegin; 1663a40ed3dSBarry Smith PetscFunctionReturn(0); 167b97cf49bSBarry Smith } 168b97cf49bSBarry Smith 1694a2ae208SSatish Balay #undef __FUNCT__ 1704a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 171dfbe8321SBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 172b97cf49bSBarry Smith { 1733eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 174dfbe8321SBarry Smith PetscErrorCode ierr; 1750f5bd95cSBarry Smith PetscTruth flag; 176b97cf49bSBarry Smith 177433994e6SBarry Smith PetscFunctionBegin; 178b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 179899cda47SBarry Smith if ((A->rmap.n != B->rmap.n) ||(a->nz != b->nz)) { 1800f5bd95cSBarry Smith flag = PETSC_FALSE; 181b97cf49bSBarry Smith } 182b97cf49bSBarry Smith 183b97cf49bSBarry Smith /* if the a->i are the same */ 184899cda47SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->rmap.n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 185b97cf49bSBarry Smith 186b97cf49bSBarry Smith /* if a->j are the same */ 187b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 188b97cf49bSBarry Smith 189ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 191b97cf49bSBarry Smith } 192b97cf49bSBarry Smith 1934a2ae208SSatish Balay #undef __FUNCT__ 1944a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 195b24ad042SBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 1966cd91f64SBarry Smith { 197dfbe8321SBarry Smith PetscErrorCode ierr; 198b24ad042SBarry Smith PetscMPIInt size; 199b24ad042SBarry Smith PetscInt i; 2006cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 2016cd91f64SBarry Smith 2026cd91f64SBarry Smith PetscFunctionBegin; 2036cd91f64SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 2046cd91f64SBarry Smith if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 205899cda47SBarry 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" 220b24ad042SBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 221d42735a1SBarry Smith { 222b24ad042SBarry Smith PetscInt i; 223d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 224d42735a1SBarry Smith 225d42735a1SBarry Smith PetscFunctionBegin; 226e005ede5SBarry Smith if (ia && a->i != *ia) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 227e005ede5SBarry Smith if (ja && a->j != *ja) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 228d42735a1SBarry Smith if (oshift) { 229d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 230d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 231d42735a1SBarry Smith (*ja)[i]--; 232d42735a1SBarry Smith } 233d42735a1SBarry Smith } 2346cd91f64SBarry Smith PetscFunctionReturn(0); 2356cd91f64SBarry Smith } 236b97cf49bSBarry Smith 23717667f90SBarry Smith #undef __FUNCT__ 23817667f90SBarry Smith #define __FUNCT__ "MatConvertFrom_MPIAdj" 23917667f90SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 24017667f90SBarry Smith { 24117667f90SBarry Smith Mat B; 24217667f90SBarry Smith PetscErrorCode ierr; 24317667f90SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 24417667f90SBarry Smith const PetscInt *rj; 24517667f90SBarry Smith const PetscScalar *ra; 24617667f90SBarry Smith MPI_Comm comm; 24717667f90SBarry Smith 24817667f90SBarry Smith PetscFunctionBegin; 24917667f90SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 25017667f90SBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 25117667f90SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 25217667f90SBarry Smith 25317667f90SBarry Smith /* count the number of nonzeros per row */ 25417667f90SBarry Smith for (i=0; i<m; i++) { 25517667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 25617667f90SBarry Smith for (j=0; j<len; j++) { 25717667f90SBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 25817667f90SBarry Smith } 25917667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 26017667f90SBarry Smith nzeros += len; 26117667f90SBarry Smith } 26217667f90SBarry Smith 26317667f90SBarry Smith /* malloc space for nonzeros */ 26417667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr); 26517667f90SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr); 26617667f90SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr); 26717667f90SBarry Smith 26817667f90SBarry Smith nzeros = 0; 26917667f90SBarry Smith ia[0] = 0; 27017667f90SBarry Smith for (i=0; i<m; i++) { 27117667f90SBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 27217667f90SBarry Smith cnt = 0; 27317667f90SBarry Smith for (j=0; j<len; j++) { 27417667f90SBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 27517667f90SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 27617667f90SBarry Smith ja[nzeros+cnt++] = rj[j]; 27717667f90SBarry Smith } 27817667f90SBarry Smith } 27917667f90SBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 28017667f90SBarry Smith nzeros += cnt; 28117667f90SBarry Smith ia[i+1] = nzeros; 28217667f90SBarry Smith } 28317667f90SBarry Smith 28417667f90SBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 28517667f90SBarry Smith ierr = MatCreate(comm,&B);CHKERRQ(ierr); 28617667f90SBarry Smith ierr = MatSetSizes(B,m,N,PETSC_DETERMINE,N);CHKERRQ(ierr); 28717667f90SBarry Smith ierr = MatSetType(B,type);CHKERRQ(ierr); 28817667f90SBarry Smith ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 28917667f90SBarry Smith 29017667f90SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 29117667f90SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 29217667f90SBarry Smith } else { 29317667f90SBarry Smith *newmat = B; 29417667f90SBarry Smith } 29517667f90SBarry Smith PetscFunctionReturn(0); 29617667f90SBarry Smith } 29717667f90SBarry Smith 298b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2990b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 3003eda8832SBarry Smith MatGetRow_MPIAdj, 3013eda8832SBarry Smith MatRestoreRow_MPIAdj, 302b97cf49bSBarry Smith 0, 30397304618SKris Buschelman /* 4*/ 0, 304b97cf49bSBarry Smith 0, 305b97cf49bSBarry Smith 0, 306b97cf49bSBarry Smith 0, 3070b3b1487SBarry Smith 0, 3080b3b1487SBarry Smith 0, 30997304618SKris Buschelman /*10*/ 0, 3100b3b1487SBarry Smith 0, 3110b3b1487SBarry Smith 0, 3120b3b1487SBarry Smith 0, 3130b3b1487SBarry Smith 0, 31497304618SKris Buschelman /*15*/ 0, 3153eda8832SBarry Smith MatEqual_MPIAdj, 3160b3b1487SBarry Smith 0, 3170b3b1487SBarry Smith 0, 3180b3b1487SBarry Smith 0, 31997304618SKris Buschelman /*20*/ 0, 3200b3b1487SBarry Smith 0, 3210b3b1487SBarry Smith 0, 3223eda8832SBarry Smith MatSetOption_MPIAdj, 3230b3b1487SBarry Smith 0, 32497304618SKris Buschelman /*25*/ 0, 3250b3b1487SBarry Smith 0, 3260b3b1487SBarry Smith 0, 3270b3b1487SBarry Smith 0, 3280b3b1487SBarry Smith 0, 32997304618SKris Buschelman /*30*/ 0, 3300b3b1487SBarry Smith 0, 331273d9f13SBarry Smith 0, 332273d9f13SBarry Smith 0, 3330b3b1487SBarry Smith 0, 33497304618SKris Buschelman /*35*/ 0, 3350b3b1487SBarry Smith 0, 3360b3b1487SBarry Smith 0, 3370b3b1487SBarry Smith 0, 3380b3b1487SBarry Smith 0, 33997304618SKris Buschelman /*40*/ 0, 3400b3b1487SBarry Smith 0, 3410b3b1487SBarry Smith 0, 3420b3b1487SBarry Smith 0, 3430b3b1487SBarry Smith 0, 34497304618SKris Buschelman /*45*/ 0, 3450b3b1487SBarry Smith 0, 3460b3b1487SBarry Smith 0, 3470b3b1487SBarry Smith 0, 3480b3b1487SBarry Smith 0, 349521d7252SBarry Smith /*50*/ 0, 3506cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 351d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 352b97cf49bSBarry Smith 0, 353b97cf49bSBarry Smith 0, 35497304618SKris Buschelman /*55*/ 0, 355b97cf49bSBarry Smith 0, 3560752156aSBarry Smith 0, 3570752156aSBarry Smith 0, 3580b3b1487SBarry Smith 0, 35997304618SKris Buschelman /*60*/ 0, 360b9b97703SBarry Smith MatDestroy_MPIAdj, 361b9b97703SBarry Smith MatView_MPIAdj, 36217667f90SBarry Smith MatConvertFrom_MPIAdj, 36397304618SKris Buschelman 0, 36497304618SKris Buschelman /*65*/ 0, 36597304618SKris Buschelman 0, 36697304618SKris Buschelman 0, 36797304618SKris Buschelman 0, 36897304618SKris Buschelman 0, 36997304618SKris Buschelman /*70*/ 0, 37097304618SKris Buschelman 0, 37197304618SKris Buschelman 0, 37297304618SKris Buschelman 0, 37397304618SKris Buschelman 0, 37497304618SKris Buschelman /*75*/ 0, 37597304618SKris Buschelman 0, 37697304618SKris Buschelman 0, 37797304618SKris Buschelman 0, 37897304618SKris Buschelman 0, 37997304618SKris Buschelman /*80*/ 0, 38097304618SKris Buschelman 0, 38197304618SKris Buschelman 0, 38297304618SKris Buschelman 0, 383865e5f61SKris Buschelman 0, 384865e5f61SKris Buschelman /*85*/ 0, 385865e5f61SKris Buschelman 0, 386865e5f61SKris Buschelman 0, 387865e5f61SKris Buschelman 0, 388865e5f61SKris Buschelman 0, 389865e5f61SKris Buschelman /*90*/ 0, 390865e5f61SKris Buschelman 0, 391865e5f61SKris Buschelman 0, 392865e5f61SKris Buschelman 0, 393865e5f61SKris Buschelman 0, 394865e5f61SKris Buschelman /*95*/ 0, 395865e5f61SKris Buschelman 0, 396865e5f61SKris Buschelman 0, 397865e5f61SKris Buschelman 0}; 398b97cf49bSBarry Smith 399a23d5eceSKris Buschelman EXTERN_C_BEGIN 400a23d5eceSKris Buschelman #undef __FUNCT__ 401a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 402be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 403a23d5eceSKris Buschelman { 404a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 405dfbe8321SBarry Smith PetscErrorCode ierr; 4062515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 407b24ad042SBarry Smith PetscInt ii; 408a23d5eceSKris Buschelman #endif 409a23d5eceSKris Buschelman 410a23d5eceSKris Buschelman PetscFunctionBegin; 411a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 4122515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 41377431f27SBarry Smith if (i[0] != 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]); 414899cda47SBarry Smith for (ii=1; ii<B->rmap.n; ii++) { 415a23d5eceSKris Buschelman if (i[ii] < 0 || i[ii] < i[ii-1]) { 41677431f27SBarry Smith SETERRQ4(PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]); 417a23d5eceSKris Buschelman } 418a23d5eceSKris Buschelman } 419899cda47SBarry Smith for (ii=0; ii<i[B->rmap.n]; ii++) { 420899cda47SBarry Smith if (j[ii] < 0 || j[ii] >= B->cmap.N) { 42177431f27SBarry Smith SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]); 422a23d5eceSKris Buschelman } 423a23d5eceSKris Buschelman } 424a23d5eceSKris Buschelman #endif 425a23d5eceSKris Buschelman 426a23d5eceSKris Buschelman b->j = j; 427a23d5eceSKris Buschelman b->i = i; 428a23d5eceSKris Buschelman b->values = values; 429a23d5eceSKris Buschelman 430899cda47SBarry Smith b->nz = i[B->rmap.n]; 431a23d5eceSKris Buschelman b->diag = 0; 432a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 433a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 434a23d5eceSKris Buschelman 435a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 436a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 437a23d5eceSKris Buschelman PetscFunctionReturn(0); 438a23d5eceSKris Buschelman } 439a23d5eceSKris Buschelman EXTERN_C_END 440b97cf49bSBarry Smith 4410bad9183SKris Buschelman /*MC 442fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 4430bad9183SKris Buschelman intended for use constructing orderings and partitionings. 4440bad9183SKris Buschelman 4450bad9183SKris Buschelman Level: beginner 4460bad9183SKris Buschelman 4470bad9183SKris Buschelman .seealso: MatCreateMPIAdj 4480bad9183SKris Buschelman M*/ 4490bad9183SKris Buschelman 450273d9f13SBarry Smith EXTERN_C_BEGIN 4514a2ae208SSatish Balay #undef __FUNCT__ 4524a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 453be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIAdj(Mat B) 454273d9f13SBarry Smith { 455273d9f13SBarry Smith Mat_MPIAdj *b; 4566849ba73SBarry Smith PetscErrorCode ierr; 457b24ad042SBarry Smith PetscMPIInt size,rank; 458273d9f13SBarry Smith 459273d9f13SBarry Smith PetscFunctionBegin; 460273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 461273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 462273d9f13SBarry Smith 463b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr); 464b0a32e0cSBarry Smith B->data = (void*)b; 465273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 466273d9f13SBarry Smith B->factor = 0; 467273d9f13SBarry Smith B->mapping = 0; 468273d9f13SBarry Smith B->assembled = PETSC_FALSE; 469273d9f13SBarry Smith 470899cda47SBarry Smith ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr); 471e2ee2a47SBarry Smith if (B->cmap.n < 0) B->cmap.n = B->cmap.N; 472e2ee2a47SBarry Smith if (B->cmap.N < 0) B->cmap.N = B->cmap.n; 473273d9f13SBarry Smith 474a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C", 475a23d5eceSKris Buschelman "MatMPIAdjSetPreallocation_MPIAdj", 476a23d5eceSKris Buschelman MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 47717667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 478273d9f13SBarry Smith PetscFunctionReturn(0); 479273d9f13SBarry Smith } 480273d9f13SBarry Smith EXTERN_C_END 481273d9f13SBarry Smith 4824a2ae208SSatish Balay #undef __FUNCT__ 4834a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 484a23d5eceSKris Buschelman /*@C 485a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 486a23d5eceSKris Buschelman 487a23d5eceSKris Buschelman Collective on MPI_Comm 488a23d5eceSKris Buschelman 489a23d5eceSKris Buschelman Input Parameters: 490a23d5eceSKris Buschelman + A - the matrix 491a23d5eceSKris Buschelman . i - the indices into j for the start of each row 492a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 493a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 494a23d5eceSKris Buschelman - values - [optional] edge weights 495a23d5eceSKris Buschelman 496a23d5eceSKris Buschelman Level: intermediate 497a23d5eceSKris Buschelman 498a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 499a23d5eceSKris Buschelman @*/ 500be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 501273d9f13SBarry Smith { 502b24ad042SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*,PetscInt*); 503273d9f13SBarry Smith 504273d9f13SBarry Smith PetscFunctionBegin; 505a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 506a23d5eceSKris Buschelman if (f) { 507a23d5eceSKris Buschelman ierr = (*f)(B,i,j,values);CHKERRQ(ierr); 508273d9f13SBarry Smith } 509273d9f13SBarry Smith PetscFunctionReturn(0); 510273d9f13SBarry Smith } 511273d9f13SBarry Smith 5124a2ae208SSatish Balay #undef __FUNCT__ 5134a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 514b97cf49bSBarry Smith /*@C 5153eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 516b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 517b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 518b97cf49bSBarry Smith 519ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 520ef5ee4d1SLois Curfman McInnes 521b97cf49bSBarry Smith Input Parameters: 522c2e958feSBarry Smith + comm - MPI communicator 5230752156aSBarry Smith . m - number of local rows 524b97cf49bSBarry Smith . n - number of columns 525b97cf49bSBarry Smith . i - the indices into j for the start of each row 526329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 527ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 528329f5518SBarry Smith - values -[optional] edge weights 529b97cf49bSBarry Smith 530b97cf49bSBarry Smith Output Parameter: 531b97cf49bSBarry Smith . A - the matrix 532b97cf49bSBarry Smith 53315091d37SBarry Smith Level: intermediate 53415091d37SBarry Smith 5354bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 5364bc6d8c0SBarry Smith MatSetValues(). 537c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 538ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 5391198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 540327e1a73SBarry Smith Should not include the matrix diagonals. 541b97cf49bSBarry Smith 542e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 543ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 544ededeb1aSvictorle 545ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 546b97cf49bSBarry Smith 547e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 548b97cf49bSBarry Smith @*/ 549be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 550b97cf49bSBarry Smith { 551dfbe8321SBarry Smith PetscErrorCode ierr; 552b97cf49bSBarry Smith 553433994e6SBarry Smith PetscFunctionBegin; 554f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 555e2ee2a47SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,n);CHKERRQ(ierr); 556273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 557273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 5583a40ed3dSBarry Smith PetscFunctionReturn(0); 559b97cf49bSBarry Smith } 560b97cf49bSBarry Smith 561c2e958feSBarry Smith 562433994e6SBarry Smith 563433994e6SBarry Smith 564433994e6SBarry Smith 565433994e6SBarry Smith 566433994e6SBarry Smith 567